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