da_get_innov_vector_pilot.inc

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