da_get_innov_vector_synop.inc

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