da_get_innov_vector_pilot.inc

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