updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_ssmi / da_get_innov_vector_ssmt2.inc
blobcf42e29bbede23a333a2d4c9e288aa998fb64ab6
1 subroutine da_get_innov_vector_ssmt2 (it, num_qcstat_conv,grid, ob, iv)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
7    implicit none
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    integer,          intent(inout) :: num_qcstat_conv(:,:,:,:)
15    integer :: n        ! Loop counter.
16    integer :: i, j, k  ! Index dimension.
17    integer :: num_levs ! Number of obs levels.
18    real    :: dx, dxm  ! Interpolation weights.
19    real    :: dy, dym  ! Interpolation weights.
20    real, allocatable :: model_rh(:,:) ! Model value rh at ob location.
22    real    :: v_h(kms:kme)      ! Model value h at ob hor. location.
23    real    :: v_p(kms:kme)      ! Model value p at ob hor. location.
25    if (trace_use_dull) call da_trace_entry("da_get_innov_vector_ssmt2")
27    allocate (model_rh(1:iv%info(ssmt2)%max_lev,iv%info(ssmt2)%n1:iv%info(ssmt2)%n2))
28    model_rh(:,:) = 0.0
30    if ( it > 1 )then
31       do n=iv%info(ssmt2)%n1,iv%info(ssmt2)%n2
32          do k = 1, iv%info(ssmt1)%levels(n)
33             if(iv % ssmt2(n) % rh(k) % qc == fails_error_max)iv % ssmt2(n) % rh(k) % qc = 0
34          end do
35       end do
36    end if
38    do n=iv%info(ssmt2)%n1,iv%info(ssmt2)%n2
39       num_levs = iv%info(ssmt2)%levels(n)
41       if (num_levs < 1) cycle
43       ! [1.1] Get horizontal interpolation weights:
45       i   = iv%info(ssmt2)%i(1,n)
46       j   = iv%info(ssmt2)%j(1,n)
47       dx  = iv%info(ssmt2)%dx(1,n)
48       dy  = iv%info(ssmt2)%dy(1,n)
49       dxm = iv%info(ssmt2)%dxm(1,n)
50       dym = iv%info(ssmt2)%dym(1,n)
52       do k=kts,kte
53          v_h(k) = dym*(dxm*grid%xb%h(i,j  ,k) + dx*grid%xb%h(i+1,j  ,k)) &
54                 + dy *(dxm*grid%xb%h(i,j+1,k) + dx*grid%xb%h(i+1,j+1,k))
55          v_p(k) = dym*(dxm*grid%xb%p(i,j  ,k) + dx*grid%xb%p(i+1,j  ,k)) &
56                 + dy *(dxm*grid%xb%p(i,j+1,k) + dx*grid%xb%p(i+1,j+1,k))
57       end do
59       num_levs=0
60       do k=1, iv%info(ssmt2)%levels(n)
61          if (iv % ssmt2(n) % h(k) > 0.0) then
62             call da_to_zk(iv % ssmt2(n) % h(k), v_h, v_interp_h, iv%info(ssmt2)%zk(k,n))
63          else if (iv % ssmt2(n) % p(k) > 1.0) then
64             call da_to_zk(iv % ssmt2(n) % p(k), v_p, v_interp_p, iv%info(ssmt2)%zk(k,n))
65          end if
67          if (iv%info(ssmt2)%zk(k,n) > 0.0) then
68             num_levs=num_levs+1
69             iv%info(ssmt2)%zk(num_levs,n)=iv%info(ssmt2)%zk(k,n)
71             ob % ssmt2(n) % rh(num_levs) = ob % ssmt2(n) % rh(k)
73             iv % ssmt2(n) % rh(num_levs) % qc = iv % ssmt2(n) % rh(k) % qc
75             iv % ssmt2(n) % rh(num_levs) % error = iv % ssmt2(n) % rh(k) % error
76          end if
77       end do
79       iv%info(ssmt2)%levels(n) = num_levs
80    end do
82    call da_convert_zk (iv%info(ssmt2))
84    call da_interp_lin_3d (grid%xb%rh, iv%info(ssmt2), model_rh)
86    do n=iv%info(ssmt2)%n1, iv%info(ssmt2)%n2
87       do k = 1, iv%info(ssmt2)%levels(n)
88          iv % ssmt2(n) % rh(k) % inv = 0.0
90          !-----------------------------------------------------------------
91          ! [3.0] Interpolation:
92          !-----------------------------------------------------------------
94          if (ob % ssmt2(n) % rh(k) > missing_r .AND. &
95              iv % ssmt2(n) % rh(k) % qc >= obs_qc_pointer) then
96             iv % ssmt2(n) % rh(k) % inv = ob % ssmt2(n) % rh(k) - model_rh(k,n)
97          end if
99          ! write(122,'(2i4,i8,5f15.5)')n, k, iv%ssmt2(n)%height_qc(k), &
100          ! iv%ssmt2(n)%info%lat, iv%ssmt2(n)%info%lon, &
101          ! iv%ssmt2(n)%h(k), &
102          ! ob % ssmt2(n) % rh(k), model_rh(k,n)
104       end do
105    end do
107    !------------------------------------------------------------------------
108    ! [5.0] Perform optional maximum error check:
109    !------------------------------------------------------------------------
111    if (check_max_iv) call da_check_max_iv_ssmt2(iv, it, num_qcstat_conv)           
113    deallocate (model_rh)
114    
115    if (trace_use_dull) call da_trace_exit("da_get_innov_vector_ssmt2")
117 end subroutine da_get_innov_vector_ssmt2