1 subroutine da_transform_xtoy_buoy (grid, iv, y)
3 !--------------------------------------------------------------------------
5 ! Updated for Analysis on Arakawa-C grid
6 ! Author: Syed RH Rizvi, MMM/ESSL/NCAR, Date: 10/22/2008
7 !--------------------------------------------------------------------------
11 type (domain), intent(inout) :: grid
12 type (iv_type), intent(in) :: iv ! Innovation vector (O-B).
13 type (y_type), intent(inout) :: y ! y = h (grid%xa) (linear)
15 integer :: n ! Loop counter.
17 real, allocatable :: model_u(:,:)
18 real, allocatable :: model_v(:,:)
19 real, allocatable :: model_t(:,:)
20 real, allocatable :: model_q(:,:)
21 real, allocatable :: model_psfc(:)
23 real, allocatable :: ub(:,:)
24 real, allocatable :: vb(:,:)
26 if (trace_use_dull) call da_trace_entry("da_transform_xtoy_buoy")
28 if (sfc_assi_options == sfc_assi_options_1) then
29 allocate (model_u(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
30 allocate (model_v(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
31 allocate (model_t(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
32 allocate (model_q(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
33 allocate (model_psfc(iv%info(buoy)%n1:iv%info(buoy)%n2))
35 allocate (ub(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
36 allocate (vb(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
38 ! [1.2] Interpolate horizontally:
40 call da_interp_lin_3d (grid%xa%u, iv%info(buoy), model_u,'u')
41 call da_interp_lin_3d (grid%xa%v, iv%info(buoy), model_v,'v')
43 call da_interp_lin_3d (grid%xa%u, iv%info(buoy), model_u)
44 call da_interp_lin_3d (grid%xa%v, iv%info(buoy), model_v)
46 call da_interp_lin_3d (grid%xa%t, iv%info(buoy), model_t)
47 call da_interp_lin_3d (grid%xa%q, iv%info(buoy), model_q)
49 call da_interp_lin_2d (grid%xa%psfc, iv%info(buoy), 1, model_psfc)
51 call da_interp_lin_3d (grid%xb%u, iv%info(buoy), ub)
52 call da_interp_lin_3d (grid%xb%v, iv%info(buoy), vb)
54 do n=iv%info(buoy)%n1,iv%info(buoy)%n2
56 call da_uv_to_sd_lin(y%buoy(n)%u,y%buoy(n)%v,model_u(1,n),model_v(1,n),ub(1,n),vb(1,n))
58 y%buoy(n)%u = model_u(1,n)
59 y%buoy(n)%v = model_v(1,n)
61 y%buoy(n)%t = model_t(1,n)
62 y%buoy(n)%q = model_q(1,n)
63 y%buoy(n)%p = model_psfc(n)
69 deallocate (model_psfc)
73 else if (sfc_assi_options == sfc_assi_options_2) then
74 ! [2.0] Surface assimilation approach 2
75 call da_transform_xtopsfc(grid,iv,buoy,iv%buoy(:),y%buoy(:))
78 if (trace_use_dull) call da_trace_exit("da_transform_xtoy_buoy")
80 end subroutine da_transform_xtoy_buoy