da_get_innov_vector_ssmi_tb.inc

References to this file elsewhere.
1 subroutine da_get_innov_vector_ssmi_tb( it, xb, ob, iv)
2 
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
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(y_type),  intent(in)    :: ob         ! Observation structure.
12    type(ob_type), intent(inout) :: iv         ! O-B structure.
13 
14    integer                   :: n           ! Loop counter.
15    integer                   :: i, j        ! Index dimension.
16    real                      :: dx, dxm     ! Interpolation weights.
17    real                      :: dy, dym     ! Interpolation weights.
18    real                      :: model_tb19h ! Model value tb19h at oblocation.
19    real                      :: model_tb19v ! Model value tb19v at oblocation.
20    real                      :: model_tb22v ! Model value tb22v at oblocation.
21    real                      :: model_tb37h ! Model value tb37h at oblocation.
22    real                      :: model_tb37v ! Model value tb37v at oblocation.
23    real                      :: model_tb85h ! Model value tb85h at oblocation.
24    real                      :: model_tb85v ! Model value tb85v at oblocation.
25 
26    integer                   :: itb19v, itb19vf, &
27                                    itb19h, itb19hf, &
28                                    itb22v, itb22vf, &
29                                    itb37v, itb37vf, &
30                                    itb37h, itb37hf, &
31                                    itb85v, itb85vf, &
32                                    itb85h, itb85hf
33 
34    if (iv % num_ssmi_tb > 0) then
35 
36       itb19v = 0
37       itb19vf = 0
38       itb19h = 0
39       itb19hf = 0
40       itb22v = 0
41       itb22vf = 0
42       itb37v = 0
43       itb37vf = 0
44       itb37h = 0
45       itb37hf = 0
46       itb85v = 0
47       itb85vf = 0
48       itb85h = 0
49       itb85hf = 0
50 
51       do n=1, iv % num_ssmi_tb
52 
53          ! COMPUTE inNOVATION VECTOR
54          ! =========================
55 
56          !  Obs coordinates on model grid
57 
58          ! TB
59 
60          i = iv%ssmi_tb(n)%loc%i
61          j = iv%ssmi_tb(n)%loc%j
62          dx = iv%ssmi_tb(n)%loc%dx
63          dy = iv%ssmi_tb(n)%loc%dy
64          dxm = iv%ssmi_tb(n)%loc%dxm
65          dym = iv%ssmi_tb(n)%loc%dym
66 
67          ! Tb19h
68 
69          if (abs(ob % ssmi_tb(n) % tb19h - missing_r) > 1.) then
70             model_tb19h = dym*(dxm*xb%tb19h(i,j)   + dx*xb%tb19h(i+1,j)) &
71                          + dy *(dxm*xb%tb19h(i,j+1) + dx*xb%tb19h(i+1,j+1))
72             iv % ssmi_tb(n) % tb19h % inv = ob % ssmi_tb(n) % tb19h - &
73                model_tb19h
74          else
75             iv % ssmi_tb(n) % tb19h % inv = 0.
76          end if
77 
78          ! Tb19v
79 
80          if (abs(ob % ssmi_tb(n) % tb19v - missing_r) > 1.) then
81             model_tb19v = dym*(dxm*xb%tb19v(i,j)   + dx *xb%tb19v(i+1,j)) &
82                          + dy *(dxm*xb%tb19v(i,j+1) + dx *xb%tb19v(i+1,j+1))
83             iv % ssmi_tb(n) % tb19v % inv = ob % ssmi_tb(n) % tb19v - &
84                model_tb19v
85          else
86             iv % ssmi_tb(n) % tb19v % inv = 0.
87          end if
88 
89         ! Tb19v
90 
91          if (abs(ob % ssmi_tb(n) % tb22v - missing_r) > 1.) then
92             model_tb22v = dym*(dxm*xb%tb22v(i,j) + dx *xb%tb22v(i+1,j)) &
93                          + dy *(dxm*xb%tb22v(i,j+1) + dx *xb%tb22v(i+1,j+1))
94             iv % ssmi_tb(n) % tb22v % inv = ob % ssmi_tb(n) % tb22v - &
95                model_tb22v
96          else
97             iv % ssmi_tb(n) % tb22v % inv = 0.0
98          end if
99 
100          ! Tb37h
101 
102          if (abs(ob % ssmi_tb(n) % tb37h - missing_r) > 1.) then
103             model_tb37h = dym*(dxm*xb%tb37h(i,j)  + dx *xb%tb37h(i+1,j)) &
104                          + dy *(dxm*xb%tb37h(i,j+1) + dx *xb%tb37h(i+1,j+1))
105             iv % ssmi_tb(n) % tb37h % inv = ob % ssmi_tb(n) % tb37h - &
106                model_tb37h
107          else
108             iv % ssmi_tb(n) % tb37h % inv = 0.
109          end if
110 
111          ! Tb37v
112 
113          if (abs(ob % ssmi_tb(n) % tb37v - missing_r) > 1.) then
114             model_tb37v = dym*(dxm*xb%tb37v(i,j)  + dx *xb%tb37v(i+1,j)) &
115                          + dy *(dxm*xb%tb37v(i,j+1) + dx *xb%tb37v(i+1,j+1))
116              iv % ssmi_tb(n) % tb37v % inv = ob % ssmi_tb(n) % tb37v - &
117                 model_tb37v
118          else
119              iv % ssmi_tb(n) % tb37v % inv = 0.
120          end if
121 
122          ! Tb85h
123 
124          if (abs(ob % ssmi_tb(n) % tb85h - missing_r) > 1.) then
125              model_tb85h = dym*(dxm*xb%tb85h(i,j) + dx *xb%tb85h(i+1,j)) &
126                          + dy *(dxm*xb%tb85h(i,j+1) + dx *xb%tb85h(i+1,j+1))
127             iv % ssmi_tb(n) % tb85h % inv = ob % ssmi_tb(n) % tb85h - &
128                model_tb85h
129          else
130             iv % ssmi_tb(n) % tb85h % inv = 0.
131          end if
132 
133          ! Tb85v
134 
135          if (abs(ob % ssmi_tb(n) % tb85v - missing_r) > 1.) then
136              model_tb85v = dym*(dxm*xb%tb85v(i,j) + dx *xb%tb85v(i+1,j)) &
137                          + dy *(dxm*xb%tb85v(i,j+1) + dx *xb%tb85v(i+1,j+1))
138             iv % ssmi_tb(n) % tb85v % inv = ob % ssmi_tb(n) % tb85v -  &
139                model_tb85v
140          else
141             iv % ssmi_tb(n) % tb85v % inv = 0.
142          end if
143 
144 
145          !----------------------------------------------------------------
146          !     Perform optional maximum error check:
147          !----------------------------------------------------------------
148 
149          if (check_max_iv) then
150             call da_check_max_iv_ssmi_tb(it, iv % ssmi_tb(n), &
151                                    itb19v, itb19vf, itb19h, itb19hf, &
152                                    itb22v, itb22vf,                  &
153                                    itb37v, itb37vf, itb37h, itb37hf, &
154                                    itb85v, itb85vf, itb85h, itb85hf )
155          end if
156       end do
157 
158       if (rootproc .and. check_max_iv_print) then
159          write(unit = check_max_iv_unit, fmt ='(A,i5,A)')'For outer iteration ', &
160             it,'Total Rejections for SSMI(Tbs) follows:'
161          write(unit = check_max_iv_unit, fmt = '(2(A,I6))') &
162             'Number of failed ssmi tb19v observations: ',itb19v, ' on ',itb19v, &
163             'Number of failed ssmi tb19h observations: ',itb19h, ' on ',itb19h, &
164             'Number of failed ssmi tb22v observations: ',itb22v, ' on ',itb22v, &
165             'Number of failed ssmi tb37v observations: ',itb37v, ' on ',itb37v, &
166             'Number of failed ssmi tb37h observations: ',itb37h, ' on ',itb37h, &
167             'Number of failed ssmi tb85v observations: ',itb85v, ' on ',itb85v, &
168             'Number of failed ssmi tb85h observations: ',itb85h, ' on ',itb85h
169       end if
170    end if
171 
172 end subroutine da_get_innov_vector_ssmi_tb
173 
174