da_get_innov_vector_sonde_sfc.inc
References to this file elsewhere.
1 subroutine da_get_innov_vector_sonde_sfc( it, xb, xp, ob, iv)
2
3 !-----------------------------------------------------------------------
4 ! Purpose: TBD
5 ! Update :
6 ! 01/24/2007 Syed RH Rizvi
7 ! Updated for "VERIFY"
8 !-----------------------------------------------------------------------
9
10 implicit none
11
12 integer, intent(in) :: it ! External iteration.
13 type(xb_type), intent(in) :: xb ! first guess state.
14 type(xpose_type), intent(in) :: xp ! Domain decomposition vars.
15 type(y_type), intent(inout) :: ob ! Observation structure.
16 type(ob_type), intent(inout) :: iv ! O-B structure.
17
18 integer :: n ! Loop counter.
19 integer :: i, j, k ! Index dimension.
20 real :: dx, dxm ! Interpolation weights.
21 real :: dy, dym ! Interpolation weights.
22 real :: model_u ! Model value u at oblocation.
23 real :: model_v ! Model value v at oblocation.
24 real :: model_t ! Model value t at oblocation.
25 real :: model_p ! Model value p at oblocation.
26 real :: model_q ! Model value q at oblocation.
27
28 real, dimension(xp%kms:xp%kme) :: v_h ! Model value h at ob hor. location.
29 real, dimension(xp%kms:xp%kme) :: v_p ! Model value p at ob hor. location.
30
31 real :: hd, psfcm
32
33 real :: hsm , ho, to, qo
34
35 integer :: itu,ituf,itvv,itvvf,itp,itpf,itt,ittf,itqv,itqvf
36
37 if (iv % num_sound < 1) return
38
39 itu = 0; itvv = 0; itp = 0; itt = 0; itqv = 0;
40 ituf = 0; itvvf = 0; itpf = 0; ittf = 0; itqvf = 0;
41
42 do n=iv%ob_numb(iv%current_ob_time-1)%sound + 1, &
43 iv%ob_numb(iv%current_ob_time)%sound
44
45 ! [1.1] Get horizontal interpolation weights:
46
47 i = iv%sonde_sfc(n)%loc%i
48 j = iv%sonde_sfc(n)%loc%j
49 dx = iv%sonde_sfc(n)%loc%dx
50 dy = iv%sonde_sfc(n)%loc%dy
51 dxm = iv%sonde_sfc(n)%loc%dxm
52 dym = iv%sonde_sfc(n)%loc%dym
53
54 !-----------------Surface correction
55
56 iv%sonde_sfc(n)%p%inv = ob%sonde_sfc(n)%p
57 iv%sonde_sfc(n)%t%inv = ob%sonde_sfc(n)%t
58 iv%sonde_sfc(n)%q%inv = ob%sonde_sfc(n)%q
59 iv%sonde_sfc(n)%u%inv = ob%sonde_sfc(n)%u
60 iv%sonde_sfc(n)%v%inv = ob%sonde_sfc(n)%v
61
62 if (sfc_assi_options == 1) then
63 iv%sonde_sfc(n)%zk=missing_r
64
65 if (iv % sonde_sfc(n) % h > missing_r) then
66 do k=xp%kts,xp%kte
67 v_h(k) = dym*(dxm*xb%h(i,j ,k) + dx*xb%h(i+1,j ,k)) &
68 + dy *(dxm*xb%h(i,j+1,k) + dx*xb%h(i+1,j+1,k))
69 end do
70
71 hd = v_h(xp%kts) - iv % sonde_sfc(n) % h
72
73 if (abs(hd) <= Max_StHeight_Diff .or. anal_type_verify ) then
74 if (iv % sonde_sfc(n) % h < v_h(xp%kts)) then
75 iv%sonde_sfc(n)%zk = 1.0+1.0e-6
76 call da_obs_sfc_correction(iv%sonde_sfc(n), xb)
77
78 ! To keep the original "ob" with no change for multiple
79 ! outer-loops use:
80
81 ! ob%sonde_sfc(n)%p = iv%sonde_sfc(n)%p%inv
82 ! ob%sonde_sfc(n)%t = iv%sonde_sfc(n)%t%inv
83 ! ob%sonde_sfc(n)%q = iv%sonde_sfc(n)%q%inv
84 ! ob%sonde_sfc(n)%u = iv%sonde_sfc(n)%u%inv
85 ! ob%sonde_sfc(n)%v = iv%sonde_sfc(n)%v%inv
86 else
87 call da_to_zk(iv % sonde_sfc(n) % h, v_h, xp, v_interp_h, &
88 iv%sonde_sfc(n)%zk)
89 end if
90 else
91 iv%sonde_sfc(n)%zk = missing_r
92 end if
93 else if (ob % sonde_sfc(n) % p > 1.0) then
94 do k=xp%kts,xp%kte
95 v_p(k) = dym*(dxm*xb%p(i,j ,k) + dx*xb%p(i+1,j ,k)) &
96 + dy *(dxm*xb%p(i,j+1,k) + dx*xb%p(i+1,j+1,k))
97 end do
98
99 call da_to_zk(ob % sonde_sfc(n) % p, v_p, xp, v_interp_p, iv%sonde_sfc(n)%zk)
100
101 if (iv%sonde_sfc(n)%zk < 0.0 .and. .not.anal_type_verify ) then
102 iv % sonde_sfc(n) % p % inv = missing_r
103 iv % sonde_sfc(n) % p % qc = missing
104 iv%sonde_sfc(n)%zk = 1.0+1.0e-6
105 end if
106 end if
107
108 !------------------------------------------------------------------------
109 ! [2.0] Initialise components of innovation vector:
110 !------------------------------------------------------------------------
111
112 if (iv%sonde_sfc(n)%zk < 0.0 .and. .not.anal_type_verify ) then
113 iv % sonde_sfc(n) % u % qc = missing
114 iv % sonde_sfc(n) % v % qc = missing
115 iv % sonde_sfc(n) % t % qc = missing
116 iv % sonde_sfc(n) % q % qc = missing
117 iv % sonde_sfc(n) % p % qc = missing
118 else
119 !------------[1.2] Interpolate horizontally:
120 call da_interp_obs_lin_2d( xb % u, xp, i, j, dx, dy, dxm, dym, &
121 model_u, iv%sonde_sfc(n)%zk)
122 call da_interp_obs_lin_2d( xb % v, xp, i, j, dx, dy, dxm, dym, &
123 model_v, iv%sonde_sfc(n)%zk)
124 call da_interp_obs_lin_2d( xb % t, xp, i, j, dx, dy, dxm, dym, &
125 model_t, iv%sonde_sfc(n)%zk)
126 call da_interp_obs_lin_2d( xb % q, xp, i, j, dx, dy, dxm, dym, &
127 model_q, iv%sonde_sfc(n)%zk)
128 call da_interp_obs_lin_2d( xb % p, xp, i, j, dx, dy, dxm, dym, &
129 model_p, iv%sonde_sfc(n)%zk)
130 end if
131 else if (sfc_assi_options == 2) then
132
133 !------- Surface data assimilation approca 2
134 !-------------------------------------------
135
136 ! 1.2.1 Surface assimilation approach 2(10-m u, v, 2-m t, q, and sfc_p)
137
138 call da_interp_lin_2d( xb % u10, xp%ims, xp%ime, xp%jms, xp%jme, &
139 i, j, dx, dy, dxm, dym, &
140 model_u)
141 call da_interp_lin_2d( xb % v10, xp%ims, xp%ime, xp%jms, xp%jme, &
142 i, j, dx, dy, dxm, dym, &
143 model_v)
144 call da_interp_lin_2d( xb % t2, xp%ims, xp%ime, xp%jms, xp%jme, &
145 i, j, dx, dy, dxm, dym, &
146 model_t)
147 call da_interp_lin_2d( xb % q2, xp%ims, xp%ime, xp%jms, xp%jme, &
148 i, j, dx, dy, dxm, dym, &
149 model_q)
150 call da_interp_lin_2d( xb % psfc, xp%ims, xp%ime, xp%jms, xp%jme, &
151 i, j, dx, dy, dxm, dym, &
152 model_p)
153
154 if (iv%sonde_sfc(n)%p%qc >= 0) then
155
156 !.......... model surface p, t, q, h at observed site:
157
158 call da_interp_lin_2d( xb % terr, xp%ims, xp%ime, xp%jms, xp%jme, &
159 i, j, dx, dy, dxm, dym, &
160 hsm)
161
162 ho = iv%sonde_sfc(n)%h
163 to = -888888.0
164 qo = -888888.0
165
166 if (iv%sonde_sfc(n)%t%qc >= 0 .and. &
167 iv%sonde_sfc(n)%q%qc >= 0) then
168 to = ob%sonde_sfc(n)%t
169 qo = ob%sonde_sfc(n)%q
170 call da_sfc_pre(psfcm, model_p, model_t, model_q, &
171 hsm, ho, to, qo)
172 else if (iv%sonde_sfc(n)%t%qc >= 0 .and. &
173 iv%sonde_sfc(n)%q%qc < 0) then
174 to = ob%sonde_sfc(n)%t
175 call da_sfc_pre(psfcm, model_p, model_t, model_q, hsm, ho,to)
176 else
177 call da_sfc_pre(psfcm, model_p, model_t, model_q, hsm, ho)
178 end if
179
180 !.......... Pressure at the observed height:
181 model_p = psfcm
182 end if
183 end if
184
185 !-----------------------------------------------------------------------
186 ! [3.0] Fast interpolation:
187 !-----------------------------------------------------------------------
188
189 if (ob % sonde_sfc(n) % u > missing_r .AND. &
190 iv % sonde_sfc(n) % u % qc >= obs_qc_pointer) then
191 iv % sonde_sfc(n) % u % inv = iv%sonde_sfc(n)%u%inv - model_u
192 else
193 iv % sonde_sfc(n) % u % inv = 0.0
194 end if
195
196 if (ob % sonde_sfc(n) % v > missing_r .AND. &
197 iv % sonde_sfc(n) % v % qc >= obs_qc_pointer) then
198 iv % sonde_sfc(n) % v % inv = iv%sonde_sfc(n)%v%inv - model_v
199 else
200 iv % sonde_sfc(n) % v % inv = 0.0
201 end if
202
203 if (ob % sonde_sfc(n) % p > 0.0 .AND. &
204 iv % sonde_sfc(n) % p % qc >= obs_qc_pointer) then
205 iv % sonde_sfc(n) % p % inv = iv%sonde_sfc(n)%p%inv - model_p
206 else
207 iv % sonde_sfc(n) % p % inv = 0.0
208 end if
209
210 if (ob % sonde_sfc(n) % t > 0.0 .AND. &
211 iv % sonde_sfc(n) % t % qc >= obs_qc_pointer) then
212 iv % sonde_sfc(n) % t % inv = iv%sonde_sfc(n)%t%inv - model_t
213 else
214 iv % sonde_sfc(n) % t % inv = 0.0
215 end if
216
217 if (ob % sonde_sfc(n) % q > 0.0 .AND. &
218 iv % sonde_sfc(n) % q % qc >= obs_qc_pointer) then
219 iv % sonde_sfc(n) % q % inv = iv%sonde_sfc(n)%q%inv - model_q
220 else
221 iv % sonde_sfc(n) % q % inv = 0.0
222 end if
223
224 !-----------------------------------------------------------------------
225 ! [5.0] Perform optional maximum error check:
226 !-----------------------------------------------------------------------
227
228 if (check_max_iv) then
229 call da_check_max_iv_sonde_sfc(it, iv % sonde_sfc(n), &
230 itu,ituf,itvv,itvvf,itp,itpf,itt,ittf,itqv,itqvf)
231 end if
232 end do
233
234 if (rootproc .and. check_max_iv_print) then
235 write(unit= check_max_iv_unit, fmt ='(A,i5,A)')'For outer iteration ',it, &
236 ', Total Rejections for Sonde_Sfc follows:'
237
238 write(unit = check_max_iv_unit, fmt = '(/,10(2(A,I6),/))') &
239 'Number of failed u-wind observations: ',ituf, ' on ',itu, &
240 'Number of failed v-wind observations: ',itvvf,' on ',itvv, &
241 'Number of failed pressure observations: ',itpf, ' on ',itp, &
242 'Number of failed temperature observations:',ittf, ' on ',itt, &
243 'Number of failed mixing ratio observations:',itqvf,' on ',itqv, &
244 'Finally Total Sonde_Sfc rejections ',ituf+itvvf+itpf+ittf+itqvf,' on ',&
245 itu +itvv +itp +itt +itqv
246 end if
247
248 end subroutine da_get_innov_vector_sonde_sfc
249
250