da_get_innov_vector_profiler.inc

References to this file elsewhere.
1 subroutine da_get_innov_vector_profiler( 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    integer :: itu,ituf,itvv,itvvf
25    
26    if (trace_use_dull) call da_trace_entry("da_get_innov_vector_profiler")
27    
28    itu   = 0; itvv    = 0;
29    ituf  = 0; itvvf   = 0;
30 
31    allocate (model_u(iv%info(profiler)%max_lev,iv%info(profiler)%n1:iv%info(profiler)%n2))
32    allocate (model_v(iv%info(profiler)%max_lev,iv%info(profiler)%n1:iv%info(profiler)%n2))
33 
34    model_u(:,:) = 0.0
35    model_v(:,:) = 0.0
36 
37    do n=iv%info(profiler)%n1,iv%info(profiler)%n2
38 
39       ! [1.3] Get horizontal interpolation weights:
40 
41       i   = iv%info(profiler)%i(1,n)
42       j   = iv%info(profiler)%j(1,n)
43       dx  = iv%info(profiler)%dx(1,n)
44       dy  = iv%info(profiler)%dy(1,n)
45       dxm = iv%info(profiler)%dxm(1,n)
46       dym = iv%info(profiler)%dym(1,n)
47 
48       do k=kts,kte
49          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))
50       end do
51 
52       do k=1, iv%info(profiler)%levels(n)
53          if (iv % profiler(n) % p(k) > 1.0) then
54             call da_to_zk(iv % profiler(n) % p(k), v_p, v_interp_p, iv%info(profiler)%zk(k,n))
55          end if
56 
57          if (iv%info(profiler)%zk(k,n) < 0.0 .and.  .not.anal_type_verify) then
58             iv % profiler(n) % u(k) % qc = missing
59             iv % profiler(n) % v(k) % qc = missing
60          end if
61       end do
62    end do
63 
64    call da_convert_zk (iv%info(profiler))
65 
66    ! [1.4] Interpolate horizontally:
67    call da_interp_lin_3d (grid%xb%u, iv%info(profiler), model_u)
68    call da_interp_lin_3d (grid%xb%v, iv%info(profiler), model_v)
69 
70 
71    do n=iv%info(profiler)%n1,iv%info(profiler)%n2
72       !------------------------------------------------------------------------
73       ! [2.0] Initialise components of innovation vector:
74       !------------------------------------------------------------------------
75 
76       do k = 1, iv%info(profiler)%levels(n)
77          iv % profiler(n) % u(k) % inv = 0.0
78          iv % profiler(n) % v(k) % inv = 0.0
79 
80          !----------------------------------------------------------------
81          ! [4.0] Fast interpolation:
82          !----------------------------------------------------------------
83 
84          if (ob % profiler(n) % u(k) > missing_r .AND. &
85              iv % profiler(n) % u(k) % qc >= obs_qc_pointer) then
86             iv % profiler(n) % u(k) % inv = ob % profiler(n) % u(k) - model_u(k,n)
87          end if
88 
89          if (ob % profiler(n) % v(k) > missing_r .AND. &
90              iv % profiler(n) % v(k) % qc >= obs_qc_pointer) then
91             iv % profiler(n) % v(k) % inv = ob % profiler(n) % v(k) - model_v(k,n)
92          end if
93       end do
94 
95       !------------------------------------------------------------------
96       ! [5.0] Perform optional maximum error check:
97       !------------------------------------------------------------------
98 
99    end do
100 
101    if (check_max_iv) then  
102       call da_check_max_iv_profiler(iv, it, itu,ituf,itvv,itvvf)
103    end if
104 
105    if (rootproc .and. check_max_iv_print) then
106       write(unit = check_max_iv_unit, fmt ='(A,i5,A)')&
107          'For outer iteration ',it, ', Total Rejections for Profiler 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 Profiler rejections ',ituf+itvvf,' on ',itu +itvv
112    end if
113 
114    deallocate (model_u)
115    deallocate (model_v)
116    
117    if (trace_use_dull) call da_trace_exit("da_get_innov_vector_profiler")
118 
119 end subroutine da_get_innov_vector_profiler
120 
121