da_get_innov_vector_ssmt2.inc

References to this file elsewhere.
1 subroutine da_get_innov_vector_ssmt2( it, xb, xp, 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(xpose_type), intent(in)  :: xp       ! Domain decomposition vars.
12    type(y_type),  intent(inout)  :: 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    real                           :: dx, dxm  ! Interpolation weights.
19    real                           :: dy, dym  ! Interpolation weights.
20    real, dimension(1:max_ob_levels):: model_rh ! Model value rh at ob location.
21 
22    real                            :: zk       ! Interpolation vertical coordinator.
23 
24    real, dimension(xp%kms:xp%kme)  :: v_h      ! Model value h at ob hor. location.
25    real, dimension(xp%kms:xp%kme)  :: v_p      ! Model value p at ob hor. location.
26 
27    integer                         :: irh, irhf
28    
29    if (iv % num_ssmt2 < 1) return
30    
31    if (trace_use) call da_trace_entry("da_get_innov_vector_ssmt2")
32    
33    irh = 0 ; irhf = 0
34    do n=1, iv % num_ssmt2
35       num_levs = iv % ssmt2(n) % info % levels
36 
37       if (num_levs < 1) cycle
38 
39       model_rh(:) = 0.0
40 
41       ! [1.1] Get horizontal interpolation weights:
42 
43       i = iv%ssmt2(n)%loc%i
44       j = iv%ssmt2(n)%loc%j
45       dx = iv%ssmt2(n)%loc%dx
46       dy = iv%ssmt2(n)%loc%dy
47       dxm = iv%ssmt2(n)%loc%dxm
48       dym = iv%ssmt2(n)%loc%dym
49 
50       do k=xp%kts,xp%kte
51          v_h(k) = dym*(dxm*xb%h(i,j  ,k) + dx*xb%h(i+1,j  ,k)) &
52                 + dy *(dxm*xb%h(i,j+1,k) + dx*xb%h(i+1,j+1,k))
53          v_p(k) = dym*(dxm*xb%p(i,j  ,k) + dx*xb%p(i+1,j  ,k)) &
54                 + dy *(dxm*xb%p(i,j+1,k) + dx*xb%p(i+1,j+1,k))
55       end do
56 
57       num_levs=0
58       do k=1, iv % ssmt2(n) % info % levels
59 
60          zk=missing_r
61 
62          if (iv % ssmt2(n) % h(k) > 0.0) then
63             call da_to_zk(iv % ssmt2(n) % h(k), v_h, xp, v_interp_h, zk)
64          else if (iv % ssmt2(n) % p(k) > 1.0) then
65             call da_to_zk(iv % ssmt2(n) % p(k), v_p, xp, v_interp_p, zk)
66          end if
67 
68          if (zk > 0.0) then
69             num_levs=num_levs+1
70             iv%ssmt2(n)%zk(num_levs)=zk
71 
72             ob % ssmt2(n) % rh(num_levs) = ob % ssmt2(n) % rh(k)
73 
74             iv % ssmt2(n) % rh(num_levs) % qc = iv % ssmt2(n) % rh(k) % qc
75 
76             iv % ssmt2(n) % rh(num_levs) % error = iv % ssmt2(n) % rh(k) % error
77          end if
78       end do
79 
80       iv % ssmt2(n) % info % levels = num_levs
81 
82       ! [1.2] Interpolate horizontally to ob:
83 
84       call da_interp_lin_3d( xb % rh, xp, i, j, dx, dy, dxm, dym, &
85                           model_rh, max_ob_levels, iv%ssmt2(n)%zk, num_levs)
86 
87       !--------------------------------------------------------------------
88       ! [2.0] Initialise components of innovation vector:
89       !--------------------------------------------------------------------
90 
91       do k = 1, iv % ssmt2(n) % info % levels
92 
93          iv % ssmt2(n) % rh(k) % inv = 0.0
94 
95          !-----------------------------------------------------------------
96          ! [3.0] Interpolation:
97          !-----------------------------------------------------------------
98 
99          if (ob % ssmt2(n) % rh(k) > missing_r .AND. &
100             iv % ssmt2(n) % rh(k) % qc >= obs_qc_pointer) then
101 
102             iv % ssmt2(n) % rh(k) % inv = ob % ssmt2(n) % rh(k) - &
103                                          model_rh(k)
104          end if
105 
106          ! write(122,'(2i4,i8,5f15.5)')n, k, iv%ssmt2(n)%height_qc(k), &
107          ! iv%ssmt2(n)%info%lat, iv%ssmt2(n)%info%lon, &
108          ! iv%ssmt2(n)%h(k), &
109          ! ob % ssmt2(n) % rh(k), model_rh(k)
110 
111       end do
112 
113       !------------------------------------------------------------------------
114       ! [5.0] Perform optional maximum error check:
115       !------------------------------------------------------------------------
116 
117       if (check_max_iv) call da_check_max_iv_ssmt2(it, iv % ssmt2(n), irh,irhf)
118    end do
119 
120    if (rootproc .and. check_max_iv_print) then
121       write(unit = check_max_iv_unit, fmt ='(A,i5,A)')    &
122          'For outer iteration ',it,' Total Rejections for SSMI(T2) follows:'
123       write(unit = check_max_iv_unit, fmt = '(/,2(2(A,I6),/))') &
124          'Number of failed RH(T2) observations:    ',irhf, ' on ',irh
125    end if
126    
127    if (trace_use) call da_trace_exit("da_get_innov_vector_ssmt2")
128 
129 end subroutine da_get_innov_vector_ssmt2
130 
131