1 subroutine da_get_innov_vector_buoy( it,num_qcstat_conv, grid, ob, iv)
3 !-----------------------------------------------------------------------
5 ! Updated for Analysis on Arakawa-C grid
6 ! Author: Syed RH Rizvi, MMM/ESSL/NCAR, Date: 10/22/2008
7 !-----------------------------------------------------------------------
11 integer, intent(in) :: it ! External iteration.
12 type(domain), intent(in) :: grid ! first guess state.
13 type(y_type), intent(inout) :: ob ! Observation structure.
14 type(iv_type), intent(inout) :: iv ! O-B structure.
15 integer, intent(inout) :: num_qcstat_conv(:,:,:,:)
17 integer :: n ! Loop counter.
18 integer :: i, j, k ! Index dimension.
19 real :: dx, dxm ! Interpolation weights.
20 real :: dy, dym ! Interpolation weights.
21 real, allocatable :: model_u(:,:) ! Model value u at oblocation.
22 real, allocatable :: model_v(:,:) ! Model value v at oblocation.
23 real, allocatable :: model_t(:,:) ! Model value t at oblocation.
24 real, allocatable :: model_p(:,:) ! Model value p at oblocation.
25 real, allocatable :: model_q(:,:) ! Model value q at oblocation.
26 real, allocatable :: model_hsm(:,:)
28 real :: v_h(kms:kme) ! Model value h at ob hor. location.
29 real :: v_p(kms:kme) ! Model value p at ob hor. location.
35 if (trace_use_dull) call da_trace_entry("da_get_innov_vector_buoy")
37 allocate (model_u(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
38 allocate (model_v(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
39 allocate (model_t(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
40 allocate (model_p(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
41 allocate (model_q(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
42 allocate (model_hsm(1,iv%info(buoy)%n1:iv%info(buoy)%n2))
45 do n=iv%info(buoy)%n1,iv%info(buoy)%n2
46 if (iv%buoy(n)%u%qc == fails_error_max) iv%buoy(n)%u%qc = 0
47 if (iv%buoy(n)%v%qc == fails_error_max) iv%buoy(n)%v%qc = 0
48 if (iv%buoy(n)%t%qc == fails_error_max) iv%buoy(n)%t%qc = 0
49 if (iv%buoy(n)%p%qc == fails_error_max) iv%buoy(n)%p%qc = 0
50 if (iv%buoy(n)%q%qc == fails_error_max) iv%buoy(n)%q%qc = 0
54 if (sfc_assi_options == sfc_assi_options_1) then
55 do n=iv%info(buoy)%n1,iv%info(buoy)%n2
56 ! [1.1] Get horizontal interpolation weights:
58 i = iv%info(buoy)%i(1,n)
59 j = iv%info(buoy)%j(1,n)
60 dx = iv%info(buoy)%dx(1,n)
61 dy = iv%info(buoy)%dy(1,n)
62 dxm = iv%info(buoy)%dxm(1,n)
63 dym = iv%info(buoy)%dym(1,n)
67 iv%buoy(n)%p%inv = ob%buoy(n)%p
68 iv%buoy(n)%t%inv = ob%buoy(n)%t
69 iv%buoy(n)%q%inv = ob%buoy(n)%q
70 iv%buoy(n)%u%inv = ob%buoy(n)%u
71 iv%buoy(n)%v%inv = ob%buoy(n)%v
73 if (iv % buoy(n) % h > missing_r) then
75 v_h(k) = dym*(dxm*grid%xb%h(i,j ,k) + dx*grid%xb%h(i+1,j ,k)) &
76 + dy *(dxm*grid%xb%h(i,j+1,k) + dx*grid%xb%h(i+1,j+1,k))
79 hd = v_h(kts) - iv % buoy(n) % h
80 if (abs(hd) <= Max_StHeight_Diff .or. anal_type_verify) then
81 if (iv % buoy(n) % h < v_h(kts)) then
82 iv%info(buoy)%zk(:,n) = 1.0+1.0e-6
83 call da_obs_sfc_correction(iv%info(buoy), iv%buoy(n), n, grid%xb)
86 call da_to_zk(iv % buoy(n) % h, v_h, v_interp_h, iv%info(buoy)%zk(1,n))
89 else if (ob % buoy(n) % p > 1.0) then
91 v_p(k) = dym*(dxm*grid%xb%p(i,j ,k) + dx*grid%xb%p(i+1,j ,k)) &
92 + dy *(dxm*grid%xb%p(i,j+1,k) + dx*grid%xb%p(i+1,j+1,k))
95 call da_to_zk(ob % buoy(n) % p, v_p, v_interp_p, iv%info(buoy)%zk(1,n))
97 if (iv%info(buoy)%zk(1,n) < 0.0 .and. .not.anal_type_verify) then
98 iv % buoy(n) % p % inv = missing_r
99 iv % buoy(n) % p % qc = missing_data
100 iv%info(buoy)%zk(:,n) = 1.0+1.0e-6
105 call da_convert_zk (iv%info(buoy))
107 if (.not.anal_type_verify) then
108 do n=iv%info(buoy)%n1,iv%info(buoy)%n2
109 if (iv%info(buoy)%zk(1,n) < 0.0) then
110 iv % buoy(n) % u % qc = missing_data
111 iv % buoy(n) % v % qc = missing_data
112 iv % buoy(n) % t % qc = missing_data
113 iv % buoy(n) % q % qc = missing_data
114 iv % buoy(n) % p % qc = missing_data
119 ! [1.2] Interpolate horizontally:
121 call da_interp_lin_3d (grid%xb%u, iv%info(buoy),model_u,'u')
122 call da_interp_lin_3d (grid%xb%v, iv%info(buoy),model_v,'v')
124 call da_interp_lin_3d (grid%xb%u, iv%info(buoy),model_u)
125 call da_interp_lin_3d (grid%xb%v, iv%info(buoy),model_v)
127 call da_interp_lin_3d (grid%xb%t, iv%info(buoy),model_t)
128 call da_interp_lin_3d (grid%xb%q, iv%info(buoy),model_q)
129 call da_interp_lin_3d (grid%xb%p, iv%info(buoy),model_p)
130 else if (sfc_assi_options == sfc_assi_options_2) then
132 ! Surface data assimilation approach 2
133 ! -----------------------------------
135 ! 1.2.1 Surface assmiilation approach 2(10-m u, v, 2-m t, q,
138 call da_interp_lin_2d (grid%xb%u10, iv%info(buoy), 1,model_u)
139 call da_interp_lin_2d (grid%xb%v10, iv%info(buoy), 1,model_v)
140 call da_interp_lin_2d (grid%xb%t2, iv%info(buoy), 1,model_t)
141 call da_interp_lin_2d (grid%xb%q2, iv%info(buoy), 1,model_q)
142 call da_interp_lin_2d (grid%xb%psfc, iv%info(buoy), 1,model_p)
144 do n=iv%info(buoy)%n1,iv%info(buoy)%n2
146 iv%buoy(n)%p%inv = ob%buoy(n)%p
147 iv%buoy(n)%t%inv = ob%buoy(n)%t
148 iv%buoy(n)%q%inv = ob%buoy(n)%q
149 iv%buoy(n)%u%inv = ob%buoy(n)%u
150 iv%buoy(n)%v%inv = ob%buoy(n)%v
152 if (iv%buoy(n)%p%qc >= 0) then
153 ! model surface p, t, q, h at observed site:
155 call da_interp_lin_2d_partial (grid%xb%terr, iv%info(buoy), 1, n, n, model_hsm(:,n))
161 if (iv%buoy(n)%t%qc >= 0 .and. iv%buoy(n)%q%qc >= 0) then
164 call da_sfc_pre(psfcm, model_p(1,n), model_t(1,n), model_q(1,n), model_hsm(1,n), ho, to, qo)
165 else if (iv%buoy(n)%t%qc >= 0 .and. iv%buoy(n)%q%qc < 0) then
167 call da_sfc_pre(psfcm, model_p(1,n), model_t(1,n), model_q(1,n), model_hsm(1,n), ho, to)
169 call da_sfc_pre(psfcm, model_p(1,n), model_t(1,n), model_q(1,n), model_hsm(1,n), ho)
172 ! Pressure at the observed height:
178 do n=iv%info(buoy)%n1,iv%info(buoy)%n2
180 !-----------------------------------------------------------------------
181 ! [3.0] Fast interpolation:
182 !-----------------------------------------------------------------------
183 if (ob % buoy(n) % u > missing_r .AND. iv % buoy(n) % u % qc >= obs_qc_pointer) then
184 iv % buoy(n) % u % inv = iv%buoy(n)%u%inv - model_u(1,n)
186 iv % buoy(n) % u % inv = 0.0
189 if (ob % buoy(n) % v > missing_r .AND. iv % buoy(n) % v % qc >= obs_qc_pointer) then
190 iv % buoy(n) % v % inv = iv%buoy(n)%v%inv - model_v(1,n)
192 iv % buoy(n) % v % inv = 0.0
195 if (ob % buoy(n) % p > 0.0 .AND. iv % buoy(n) % p % qc >= obs_qc_pointer) then
196 iv % buoy(n) % p % inv = iv%buoy(n)%p%inv - model_p(1,n)
198 iv % buoy(n) % p % inv = 0.0
201 if (ob % buoy(n) % t > 0.0 .AND. iv % buoy(n) % t % qc >= obs_qc_pointer) then
202 iv % buoy(n) % t % inv = iv%buoy(n)%t%inv - model_t(1,n)
204 iv % buoy(n) % t % inv = 0.0
207 if (ob % buoy(n) % q > 0.0 .AND. iv % buoy(n) % q % qc >= obs_qc_pointer) then
208 iv % buoy(n) % q % inv = iv%buoy(n)%q%inv - model_q(1,n)
210 iv % buoy(n) % q % inv = 0.0
214 ! -----------------------------------------------------------------------
215 ! [5.0] Perform optional maximum error check:
216 !-----------------------------------------------------------------------
218 if ( check_max_iv ) &
219 call da_check_max_iv_buoy(iv,ob, it, num_qcstat_conv)
226 deallocate (model_hsm)
228 if (trace_use_dull) call da_trace_exit("da_get_innov_vector_buoy")
230 end subroutine da_get_innov_vector_buoy