wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_airsr / da_get_innov_vector_airsr.inc
blob1077aa7e6375dfc423330c97dd02e7cdbb764595
1 subroutine da_get_innov_vector_airsr( it,num_qcstat_conv, grid, ob, iv)
3    !----------------------------------------------------------------------
4    ! Purpose:   a) Rcomputes innovation vecrot at
5    !                   AIRS retrieval locations 
6    !                b) Does quality control check on innovation vector
7    !                   if required.
8    !    Updated for Analysis on Arakawa-C grid
9    !    Author: Syed RH Rizvi,  MMM/ESSL/NCAR,  Date: 10/22/2008
10    !----------------------------------------------------------------------
12    implicit none
14    integer,          intent(in)    :: it       ! External iteration.
15    type(domain),     intent(in)    :: grid       ! first guess state.
16    type(y_type),     intent(inout) :: ob       ! Observation structure.
17    type(iv_type),    intent(inout) :: iv       ! O-B structure.
18    integer,          intent(inout) :: num_qcstat_conv(:,:,:,:)
21    integer :: n  ! Loop counter.
22    integer :: i, j, k  ! Index dimension.
23    integer :: num_levs ! Number of obs levels.
24    real    :: dx, dxm  ! Interpolation weights.
25    real    :: dy, dym  ! Interpolation weights.
27    ! real    :: model_h(1:max_ob_levels)  ! Model value h at ob location.
29    real    :: v_h(kms:kme)      ! Model value h at ob hor. location.
30    real    :: v_p(kms:kme)      ! Model value p at ob hor. location.
32    real, allocatable :: model_q(:,:)  ! Model value q at ob location.
33    real, allocatable :: model_t(:,:)  ! Model value t at ob location.
35    if (trace_use_dull) call da_trace_entry("da_get_innov_vector_airsr")
37    allocate (model_t(iv%info(airsr)%max_lev,iv%info(airsr)%n1:iv%info(airsr)%n2))
38    allocate (model_q(iv%info(airsr)%max_lev,iv%info(airsr)%n1:iv%info(airsr)%n2))
40    model_t(:,:) = 0.0
41    model_q(:,:) = 0.0
43    if ( it > 1 ) then
44       do n = iv%info(airsr)%n1, iv%info(airsr)%n2
45          do k = 1, iv%info(airsr)%levels(n)
46             if (iv%airsr(n)%t(k)%qc == fails_error_max) iv%airsr(n)%t(k)%qc = 0
47             if (iv%airsr(n)%q(k)%qc == fails_error_max) iv%airsr(n)%q(k)%qc = 0
48          end do
49       end do
50    end if
52    do n=iv%info(airsr)%n1, iv%info(airsr)%n2
53       num_levs = iv%info(airsr)%levels(n)
55       if (num_levs < 1) cycle
57       ! [1.1] Get horizontal interpolation weights:
59       i   = iv%info(airsr)%i(1,n)
60       j   = iv%info(airsr)%j(1,n)
61       dx  = iv%info(airsr)%dx(1,n)
62       dy  = iv%info(airsr)%dy(1,n)
63       dxm = iv%info(airsr)%dxm(1,n)
64       dym = iv%info(airsr)%dym(1,n)
66       do k=kts,kte
67          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))
68          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))
69       end do
71       do k=1, num_levs
72          if (iv%airsr(n)%p(k) > 1.0) then
73             call da_to_zk(iv%airsr(n)%p(k), v_p, v_interp_p, iv%info(airsr)%zk(k,n))
74          else if (iv%airsr(n)%h(k) > 0.0) then
75             call da_to_zk(iv%airsr(n)%h(k), v_h, v_interp_h, iv%info(airsr)%zk(k,n))
76          end if
78          if (iv%info(airsr)%zk(k,n) < 0.0 .and. .not. anal_type_verify) then
79             iv%airsr(n)%t(k)%qc = missing_data
80             iv%airsr(n)%q(k)%qc = missing_data
81          end if
82       end do
83    end do
85    call da_convert_zk (iv%info(airsr))
87    ! [1.2] Interpolate horizontally to ob:
88    call da_interp_lin_3d (grid%xb%t, iv%info(airsr), model_t)
89    call da_interp_lin_3d (grid%xb%q, iv%info(airsr), model_q)
91    do n=iv%info(airsr)%n1,iv%info(airsr)%n2
92       !------------------------------------------------------------------------
93       ! [2.0] Initialise components of innovation vector:
94       !------------------------------------------------------------------------
96       do k = 1, iv%info(airsr)%levels(n)
97          iv%airsr(n)%t(k)%inv = 0.0
98          iv%airsr(n)%q(k)%inv = 0.0
100          !------------------------------------------------------------------------
101          ! [3.0] Interpolation:
102          !------------------------------------------------------------------------
104          if (ob%airsr(n)%t(k) > missing_r .AND. iv%airsr(n)%t(k)%qc >= obs_qc_pointer) then
105            iv%airsr(n)%t(k)%inv = ob%airsr(n)%t(k) - model_t(k,n)
106          end if
108          if (ob%airsr(n)%q(k) > missing_r .AND. iv%airsr(n)%q(k)%qc >= obs_qc_pointer) then
109             iv%airsr(n)%q(k)%inv = ob%airsr(n)%q(k) - model_q(k,n)
110          end if
111       end do
112    end do
114    !------------------------------------------------------------------------
115    ! [5.0] Perform optional maximum error check:
116    !------------------------------------------------------------------------
118    if ( check_max_iv ) &
119       call da_check_max_iv_airsr(iv, it, num_qcstat_conv)     
121    deallocate (model_t)
122    deallocate (model_q)
123    
124    if (trace_use_dull) call da_trace_exit("da_get_innov_vector_airsr")
126 end subroutine da_get_innov_vector_airsr