wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_ssmi / da_get_innov_vector_ssmi_tb.inc
blob55c6fb23853ae3702dd594f1105405b8a1b8fe1e
1 subroutine da_get_innov_vector_ssmi_tb (it, grid, ob, iv)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
7    implicit none
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.
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.
26    if (trace_use) call da_trace_entry("da_get_innov_vector_ssmi_tb")
28    if ( it > 1 ) then
29       do n=iv%info(ssmi_tb)%n1,iv%info(ssmi_tb)%n2
30          if(iv % ssmi_tb(n) % tb19v % qc == fails_error_max ) iv % ssmi_tb(n) % tb19v % qc = 0
31          if(iv % ssmi_tb(n) % tb19h % qc == fails_error_max ) iv % ssmi_tb(n) % tb19h % qc = 0
32          if(iv % ssmi_tb(n) % tb22v % qc == fails_error_max ) iv % ssmi_tb(n) % tb22v % qc = 0
33          if(iv % ssmi_tb(n) % tb37v % qc == fails_error_max ) iv % ssmi_tb(n) % tb37v % qc = 0
34          if(iv % ssmi_tb(n) % tb37h % qc == fails_error_max ) iv % ssmi_tb(n) % tb37h % qc = 0
35          if(iv % ssmi_tb(n) % tb85v % qc == fails_error_max ) iv % ssmi_tb(n) % tb85v % qc = 0
36          if(iv % ssmi_tb(n) % tb85h % qc == fails_error_max ) iv % ssmi_tb(n) % tb85h % qc = 0
37       end do
38    end if
40    do n=iv%info(ssmi_tb)%n1,iv%info(ssmi_tb)%n2
41       ! compute innovation vector
42       ! =========================
44       !  Obs coordinates on model grid
46       ! TB
48       i   = iv%info(ssmi_tb)%i(1,n)
49       j   = iv%info(ssmi_tb)%j(1,n)
50       dx  = iv%info(ssmi_tb)%dx(1,n)
51       dy  = iv%info(ssmi_tb)%dy(1,n)
52       dxm = iv%info(ssmi_tb)%dxm(1,n)
53       dym = iv%info(ssmi_tb)%dym(1,n)
55       ! Tb19h
57       if (abs(ob % ssmi_tb(n) % tb19h - missing_r) > 1.0) then
58          model_tb19h = dym*(dxm*grid%xb%tb19h(i,j)   + dx*grid%xb%tb19h(i+1,j)) &
59             + dy *(dxm*grid%xb%tb19h(i,j+1) + dx*grid%xb%tb19h(i+1,j+1))
60          iv % ssmi_tb(n) % tb19h % inv = ob % ssmi_tb(n) % tb19h - &
61             model_tb19h
62       else
63          iv % ssmi_tb(n) % tb19h % inv = 0.0
64       end if
66       ! Tb19v
68       if (abs(ob % ssmi_tb(n) % tb19v - missing_r) > 1.0) then
69          model_tb19v = dym*(dxm*grid%xb%tb19v(i,j)   + dx *grid%xb%tb19v(i+1,j)) &
70             + dy *(dxm*grid%xb%tb19v(i,j+1) + dx *grid%xb%tb19v(i+1,j+1))
71          iv % ssmi_tb(n) % tb19v % inv = ob % ssmi_tb(n) % tb19v - &
72             model_tb19v
73       else
74          iv % ssmi_tb(n) % tb19v % inv = 0.0
75       end if
77      ! Tb19v
79       if (abs(ob % ssmi_tb(n) % tb22v - missing_r) > 1.0) then
80          model_tb22v = dym*(dxm*grid%xb%tb22v(i,j) + dx *grid%xb%tb22v(i+1,j)) &
81             + dy *(dxm*grid%xb%tb22v(i,j+1) + dx *grid%xb%tb22v(i+1,j+1))
82          iv % ssmi_tb(n) % tb22v % inv = ob % ssmi_tb(n) % tb22v - &
83             model_tb22v
84       else
85          iv % ssmi_tb(n) % tb22v % inv = 0.0
86       end if
88       ! Tb37h
90       if (abs(ob % ssmi_tb(n) % tb37h - missing_r) > 1.0) then
91          model_tb37h = dym*(dxm*grid%xb%tb37h(i,j)  + dx *grid%xb%tb37h(i+1,j)) &
92             + dy *(dxm*grid%xb%tb37h(i,j+1) + dx *grid%xb%tb37h(i+1,j+1))
93          iv % ssmi_tb(n) % tb37h % inv = ob % ssmi_tb(n) % tb37h - &
94             model_tb37h
95       else
96          iv % ssmi_tb(n) % tb37h % inv = 0.0
97       end if
99       ! Tb37v
101       if (abs(ob % ssmi_tb(n) % tb37v - missing_r) > 1.0) then
102          model_tb37v = dym*(dxm*grid%xb%tb37v(i,j)  + dx *grid%xb%tb37v(i+1,j)) &
103             + dy *(dxm*grid%xb%tb37v(i,j+1) + dx *grid%xb%tb37v(i+1,j+1))
104          iv % ssmi_tb(n) % tb37v % inv = ob % ssmi_tb(n) % tb37v - &
105             model_tb37v
106       else
107          iv % ssmi_tb(n) % tb37v % inv = 0.0
108       end if
110       ! Tb85h
112       if (abs(ob % ssmi_tb(n) % tb85h - missing_r) > 1.0) then
113          model_tb85h = dym*(dxm*grid%xb%tb85h(i,j) + dx *grid%xb%tb85h(i+1,j)) &
114             + dy *(dxm*grid%xb%tb85h(i,j+1) + dx *grid%xb%tb85h(i+1,j+1))
115          iv % ssmi_tb(n) % tb85h % inv = ob % ssmi_tb(n) % tb85h - &
116             model_tb85h
117       else
118          iv % ssmi_tb(n) % tb85h % inv = 0.0
119       end if
121       ! Tb85v
123       if (abs(ob % ssmi_tb(n) % tb85v - missing_r) > 1.0) then
124          model_tb85v = dym*(dxm*grid%xb%tb85v(i,j) + dx *grid%xb%tb85v(i+1,j)) &
125             + dy *(dxm*grid%xb%tb85v(i,j+1) + dx *grid%xb%tb85v(i+1,j+1))
126          iv % ssmi_tb(n) % tb85v % inv = ob % ssmi_tb(n) % tb85v -  &
127             model_tb85v
128       else
129          iv % ssmi_tb(n) % tb85v % inv = 0.0
130       end if
131    end do
133    !----------------------------------------------------------------
134    !     Perform optional maximum error check:
135    !----------------------------------------------------------------
137    if (check_max_iv) call da_check_max_iv_ssmi_tb(iv, it)  
138    
139    if (trace_use) call da_trace_exit("da_get_innov_vector_ssmi_tb")
141 end subroutine da_get_innov_vector_ssmi_tb