da_get_innov_vector_bogus.inc

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