1 subroutine da_transform_xtoy_bogus (grid, iv, y)
3 !---------------------------------------------------------------------
4 ! Purpose: the linearized bogus observation operators.
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(in) :: 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.
16 integer :: i, j, k ! Index dimension.
17 integer :: num_levs ! Number of obs levels.
20 real :: model_t(kts:kte)
21 real :: model_t9(kts:kte)
22 real :: model_q(kts:kte)
23 real :: model_q9(kts:kte)
24 real :: model_p_c(kts:kte)
25 real :: model_p_c9(kts:kte)
26 real :: psfc_loc,terr_loc,psfc_loc9
28 real, allocatable :: u(:,:)
29 real, allocatable :: v(:,:)
30 real, allocatable :: t(:,:)
31 real, allocatable :: q(:,:)
33 if (trace_use_dull) call da_trace_entry("da_transform_xtoy_bogus")
35 do n = iv%info(bogus)%n1,iv%info(bogus)%n2
36 num_levs = iv%info(bogus)%levels(n)
38 if (num_levs < 1) cycle
40 ! [1.1] Get cross pt. horizontal interpolation weights:
42 i = iv%info(bogus)%i(1,n)
43 dy = iv%info(bogus)%dy(1,n)
44 dym = iv%info(bogus)%dym(1,n)
45 j = iv%info(bogus)%j(1,n)
46 dx = iv%info(bogus)%dx(1,n)
47 dxm = iv%info(bogus)%dxm(1,n)
49 ! [1.2] Interpolate horizontally from cross points:
52 model_t9(k) = dym*(dxm*grid%xa%t(i,j,k) + dx*grid%xa%t(i+1,j,k)) &
53 + dy *(dxm*grid%xa%t(i,j+1,k) + dx*grid%xa%t(i+1,j+1,k))
54 model_t(k) = dym*(dxm*grid%xb%t(i,j,k) + dx*grid%xb%t(i+1,j,k)) &
55 + dy *(dxm*grid%xb%t(i,j+1,k) + dx*grid%xb%t(i+1,j+1,k))
56 model_q9(k) = dym*(dxm*grid%xa%q(i,j,k) + dx*grid%xa%q(i+1,j,k)) &
57 + dy *(dxm*grid%xa%q(i,j+1,k) + dx*grid%xa%q(i+1,j+1,k))
58 model_q(k) = dym*(dxm*grid%xb%q(i,j,k) + dx*grid%xb%q(i+1,j,k)) &
59 + dy *(dxm*grid%xb%q(i,j+1,k) + dx*grid%xb%q(i+1,j+1,k))
60 model_p_c9(k) = dym*(dxm*grid%xa%p(i,j,k) + dx*grid%xa%p(i+1,j,k)) &
61 + dy *(dxm*grid%xa%p(i,j+1,k) + dx*grid%xa%p(i+1,j+1,k))
62 model_p_c(k) = dym*(dxm*grid%xb%p(i,j,k) + dx*grid%xb%p(i+1,j,k)) &
63 + dy *(dxm*grid%xb%p(i,j+1,k) + dx*grid%xb%p(i+1,j+1,k))
66 terr_loc = dym*(dxm*grid%xb%terr(i,j) + dx*grid%xb%terr(i+1,j)) &
67 + dy *(dxm*grid%xb%terr(i,j+1) + dx*grid%xb%terr(i+1,j+1))
68 psfc_loc = dym*(dxm*grid%xb%psfc(i,j) + dx*grid%xb%psfc(i+1,j)) &
69 + dy *(dxm*grid%xb%psfc(i,j+1) + dx*grid%xb%psfc(i+1,j+1))
70 psfc_loc9 = dym*(dxm*grid%xa%psfc(i,j) + dx*grid%xa%psfc(i+1,j)) &
71 + dy *(dxm*grid%xa%psfc(i,j+1) + dx*grid%xa%psfc(i+1,j+1))
73 call da_tpq_to_slp_lin (model_t, model_q, model_p_c, terr_loc, psfc_loc, &
74 model_t9, model_q9, model_p_c9, psfc_loc9, y%bogus(n)%slp)
77 allocate (u(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
78 allocate (v(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
79 allocate (t(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
80 allocate (q(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
87 call da_interp_lin_3d (grid%xa%u, iv%info(bogus), u,'u')
88 call da_interp_lin_3d (grid%xa%v, iv%info(bogus), v,'v')
90 call da_interp_lin_3d (grid%xa%u, iv%info(bogus), u)
91 call da_interp_lin_3d (grid%xa%v, iv%info(bogus), v)
93 call da_interp_lin_3d (grid%xa%t, iv%info(bogus), t)
94 call da_interp_lin_3d (grid%xa%q, iv%info(bogus), q)
96 do n=iv%info(bogus)%n1,iv%info(bogus)%n2
97 y%bogus(n)%u(:) = u(1:size(y%bogus(n)%u),n)
98 y%bogus(n)%v(:) = v(1:size(y%bogus(n)%v),n)
99 y%bogus(n)%t(:) = t(1:size(y%bogus(n)%t),n)
100 y%bogus(n)%q(:) = q(1:size(y%bogus(n)%q),n)
108 if (trace_use_dull) call da_trace_exit("da_transform_xtoy_bogus")
110 end subroutine da_transform_xtoy_bogus