da_get_innov_vector_airsr.inc

References to this file elsewhere.
1 subroutine da_get_innov_vector_airsr( it, xb, xp, ob, iv)
2 
3    !----------------------------------------------------------------------
4    ! Purpose:   a) Rcomputes innovation vecrot at
5    !                   AIRS retrieval locations 
6    !                b) Does quality control check on innovation vector
7    !                   if required.
8    ! Update :
9    !     01/24/2007    Syed RH Rizvi
10    !     Updated for "VERIFY"       
11    !----------------------------------------------------------------------
12 
13    implicit none
14 
15    integer,           intent(in)    :: it       ! External iteration.
16    type(xb_type),    intent(in)    :: xb       ! first guess state.
17    type(xpose_type), intent(in)    :: xp       ! Domain decomposition vars.
18    type(y_type),     intent(inout) :: ob       ! Observation structure.
19    type(ob_type),    intent(inout) :: iv       ! O-B structure.
20 
21    integer                        :: n        ! Loop counter.
22    integer                        :: i, j, k  ! Index dimension.
23    integer                        :: num_levs ! Number of obs levels.
24    real                           :: dx, dxm  ! Interpolation weights.
25    real                           :: dy, dym  ! Interpolation weights.
26 
27    real, dimension(1:max_ob_levels) :: model_t  ! Model value t at ob location.
28    real, dimension(1:max_ob_levels) :: model_q  ! Model value q at ob location.
29    ! real, dimension(1:max_ob_levels) :: model_h  ! Model value h at ob location.
30 
31    real, dimension(xp%kms:xp%kme)   :: v_h      ! Model value h at ob hor. location.
32    real, dimension(xp%kms:xp%kme)   :: v_p      ! Model value p at ob hor. location.
33 
34    integer                      :: itt,ittf,itqv,itqvf
35 
36    if (iv%num_airsr < 1) return
37 
38    itt  = 0; itqv  = 0;
39    ittf = 0; itqvf = 0;
40 
41    do n=iv%ob_numb(iv%current_ob_time-1)%airsr + 1, iv%ob_numb(iv%current_ob_time)%airsr
42       num_levs = iv%airsr(n)%info%levels
43 
44       if (num_levs < 1) cycle
45 
46       model_t(:) = 0.0
47       model_q(:) = 0.0
48       ! model_h(:) = 0.0
49 
50       ! [1.1] Get horizontal interpolation weights:
51 
52       i = iv%airsr(n)%loc%i
53       j = iv%airsr(n)%loc%j
54       dx = iv%airsr(n)%loc%dx
55       dy = iv%airsr(n)%loc%dy
56       dxm = iv%airsr(n)%loc%dxm
57       dym = iv%airsr(n)%loc%dym
58 
59       do k=xp%kts,xp%kte
60          v_h(k) = dym*(dxm*xb%h(i,j  ,k) + dx*xb%h(i+1,j  ,k)) &
61                 + dy *(dxm*xb%h(i,j+1,k) + dx*xb%h(i+1,j+1,k))
62          v_p(k) = dym*(dxm*xb%p(i,j  ,k) + dx*xb%p(i+1,j  ,k)) &
63                 + dy *(dxm*xb%p(i,j+1,k) + dx*xb%p(i+1,j+1,k))
64       end do
65 
66       do k=1, num_levs
67          iv%airsr(n)%zk(k)=missing_r
68 
69          if (iv%airsr(n)%p(k) > 1.0) then
70             call da_to_zk(iv%airsr(n)%p(k), v_p, xp, v_interp_p, iv%airsr(n)%zk(k))
71          else if (iv%airsr(n)%h(k) > 0.0) then
72             call da_to_zk(iv%airsr(n)%h(k), v_h, xp, v_interp_h, iv%airsr(n)%zk(k))
73          end if
74 
75          if (iv%airsr(n)%zk(k) < 0.0 .and.  .not.anal_type_verify) then
76             iv%airsr(n)%t(k)%qc = missing
77             iv%airsr(n)%q(k)%qc = missing
78          end if
79       end do
80 
81       ! [1.2] Interpolate horizontally to ob:
82 
83       call da_interp_lin_3d( xb%t, xp, i, j, dx, dy, dxm, dym, &
84                            model_t, max_ob_levels, iv%airsr(n)%zk, num_levs)
85       call da_interp_lin_3d( xb%q, xp, i, j, dx, dy, dxm, dym, &
86                            model_q, max_ob_levels, iv%airsr(n)%zk, num_levs)
87 
88       ! call da_interp_lin_3d( xb%h, xp, i, j, dx, dy, dxm, dym, &
89       !   model_h, max_ob_levels, iv%airsr(n)%zk, num_levs)
90 
91       !------------------------------------------------------------------------
92       ! [2.0] Initialise components of innovation vector:
93       !------------------------------------------------------------------------
94 
95       do k = 1, iv%airsr(n)%info%levels
96          iv%airsr(n)%t(k)%inv = 0.0
97          iv%airsr(n)%q(k)%inv = 0.0
98 
99          !------------------------------------------------------------------------
100          ! [3.0] Interpolation:
101          !------------------------------------------------------------------------
102 
103          if (ob%airsr(n)%t(k) > missing_r .AND. &
104               iv%airsr(n)%t(k)%qc >= obs_qc_pointer) then
105            iv%airsr(n)%t(k)%inv = ob%airsr(n)%t(k) - model_t(k)
106          end if
107 
108          if (ob%airsr(n)%q(k) > missing_r .AND. &
109              iv%airsr(n)%q(k)%qc >= obs_qc_pointer) then
110             iv%airsr(n)%q(k)%inv = ob%airsr(n)%q(k) - model_q(k)
111          end if
112       end do
113 
114       !------------------------------------------------------------------------
115       ! [5.0] Perform optional maximum error check:
116       !------------------------------------------------------------------------
117 
118       if (check_max_iv) &    
119          call da_check_max_iv_airsr(it, iv % airsr(n),itt,ittf,itqv,itqvf)
120    end do
121 
122    if (rootproc .and. check_max_iv_print) then
123       write(unit = check_max_iv_unit, fmt ='(A,i5,A)')'For outer iteration ',it, &
124          ', Total Rejections for airsr follows:'
125 
126       write(unit = check_max_iv_unit, fmt = '(/,10(2(A,I6),/))') &
127          'Number of failed temperature observations:',ittf, ' on ',itt,   &
128          'Number of failed mixing ratio observations:',itqvf,' on ',itqv, &
129          'Finally Total airsr rejections: ',ittf+itqvf,' on ',itt +itqv
130    end if
131 
132 end subroutine da_get_innov_vector_airsr
133 
134