Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_sound / da_get_innov_vector_sound.inc
blob54cab07b7255e34c2d5414b8acba9671b59e6185
1 subroutine da_get_innov_vector_sound (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, k        ! 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)  
25    real, allocatable :: model_u(:,:)  ! Model value u at ob location.
26    real, allocatable :: model_v(:,:)  ! Model value v at ob location.
27    real, allocatable :: model_t(:,:)  ! Model value t at ob location.
28    real, allocatable :: model_q(:,:)  ! Model value q at ob location.
30    real    :: speed, direction
32    real    :: v_h(kts:kte)      ! Model value h at ob hor. location.
33    real    :: v_p(kts:kte)      ! Model value p at ob hor. location.
35    if (trace_use_dull) call da_trace_entry ("da_get_innov_vector_sound")
37    allocate (model_u(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
38    allocate (model_v(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
39    allocate (model_t(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
40    allocate (model_q(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%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(sound)%n1,iv%info(sound)%n2
49         do k=1, iv%info(sound)%levels(n)
50            if (iv%sound(n)%u(k)%qc == fails_error_max) iv%sound(n)%u(k)%qc = 0
51            if (iv%sound(n)%v(k)%qc == fails_error_max) iv%sound(n)%v(k)%qc = 0
52            if (iv%sound(n)%t(k)%qc == fails_error_max) iv%sound(n)%t(k)%qc = 0
53            if (iv%sound(n)%q(k)%qc == fails_error_max) iv%sound(n)%q(k)%qc = 0
54         end do
55      end do
56    end if
58    do n=iv%info(sound)%n1, iv%info(sound)%n2
59       if (iv%info(sound)%levels(n) < 1) cycle
61       ! [1.1] Get horizontal interpolation weights:
63       if (position_lev_dependant) then
64          i(:)   = iv%info(sound)%i(:,n)
65          j(:)   = iv%info(sound)%j(:,n)
66          dx(:)  = iv%info(sound)%dx(:,n)
67          dy(:)  = iv%info(sound)%dy(:,n)
68          dxm(:) = iv%info(sound)%dxm(:,n)
69          dym(:) = iv%info(sound)%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(sound)%i(1,n)
78          j(1)   = iv%info(sound)%j(1,n)
79          dx(1)  = iv%info(sound)%dx(1,n)
80          dy(1)  = iv%info(sound)%dy(1,n)
81          dxm(1) = iv%info(sound)%dxm(1,n)
82          dym(1) = iv%info(sound)%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(sound)%levels(n)
91          if (iv%sound(n)%p(k) > 1.0) then
92             call da_to_zk (iv%sound(n)%p(k), v_p, v_interp_p, iv%info(sound)%zk(k,n))
93          else if (iv%sound(n)%h(k) > 0.0) then
94             call da_to_zk (iv%sound(n)%h(k), v_h, v_interp_h, iv%info(sound)%zk(k,n))
95          end if
97       end do
99    end do
101    call da_convert_zk (iv%info(sound))
103    if (.not. anal_type_verify) then
104       do n=iv%info(sound)%n1,iv%info(sound)%n2
105          do k=1, iv%info(sound)%levels(n)
106             if (iv%info(sound)%zk(k,n) < 0.0) then
107                iv%sound(n)%u(k)%qc = missing_data
108                iv%sound(n)%v(k)%qc = missing_data
109                iv%sound(n)%t(k)%qc = missing_data
110                iv%sound(n)%q(k)%qc = missing_data
111             end if
112          end do
113       end do
114    end if
116    ! [1.2] Interpolate horizontally to ob:
118 #ifdef A2C
119    call da_interp_lin_3d (grid%xb%u, iv%info(sound), model_u,'u')
120    call da_interp_lin_3d (grid%xb%v, iv%info(sound), model_v,'v')
121 #else
122    call da_interp_lin_3d (grid%xb%u, iv%info(sound), model_u)
123    call da_interp_lin_3d (grid%xb%v, iv%info(sound), model_v)
124 #endif
125    call da_interp_lin_3d (grid%xb%t, iv%info(sound), model_t)
126    call da_interp_lin_3d (grid%xb%q, iv%info(sound), model_q)
128    do n=iv%info(sound)%n1, iv%info(sound)%n2
129       !----------------------------------------------------------------------
130       ! [2.0] Initialise components of innovation vector:
131       !----------------------------------------------------------------------
133       do k = 1, iv%info(sound)%levels(n)
134          iv%sound(n)%u(k)%inv = 0.0
135          iv%sound(n)%v(k)%inv = 0.0
136          iv%sound(n)%t(k)%inv = 0.0
137          iv%sound(n)%q(k)%inv = 0.0
138          !-------------------------------------------------------------------
139          ! [3.0] Interpolation:
140          !-------------------------------------------------------------------
142           if (wind_sd_sound) then
143               call da_ffdduv_model (speed,direction,model_u(k,n), model_v(k,n), convert_uv2fd)
145               if (ob%sound(n)%u(k) > missing_r .AND. iv%sound(n)%u(k)%qc >= obs_qc_pointer) then
146                   iv%sound(n)%u(k)%inv = ob%sound(n)%u(k) - speed
147               end if
149               if (ob%sound(n)%v(k) > missing_r .AND. iv%sound(n)%v(k)%qc >= obs_qc_pointer) then
150                   iv%sound(n)%v(k)%inv = ob%sound(n)%v(k) - direction
151                   if (iv%sound(n)%v(k)%inv > 180.0 ) iv%sound(n)%v(k)%inv = iv%sound(n)%v(k)%inv - 360.0
152                   if (iv%sound(n)%v(k)%inv < -180.0 ) iv%sound(n)%v(k)%inv = iv%sound(n)%v(k)%inv + 360.0
153               end if
154            else
155              if (ob%sound(n)%u(k) > missing_r .AND. iv%sound(n)%u(k)%qc >= obs_qc_pointer) then
156                  iv%sound(n)%u(k)%inv = ob%sound(n)%u(k) - model_u(k,n)
157              end if
159              if (ob%sound(n)%v(k) > missing_r .AND. iv%sound(n)%v(k)%qc >= obs_qc_pointer) then
160                  iv%sound(n)%v(k)%inv = ob%sound(n)%v(k) - model_v(k,n)
161              end if
162           end if
164          if (ob%sound(n)%t(k) > missing_r .AND. iv%sound(n)%t(k)%qc >= obs_qc_pointer) then
165             iv%sound(n)%t(k)%inv = ob%sound(n)%t(k) - model_t(k,n)
166          end if
168          if (ob%sound(n)%q(k) > missing_r .AND. iv%sound(n)%q(k)%qc >= obs_qc_pointer) then
169             iv%sound(n)%q(k)%inv = ob%sound(n)%q(k) - model_q(k,n)
170          end if
171       end do
172    end do
174    !----------------------------------------------------------------------
175    ! [5.0] Perform optional maximum error check:
176    !----------------------------------------------------------------------
178    if ( check_max_iv ) &
179       call da_check_max_iv_sound (iv, it, num_qcstat_conv)
181    if (check_buddy) call da_check_buddy_sound(iv, grid%dx, it)
183    deallocate (model_u)
184    deallocate (model_v)
185    deallocate (model_t)
186    deallocate (model_q)
187    
188    if (trace_use_dull) call da_trace_exit ("da_get_innov_vector_sound")
190 end subroutine da_get_innov_vector_sound