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