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