da_get_innov_vector_polaramv.inc

References to this file elsewhere.
1 subroutine da_get_innov_vector_polaramv( it, xb, xp, ob, iv)
2 
3    !----------------------------------------------------------------------
4    ! Purpose: Calculates innovation vector, does QC for  Polar AMV's             
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(in)    :: 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    real, dimension(1:max_ob_levels) :: model_u  ! Model value u at ob location.
25    real, dimension(1:max_ob_levels) :: model_v  ! Model value v at ob location.
26 
27    real, dimension(xp%kms:xp%kme) :: v_p      ! Model value p at ob hor. location.
28 
29    integer           :: itu,ituf,itvv,itvvf
30 
31    if (iv % num_polaramv < 1) return
32    
33    if (trace_use) call da_trace_entry("da_get_innov_vector_polaramv")
34 
35    itu   = 0; itvv    = 0;
36    ituf  = 0; itvvf   = 0;
37 
38    do n=iv%ob_numb(iv%current_ob_time-1)%polaramv + 1, iv%ob_numb(iv%current_ob_time)%polaramv
39       ! [1.3] Get horizontal interpolation weights:
40 
41       num_levs = iv % polaramv(n) % info % levels
42       if (num_levs < 1) cycle
43 
44       model_u(:) = 0.0
45       model_v(:) = 0.0
46 
47       i = iv%polaramv(n)%loc%i
48       j = iv%polaramv(n)%loc%j
49       dx = iv%polaramv(n)%loc%dx
50       dy = iv%polaramv(n)%loc%dy
51       dxm = iv%polaramv(n)%loc%dxm
52       dym = iv%polaramv(n)%loc%dym
53 
54       do k=xp%kts,xp%kte
55          v_p(k) = dym*(dxm*xb%p(i,j  ,k) + dx*xb%p(i+1,j  ,k)) &
56                 + dy *(dxm*xb%p(i,j+1,k) + dx*xb%p(i+1,j+1,k))
57       end do
58 
59       do k=1, iv % polaramv(n) % info % levels
60          iv%polaramv(n)%zk(k)=missing_r
61 
62          if (iv % polaramv(n) % p(k) > 1.0) then
63             call da_to_zk(iv % polaramv(n) % p(k), v_p, xp, v_interp_p, iv%polaramv(n)%zk(k))
64          end if
65 
66          if (iv%polaramv(n)%zk(k) < 0.0 .and.  .not.anal_type_verify) then
67             iv % polaramv(n) % u(k) % qc = missing
68             iv % polaramv(n) % v(k) % qc = missing
69          end if
70       end do
71 
72       call da_interp_lin_3d( xb % u, xp, i, j, dx, dy, dxm, dym, &
73                           model_u, max_ob_levels, iv%polaramv(n)%zk, num_levs)
74       call da_interp_lin_3d( xb % v, xp, i, j, dx, dy, dxm, dym, &
75                           model_v, max_ob_levels, iv%polaramv(n)%zk, num_levs)
76 
77       do k = 1, iv % polaramv(n) % info % levels
78          iv % polaramv(n) % u(k) % inv = 0.0
79          iv % polaramv(n) % v(k) % inv = 0.0
80 
81 
82          if (ob % polaramv(n) % u(k) > missing_r .AND. &
83               iv % polaramv(n) % u(k) % qc >= obs_qc_pointer) then
84 
85               iv % polaramv(n) % u(k) % inv = ob % polaramv(n) % u(k) - &
86                                            model_u(k)
87          end if
88 
89          if (ob % polaramv(n) % v(k) > missing_r .AND. &
90               iv % polaramv(n) % v(k) % qc >= obs_qc_pointer) then
91 
92               iv % polaramv(n) % v(k) % inv = ob % polaramv(n) % v(k) - &
93                                            model_v(k)
94          end if
95       end do
96 
97       !------------------------------------------------------------------------
98       ! Perform optional maximum error check:
99       !------------------------------------------------------------------------
100 
101       if (check_max_iv) call da_check_max_iv_polaramv(it, iv % polaramv(n), &
102                              itu,ituf,itvv,itvvf)
103    end do
104 
105    if (rootproc .and. check_max_iv_print) then
106       write(unit = check_max_iv_unit, fmt ='(A,i5,A)')'For outer iteration ',it, &
107          ', Total Rejections for Polar AMVs follows:'
108       write(unit = check_max_iv_unit, fmt = '(/,10(2(A,I6),/))') &
109          'Number of failed u-wind observations:     ',ituf, ' on ',itu,   &
110          'Number of failed v-wind observations:     ',itvvf,' on ',itvv,   &
111          'Finally Total Polar AMVs rejections ',ituf+itvvf,' on ',itu+itvv
112    end if
113    
114    if (trace_use) call da_trace_exit("da_get_innov_vector_polaramv")
115 
116 end subroutine da_get_innov_vector_polaramv
117 
118