da_get_innov_vector_bogus.inc

References to this file elsewhere.
1 subroutine da_get_innov_vector_bogus( it, xb, xp, ob, iv)
2 
3    !------------------------------------------------------------------------------
4    ! Purpose: calculate the innovations for the bogus data.
5    !------------------------------------------------------------------------------
6 
7    implicit none
8 
9    integer, intent(in)            :: it       ! External iteration.
10    type(xb_type), intent(in)     :: xb       ! first guess state.
11    type(xpose_type), intent(in)  :: xp       ! Domain decomposition vars.
12    type(y_type),  intent(in)     :: ob       ! Observation structure.
13    type(ob_type), intent(inout)  :: iv       ! O-B structure.
14 
15    integer                        :: n        ! Loop counter.
16    integer                        :: i, j, k  ! Index dimension.
17    integer                        :: num_levs ! Number of obs levels.
18 
19    real                           :: dx, dxm  ! Interpolation weights.
20    real                           :: dy, dym  ! Interpolation weights.
21    real, dimension(xp%kms:xp%kme) :: v_h      ! Model value h at ob hor. location.
22    real, dimension(xp%kms:xp%kme) :: v_p      ! Model value p at ob hor. location.
23 
24    real, dimension(1:max_ob_levels) :: model_u   ! Model value u at ob location.
25    real, dimension(1:max_ob_levels) :: model_v   ! Model value v at ob location.
26    real, dimension(1:max_ob_levels) :: model_t   ! Model value t at ob location.
27    real, dimension(1:max_ob_levels) :: model_q   ! Model value q at ob location.
28    real                             :: model_slp ! Model value slp at ob location.
29 
30    integer           :: itu,ituf,itvv,itvvf,itt,ittf,itqv,itqvf,itslp,itslpf
31 
32    if (iv % num_bogus < 1) return
33 
34    itu   = 0; itvv    = 0; itt  = 0; itqv  = 0; itslp  = 0;
35    ituf  = 0; itvvf   = 0; ittf = 0; itqvf = 0; itslpf = 0;
36 
37    do n=iv%ob_numb(iv%current_ob_time-1)%bogus+1,iv%ob_numb(iv%current_ob_time)%bogus
38       num_levs = iv % bogus(n) % info % levels
39 
40       if (num_levs < 1) cycle
41 
42       model_u(:) = 0.0
43       model_v(:) = 0.0
44       model_t(:) = 0.0
45       model_q(:) = 0.0
46       model_slp  = 0.0
47 
48       i = iv%bogus(n)%loc%i
49       j = iv%bogus(n)%loc%j
50       dx = iv%bogus(n)%loc%dx
51       dy = iv%bogus(n)%loc%dy
52       dxm = iv%bogus(n)%loc%dxm
53       dym = iv%bogus(n)%loc%dym
54 
55       do k=xp%kts,xp%kte
56          v_h(k) = dym*(dxm*xb%h(i,j  ,k) + dx*xb%h(i+1,j  ,k)) &
57                 + dy *(dxm*xb%h(i,j+1,k) + dx*xb%h(i+1,j+1,k))
58          v_p(k) = dym*(dxm*xb%p(i,j  ,k) + dx*xb%p(i+1,j  ,k)) &
59                 + dy *(dxm*xb%p(i,j+1,k) + dx*xb%p(i+1,j+1,k))
60       end do
61 
62       do k=1, iv % bogus(n) % info % levels
63         iv%bogus(n)%zk(k)=missing_r
64         if (iv % bogus(n) % p(k) > 1.0) then
65           call da_to_zk(iv % bogus(n) % p(k), v_p, xp, v_interp_p, iv%bogus(n)%zk(k))
66         else if (iv % bogus(n) % h(k) > 0.0) then
67           call da_to_zk(iv % bogus(n) % h(k), v_h, xp, v_interp_h, iv%bogus(n)%zk(k))
68         end if
69 
70         if (iv%bogus(n)%zk(k) < 0.0 .and.  .not.anal_type_verify) then
71           iv % bogus(n) % u(k) % qc = missing
72           iv % bogus(n) % v(k) % qc = missing
73           iv % bogus(n) % t(k) % qc = missing
74           iv % bogus(n) % q(k) % qc = missing
75         end if
76       end do
77 
78       ! [1.4] Interpolate horizontally:
79       call da_interp_lin_3d( xb % u, xp, i, j, dx, dy, dxm, dym, &
80                          model_u, max_ob_levels, iv%bogus(n)%zk, num_levs)
81       call da_interp_lin_3d( xb % v, xp, i, j, dx, dy, dxm, dym, &
82                          model_v, max_ob_levels, iv%bogus(n)%zk, num_levs)
83       call da_interp_lin_3d( xb % t, xp, i, j, dx, dy, dxm, dym, &
84                          model_t, max_ob_levels, iv%bogus(n)%zk, num_levs)
85       call da_interp_lin_3d( xb % q, xp, i, j, dx, dy, dxm, dym, &
86                          model_q, max_ob_levels, iv%bogus(n)%zk, num_levs)
87 
88       model_slp = dym*(dxm*xb%slp(i,j)   + dx*xb%slp(i+1,j)) &
89                     + dy *(dxm*xb%slp(i,j+1) + dx*xb%slp(i+1,j+1))
90 
91       !------------------------------------------------------------------------
92       ! [2.0] Initialise components of innovation vector:
93       !------------------------------------------------------------------------
94 
95       iv % bogus(n) % slp % inv = 0.0
96 
97       if (ABS(ob % bogus(n) % slp - missing_r) > 1. .AND. &
98            iv % bogus(n) % slp % qc >= obs_qc_pointer) then
99         iv % bogus(n) % slp % inv = ob % bogus(n) % slp - model_slp
100       end if
101 
102       do k = 1, iv % bogus(n) % info % levels
103         iv % bogus(n) % u(k) % inv = 0.0
104         iv % bogus(n) % v(k) % inv = 0.0
105         iv % bogus(n) % t(k) % inv = 0.0
106         iv % bogus(n) % q(k) % inv = 0.0
107 
108         !------------------------------------------------------------------------
109         ! [4.0] Fast interpolation:
110         !------------------------------------------------------------------------
111 
112         if (ob % bogus(n) % u(k) > missing_r .AND. &
113              iv % bogus(n) % u(k) % qc >= obs_qc_pointer) then
114           iv % bogus(n) % u(k) % inv = ob % bogus(n) % u(k) - model_u(k)
115         end if
116 
117         if (ob % bogus(n) % v(k) > missing_r .AND. &
118              iv % bogus(n) % v(k) % qc >= obs_qc_pointer) then
119           iv % bogus(n) % v(k) % inv = ob % bogus(n) % v(k) - model_v(k)
120         end if
121 
122         if (ob % bogus(n) % t(k) > missing_r .AND. &
123              iv % bogus(n) % t(k) % qc >= obs_qc_pointer) then
124 
125           ! only for global Bogus(YRG 07/15/2005):
126           if (iv % bogus(n) % info % platform(8:12) /= 'TCBOG') then
127             iv % bogus(n) % t(k) % inv = ob % bogus(n) % t(k) - model_t(k)
128           else
129             iv % bogus(n) % t(k) % inv = missing_r 
130             iv % bogus(n) % t(k) % qc  = missing_data
131           end if
132         end if
133 
134         if (ob % bogus(n) % q(k) > missing_r .AND. &
135              iv % bogus(n) % q(k) % qc >= obs_qc_pointer) then
136 
137           ! only for global Bogus(YRG 07/15/2005):
138           if (iv % bogus(n) % info % platform(8:12) /= 'TCBOG') then
139             iv % bogus(n) % q(k) % inv = ob % bogus(n) % q(k) - model_q(k)
140           else
141             iv % bogus(n) % q(k) % inv = missing_r 
142             iv % bogus(n) % q(k) % qc  = missing_data
143           end if
144         end if
145       end do
146 
147       !------------------------------------------------------------------------
148       ! [5.0] Perform optional maximum error check:
149       !------------------------------------------------------------------------
150 
151       if (check_max_iv) then
152         call da_check_max_iv_bogus(it, iv % bogus(n), &
153           itu,ituf,itvv,itvvf,itt,ittf,itqv,itqvf,itslp,itslpf)
154       end if
155    end do
156 
157    if (rootproc .and. check_max_iv_print) then
158       write(unit = check_max_iv_unit, fmt ='(A,i5,A)')'For outer iteration ',it, &
159          ', Total Rejections for Bogus follows:'
160 
161       write(unit = check_max_iv_unit, fmt = '(/,10(2(A,I6),/))') &
162         'Number of failed u-wind       observations:     ',ituf,  ' on ',itu,   &
163         'Number of failed v-wind       observations:     ',itvvf, ' on ',itvv,   &
164         'Number of failed temperature  observations:     ',ittf,  ' on ',itt,   &
165         'Number of failed mixing ratio observations:     ',itqvf, ' on ',itqv,   &
166         'Number of failed slp          observations:     ',itslpf,' on ',itslp,   &
167         'Finally Total Bogus rejections ',ituf+itvvf+ittf+itqvf+itslpf,' on ', &
168                                          itu +itvv +itt +itqv +itslp
169    end if
170 
171 end subroutine da_get_innov_vector_bogus
172 
173