Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_tamdar / da_get_innov_vector_tamdar.inc
blobf1b273af2d19cd3c5781465686bde8ff9e6e6382
1 subroutine da_get_innov_vector_tamdar (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, k        ! Loop counter.
16    integer :: i  (kms:kme)
17    integer :: j  (kms:kme)
18    real    :: dx (kms:kme)
19    real    :: dxm(kms:kme)  
20    real    :: dy (kms:kme)
21    real    :: dym(kms:kme)  
23    real, allocatable :: model_u(:,:)  ! Model value u at ob location.
24    real, allocatable :: model_v(:,:)  ! Model value v at ob location.
25    real, allocatable :: model_t(:,:)  ! Model value t at ob location.
26    real, allocatable :: model_q(:,:)  ! Model value q at ob location.
28    real    :: speed, direction
30    real    :: v_h(kts:kte)      ! Model value h at ob hor. location.
31    real    :: v_p(kts:kte)      ! Model value p at ob hor. location.
33    if (trace_use_dull) call da_trace_entry ("da_get_innov_vector_tamdar")
35    allocate (model_u(iv%info(tamdar)%max_lev,iv%info(tamdar)%n1:iv%info(tamdar)%n2))
36    allocate (model_v(iv%info(tamdar)%max_lev,iv%info(tamdar)%n1:iv%info(tamdar)%n2))
37    allocate (model_t(iv%info(tamdar)%max_lev,iv%info(tamdar)%n1:iv%info(tamdar)%n2))
38    allocate (model_q(iv%info(tamdar)%max_lev,iv%info(tamdar)%n1:iv%info(tamdar)%n2))
40    model_u(:,:) = 0.0
41    model_v(:,:) = 0.0
42    model_t(:,:) = 0.0
43    model_q(:,:) = 0.0
45    if ( it > 1 ) then
46       do n=iv%info(tamdar)%n1,iv%info(tamdar)%n2
47          do k=1, iv%info(tamdar)%levels(n)
48             if (iv%tamdar(n)%u(k)%qc == fails_error_max) iv%tamdar(n)%u(k)%qc = 0
49             if (iv%tamdar(n)%v(k)%qc == fails_error_max) iv%tamdar(n)%v(k)%qc = 0
50             if (iv%tamdar(n)%t(k)%qc == fails_error_max) iv%tamdar(n)%t(k)%qc = 0
51             if (iv%tamdar(n)%q(k)%qc == fails_error_max) iv%tamdar(n)%q(k)%qc = 0
52          end do
53       end do
54    end if
56    do n=iv%info(tamdar)%n1, iv%info(tamdar)%n2
57       if (iv%info(tamdar)%levels(n) < 1) cycle
59       ! [1.1] Get horizontal interpolation weights:
61       if (position_lev_dependant) then
62          i(:)   = iv%info(tamdar)%i(:,n)
63          j(:)   = iv%info(tamdar)%j(:,n)
64          dx(:)  = iv%info(tamdar)%dx(:,n)
65          dy(:)  = iv%info(tamdar)%dy(:,n)
66          dxm(:) = iv%info(tamdar)%dxm(:,n)
67          dym(:) = iv%info(tamdar)%dym(:,n)
68          do k=kts,kte
69             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)) &
70                + 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))
71             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)) &
72                + 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))
73          end do
74       else
75          i(1)   = iv%info(tamdar)%i(1,n)
76          j(1)   = iv%info(tamdar)%j(1,n)
77          dx(1)  = iv%info(tamdar)%dx(1,n)
78          dy(1)  = iv%info(tamdar)%dy(1,n)
79          dxm(1) = iv%info(tamdar)%dxm(1,n)
80          dym(1) = iv%info(tamdar)%dym(1,n)
82          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)) &
83                        + 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))
84          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)) &
85                        + 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))
86       end if
88       do k=1, iv%info(tamdar)%levels(n)
89          if (iv%tamdar(n)%p(k) > 1.0) then
90             call da_to_zk (iv%tamdar(n)%p(k), v_p, v_interp_p, iv%info(tamdar)%zk(k,n))
91          else if (iv%tamdar(n)%h(k) > 0.0) then
92             call da_to_zk (iv%tamdar(n)%h(k), v_h, v_interp_h, iv%info(tamdar)%zk(k,n))
93          end if
94       end do
96    end do
98    call da_convert_zk (iv%info(tamdar))
100    if (.not. anal_type_verify) then
101       do n=iv%info(tamdar)%n1,iv%info(tamdar)%n2
102          do k=1, iv%info(tamdar)%levels(n)
103             if (iv%info(tamdar)%zk(k,n) < 0.0) then
104                iv%tamdar(n)%u(k)%qc = missing_data
105                iv%tamdar(n)%v(k)%qc = missing_data
106                iv%tamdar(n)%t(k)%qc = missing_data
107                iv%tamdar(n)%q(k)%qc = missing_data
108             end if
109          end do
110       end do
111    end if
113    ! [1.2] Interpolate horizontally to ob:
114 #ifdef A2C
115    call da_interp_lin_3d (grid%xb%u, iv%info(tamdar), model_u,'u')
116    call da_interp_lin_3d (grid%xb%v, iv%info(tamdar), model_v,'v')
117 #else
118    call da_interp_lin_3d (grid%xb%u, iv%info(tamdar), model_u)
119    call da_interp_lin_3d (grid%xb%v, iv%info(tamdar), model_v)
120 #endif
121    call da_interp_lin_3d (grid%xb%t, iv%info(tamdar), model_t)
122    call da_interp_lin_3d (grid%xb%q, iv%info(tamdar), model_q)
124    do n=iv%info(tamdar)%n1, iv%info(tamdar)%n2
125       !----------------------------------------------------------------------
126       ! [2.0] Initialise components of innovation vector:
127       !----------------------------------------------------------------------
129       do k = 1, iv%info(tamdar)%levels(n)
130          iv%tamdar(n)%u(k)%inv = 0.0
131          iv%tamdar(n)%v(k)%inv = 0.0
132          iv%tamdar(n)%t(k)%inv = 0.0
133          iv%tamdar(n)%q(k)%inv = 0.0
135          !-------------------------------------------------------------------
136          ! [3.0] Interpolation:
137          !-----------------------------------------------------------------
138           if (wind_sd_tamdar) then
139               call da_ffdduv_model (speed,direction,model_u(k,n), model_v(k,n), convert_uv2fd)
141               if (ob%tamdar(n)%u(k) > missing_r .AND. iv%tamdar(n)%u(k)%qc >= obs_qc_pointer) then
142                   iv%tamdar(n)%u(k)%inv = ob%tamdar(n)%u(k) - speed
143               end if
145               if (ob%tamdar(n)%v(k) > missing_r .AND. iv%tamdar(n)%v(k)%qc >= obs_qc_pointer) then
146                   iv%tamdar(n)%v(k)%inv = ob%tamdar(n)%v(k) - direction
147                   if (iv%tamdar(n)%v(k)%inv > 180.0 ) iv%tamdar(n)%v(k)%inv = iv%tamdar(n)%v(k)%inv - 360.0
148                   if (iv%tamdar(n)%v(k)%inv < -180.0 ) iv%tamdar(n)%v(k)%inv = iv%tamdar(n)%v(k)%inv + 360.0
149               end if
150           else
151              if (ob%tamdar(n)%u(k) > missing_r .AND. iv%tamdar(n)%u(k)%qc >= obs_qc_pointer) then
152                  iv%tamdar(n)%u(k)%inv = ob%tamdar(n)%u(k) - model_u(k,n)
153              end if
155              if (ob%tamdar(n)%v(k) > missing_r .AND. iv%tamdar(n)%v(k)%qc >= obs_qc_pointer) then
156                  iv%tamdar(n)%v(k)%inv = ob%tamdar(n)%v(k) - model_v(k,n)
157              end if
158           end if
161              if (ob%tamdar(n)%t(k) > missing_r .AND. iv%tamdar(n)%t(k)%qc >= obs_qc_pointer) then
162                  iv%tamdar(n)%t(k)%inv = ob%tamdar(n)%t(k) - model_t(k,n)
163              end if
165              if (ob%tamdar(n)%q(k) > missing_r .AND. iv%tamdar(n)%q(k)%qc >= obs_qc_pointer) then
166                  iv%tamdar(n)%q(k)%inv = ob%tamdar(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) call da_check_max_iv_tamdar (iv, it, num_qcstat_conv)
177    deallocate (model_u)
178    deallocate (model_v)
179    deallocate (model_t)
180    deallocate (model_q)
181    if (trace_use_dull) call da_trace_exit ("da_get_innov_vector_tamdar")
183 end subroutine da_get_innov_vector_tamdar