Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_airep / da_get_innov_vector_airep.inc
blob8a8812ce5195299bba0d634d197987b00720701e
1 subroutine da_get_innov_vector_airep( it,num_qcstat_conv, grid, ob, iv)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD   
5    !    Updated for Analysis on Arakawa-C grid
6    !    Author: Syed RH Rizvi,  MMM/ESSL/NCAR,  Date: 10/22/2008
7    !-----------------------------------------------------------------------
9    implicit none
11    integer,          intent(in)    :: it       ! External iteration.
12    type(domain),     intent(in)    :: grid     ! first guess state.
13    type(y_type),     intent(inout) :: ob       ! Observation structure.
14    type(iv_type),    intent(inout) :: iv       ! O-B structure.
15    integer,          intent(inout) :: num_qcstat_conv(:,:,:,:)
17    integer :: n                         ! Loop counter.
18    integer :: i  (kms:kme)
19    integer :: j  (kms:kme)
20    real    :: dx (kms:kme)
21    real    :: dxm(kms:kme)  
22    real    :: dy (kms:kme)
23    real    :: dym(kms:kme)  
24    integer :: k                   ! Index dimension.
25    real    :: speed, direction
27    real, allocatable :: model_u(:,:)  ! Model value u at ob location.
28    real, allocatable :: model_v(:,:)  ! Model value v at ob location.
29    real, allocatable :: model_t(:,:)  ! Model value t at ob location.
30    real, allocatable :: model_q(:,:)  ! Model value q at ob location.
32    real    :: v_h(kms:kme)      ! Model value h at ob hor. location.
33    real    :: v_p(kms:kme)      ! Model value p at ob hor. location.
35    if (trace_use_dull) call da_trace_entry("da_get_innov_vector_airep")
37    allocate (model_u(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
38    allocate (model_v(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
39    allocate (model_t(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
40    allocate (model_q(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
42    model_u(:,:) = 0.0
43    model_v(:,:) = 0.0
44    model_t(:,:) = 0.0
45    model_q(:,:) = 0.0
47    if ( it > 1) then
48       do n=iv%info(airep)%n1, iv%info(airep)%n2
49          do k=1, iv%info(airep)%levels(n)
50             if (iv%airep(n)%u(k)%qc == fails_error_max) iv%airep(n)%u(k)%qc = 0               
51             if (iv%airep(n)%v(k)%qc == fails_error_max) iv%airep(n)%v(k)%qc = 0               
52             if (iv%airep(n)%t(k)%qc == fails_error_max) iv%airep(n)%t(k)%qc = 0               
53             if (iv%airep(n)%q(k)%qc == fails_error_max) iv%airep(n)%q(k)%qc = 0               
54          end do
55       end do
56    end if
58    do n=iv%info(airep)%n1, iv%info(airep)%n2
59       if (iv%info(airep)%levels(n) < 1) cycle
61       ! [1.1] Get horizontal interpolation weights:
63       if (position_lev_dependant) then
64          i(:)   = iv%info(airep)%i(:,n)
65          j(:)   = iv%info(airep)%j(:,n)
66          dx(:)  = iv%info(airep)%dx(:,n)
67          dy(:)  = iv%info(airep)%dy(:,n)
68          dxm(:) = iv%info(airep)%dxm(:,n)
69          dym(:) = iv%info(airep)%dym(:,n)
70          do k=kts,kte
71             v_h(k) = dym(k)*(dxm(k)*grid%xb%h(i(k),j(k),  k) + dx(k)*grid%xb%h(i(k)+1,j(k),  k)) &
72                     + dy(k)*(dxm(k)*grid%xb%h(i(k),j(k)+1,k) + dx(k)*grid%xb%h(i(k)+1,j(k)+1,k))
73             v_p(k) = dym(k)*(dxm(k)*grid%xb%p(i(k),j(k),  k) + dx(k)*grid%xb%p(i(k)+1,j(k),  k)) &
74                     + dy(k)*(dxm(k)*grid%xb%p(i(k),j(k)+1,k) + dx(k)*grid%xb%p(i(k)+1,j(k)+1,k))
75          end do
76       else
77          i(1)   = iv%info(airep)%i(1,n)
78          j(1)   = iv%info(airep)%j(1,n)
79          dx(1)  = iv%info(airep)%dx(1,n)
80          dy(1)  = iv%info(airep)%dy(1,n)
81          dxm(1) = iv%info(airep)%dxm(1,n)
82          dym(1) = iv%info(airep)%dym(1,n)
84          v_h(kts:kte) = dym(1)*(dxm(1)*grid%xb%h(i(1),j(1),  kts:kte) + dx(1)*grid%xb%h(i(1)+1,j(1),  kts:kte)) &
85                        + dy(1)*(dxm(1)*grid%xb%h(i(1),j(1)+1,kts:kte) + dx(1)*grid%xb%h(i(1)+1,j(1)+1,kts:kte))
86          v_p(kts:kte) = dym(1)*(dxm(1)*grid%xb%p(i(1),j(1),  kts:kte) + dx(1)*grid%xb%p(i(1)+1,j(1),  kts:kte)) &
87                        + dy(1)*(dxm(1)*grid%xb%p(i(1),j(1)+1,kts:kte) + dx(1)*grid%xb%p(i(1)+1,j(1)+1,kts:kte))
88       end if
90       do k=1, iv%info(airep)%levels(n)
91          if (iv%airep(n)%p(k) > 1.0) then
92             call da_to_zk (iv%airep(n)%p(k), v_p, v_interp_p, iv%info(airep)%zk(k,n))
93          else if (iv%airep(n)%h(k) > 0.0) then
94             call da_to_zk (iv%airep(n)%h(k), v_h, v_interp_h, iv%info(airep)%zk(k,n))
95          end if
96       end do
97    end do
99    call da_convert_zk (iv%info(airep))
101    if (.not. anal_type_verify) then
102       do n=iv%info(airep)%n1,iv%info(airep)%n2
103          do k=1, iv%info(airep)%levels(n)
104             if (iv%info(airep)%zk(k,n) < 0.0) then
105                iv%airep(n)%u(k)%qc = missing_data
106                iv%airep(n)%v(k)%qc = missing_data
107                iv%airep(n)%t(k)%qc = missing_data
108                iv%airep(n)%q(k)%qc = missing_data
109             end if
110          end do
111       end do
112    end if
114 #ifdef A2C
115    call da_interp_lin_3d (grid%xb%u, iv%info(airep), model_u,'u')
116    call da_interp_lin_3d (grid%xb%v, iv%info(airep), model_v,'v')
117 #else
118    call da_interp_lin_3d (grid%xb%u, iv%info(airep), model_u)
119    call da_interp_lin_3d (grid%xb%v, iv%info(airep), model_v)
120 #endif
121    call da_interp_lin_3d (grid%xb%t, iv%info(airep), model_t)
122    call da_interp_lin_3d (grid%xb%q, iv%info(airep), model_q)
124    do n=iv%info(airep)%n1, iv%info(airep)%n2
126       !-------------------------------------------------------------------
127       ! [2.0] Initialise components of innovation vector:
128       !-------------------------------------------------------------------
130       do k = 1, iv%info(airep)%levels(n)
131          iv % airep(n) % u(k) % inv = 0.0
132          iv % airep(n) % v(k) % inv = 0.0
133          iv % airep(n) % t(k) % inv = 0.0
134          iv % airep(n) % q(k) % inv = 0.0
136          !----------------------------------------------------------------
137          ! [3.0] Fast interpolation:
138          !----------------------------------------------------------------
139          if (wind_sd_airep) then
140              call da_ffdduv_model (speed,direction,model_u(k,n), model_v(k,n), convert_uv2fd)
142              if (ob%airep(n)%u(k) > missing_r .AND. iv%airep(n)%u(k)%qc >= obs_qc_pointer) then
143                  iv%airep(n)%u(k)%inv = ob%airep(n)%u(k) - speed
144              end if
146              if (ob%airep(n)%v(k) > missing_r .AND. iv%airep(n)%v(k)%qc >= obs_qc_pointer) then
147                  iv%airep(n)%v(k)%inv = ob%airep(n)%v(k) - direction
148                  if (iv%airep(n)%v(k)%inv > 180.0 ) iv%airep(n)%v(k)%inv = iv%airep(n)%v(k)%inv - 360.0
149                  if (iv%airep(n)%v(k)%inv < -180.0 ) iv%airep(n)%v(k)%inv = iv%airep(n)%v(k)%inv + 360.0
150              end if
151          else
152              if (ob % airep(n) % u(k) > missing_r .AND. iv % airep(n) % u(k) % qc >= obs_qc_pointer) then
153                  iv % airep(n) % u(k) % inv = ob % airep(n) % u(k) - model_u(k,n)
154              end if
156              if (ob % airep(n) % v(k) > missing_r .AND. iv % airep(n) % v(k) % qc >= obs_qc_pointer) then
157                  iv % airep(n) % v(k) % inv = ob % airep(n) % v(k) - model_v(k,n)
158              end if
159          endif
161          if (ob % airep(n) % t(k) > missing_r .AND. iv % airep(n) % t(k) % qc >= obs_qc_pointer) then
162             iv % airep(n) % t(k) % inv = ob % airep(n) % t(k) - model_t(k,n)
163          end if
165          if (ob % airep(n) % q(k) > missing_r .AND. iv % airep(n) % q(k) % qc >= obs_qc_pointer) then
166             iv % airep(n) % q(k) % inv = ob % airep(n) % q(k) - model_q(k,n)
167          end if
168       end do
169    end do
171    !-------------------------------------------------------------------
172    ! [5.0] Perform optional maximum error check:
173    !-------------------------------------------------------------------
175    if ( check_max_iv ) &
176       call da_check_max_iv_airep (iv, it,num_qcstat_conv)
177    
178    deallocate (model_u)
179    deallocate (model_v)
180    deallocate (model_t)
181    deallocate (model_q)
182    
183    if (trace_use_dull) call da_trace_exit("da_get_innov_vector_airep")
185 end subroutine da_get_innov_vector_airep