da_get_innov_vector_ssmt1.inc

References to this file elsewhere.
1 subroutine da_get_innov_vector_ssmt1( 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    integer :: num_levs ! Number of obs levels.
17    real    :: dx, dxm  ! Interpolation weights.
18    real    :: dy, dym  ! Interpolation weights.
19    real, allocatable :: model_t(:,:)  ! Model value t at ob location.
20 
21    real   :: v_h(kms:kme)      ! Model value h at ob hor. location.
22    real   :: v_p(kms:kme)      ! Model value p at ob hor. location.
23 
24    integer :: itt, ittf
25     
26    if (trace_use_dull) call da_trace_entry("da_get_innov_vector_ssmt1")
27 
28    itt = 0 ; ittf = 0
29 
30    allocate (model_t(1:iv%info(ssmt1)%max_lev,iv%info(ssmt1)%n1:iv%info(ssmt1)%n2))
31    model_t(:,:) = 0.0
32 
33    do n=iv%info(ssmt1)%n1,iv%info(ssmt1)%n2
34 
35       num_levs = iv%info(ssmt1)%levels(n)
36 
37       if (num_levs < 1) cycle
38 
39       ! [1.1] Get horizontal interpolation weights:
40 
41       i   = iv%info(ssmt1)%i(1,n)
42       j   = iv%info(ssmt1)%j(1,n)
43       dx  = iv%info(ssmt1)%dx(1,n)
44       dy  = iv%info(ssmt1)%dy(1,n)
45       dxm = iv%info(ssmt1)%dxm(1,n)
46       dym = iv%info(ssmt1)%dym(1,n)
47 
48       do k=kts,kte
49          v_h(k) = dym*(dxm*grid%xb%h(i,j,k)+dx*grid%xb%h(i+1,j,k)) + dy*(dxm*grid%xb%h(i,j+1,k)+dx*grid%xb%h(i+1,j+1,k))
50          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))
51       end do
52 
53       num_levs=0
54 
55       do k=1, iv%info(ssmt1)%levels(n)
56          if (iv % ssmt1(n) % h(k) > 0.0) then
57             call da_to_zk(iv % ssmt1(n) % h(k), v_h, v_interp_h, iv%info(ssmt1)%zk(k,n))
58          else if (iv % ssmt1(n) % p(k) > 1.0) then
59             call da_to_zk(iv % ssmt1(n) % p(k), v_p, v_interp_p, iv%info(ssmt1)%zk(k,n))
60          end if
61 
62          if ( iv%info(ssmt1)%zk(k,n) > 0.0) then
63             num_levs=num_levs+1
64             iv%info(ssmt1)%zk(num_levs,n)=iv%info(ssmt1)%zk(k,n)
65 
66             ob % ssmt1(n) % t(num_levs)         = ob % ssmt1(n) % t(k)
67             iv % ssmt1(n) % t(num_levs) % qc    = iv % ssmt1(n) % t(k) % qc
68             iv % ssmt1(n) % t(num_levs) % error = iv % ssmt1(n) % t(k) % error
69          end if
70       end do
71 
72       iv%info(ssmt1)%levels(n) = num_levs
73    end do
74 
75    call da_convert_zk (iv%info(ssmt1))
76 
77    ! [1.2] Interpolate horizontally to ob:
78 
79    call da_interp_lin_3d (grid%xb%t, iv%info(ssmt1), model_t)
80 
81    do n=iv%info(ssmt1)%n1,iv%info(ssmt1)%n2
82 
83       !---------------------------------------------------------------------
84       ! [2.0] Initialise components of innovation vector:
85       !---------------------------------------------------------------------
86 
87       do k = 1, iv%info(ssmt1)%levels(n)
88          iv % ssmt1(n) % t(k) % inv = 0.0
89 
90          !----------------------------------------------------------------
91          ! [3.0] Interpolation:
92          !----------------------------------------------------------------
93 
94          if (ob % ssmt1(n) % t(k) > missing_r .AND. &
95              iv % ssmt1(n) % t(k) % qc >= obs_qc_pointer) then
96             iv % ssmt1(n) % t(k) % inv = ob % ssmt1(n) % t(k) - model_t(k,n)
97          end if
98       end do
99    end do
100 
101    !------------------------------------------------------------------
102    ! [5.0] Perform optional maximum error check:
103    !------------------------------------------------------------------
104 
105    if (check_max_iv) call da_check_max_iv_ssmt1(iv, it, itt,ittf)
106 
107    if (rootproc .and. check_max_iv_print) then
108       write(unit = check_max_iv_unit, fmt ='(A,i5,A)')    &
109          'For outer iteration ',it,' Total Rejections for SSMI(T1) follows:'
110       write(unit = check_max_iv_unit, fmt = '(/,2(2(A,I6),/))') &
111          'Number of failed RH(T1) observations:    ',ittf, ' on ',itt
112    end if
113 
114    deallocate (model_t)
115     
116    if (trace_use_dull) call da_trace_exit("da_get_innov_vector_ssmt1")
117 
118 end subroutine da_get_innov_vector_ssmt1
119 
120