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