Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_vtox_transforms / da_transform_vptox_inv.inc
blob93649b675eb63a68350ecce79c89c0cd7be5c48f
1 subroutine da_transform_vptox_inv(grid, vp, be, ep)
3    !-----------------------------------------------------------------------
4    ! Purpose: Inverse of balance (physical) transform of increment
5    ! 
6    ! Author: Zhiquan (Jake) Liu,  NCAR/MMM,  2015-9
7    !-----------------------------------------------------------------------
9    implicit none
11    type (domain), intent(inout)         :: grid
12    
13    type (vp_type), intent(inout)        :: vp  ! input: full variables
14                                         ! output: unbalanced variables on model grid
15    type (be_type), intent(in), optional :: be  ! Background errors.
16    type (ep_type), intent(in), optional :: ep  ! Ensemble perturbations.
18    integer :: i, k, j, k1, ij            ! Loop counters.
19    real, allocatable          :: chi_u(:,:,:)  ! Unbalanced chi
21    if (trace_use) call da_trace_entry("da_transform_vptox_inv") 
23    !---------------------------------------------------------------------------
24    !  [1] Add flow-dependent increments in control variable space (vp):
25    !---------------------------------------------------------------------------
27    !if (be % ne > 0 .and. alphacv_method == alphacv_method_vp) then
28    !   call da_add_flow_dependence_vp(be % ne, ep, vp, its,ite, jts,jte, kts,kte)
29    !   call da_add_flow_dependence_vp_inv !!! ??
30    !end if
32    !--------------------------------------------------------------------------
33    ! [2] Impose statistical balance constraints:
34    ! Assume input vp%* is full variable, out vp% is unbalanced variables
35    ! to avoid (Psi,Chi) -> (U,V) transform, which has no exact inverse,
36    ! we need to store full variables at vp%* after each outloop.
37    ! da_transform_vptox.inc is also modified for this purpose.
38    !
39    ! for cv7, control variables are all full variables w/o multi-variate correlation.
40    ! so there is no need for balance transform and its inverse.
41    !--------------------------------------------------------------------------
43    !$OMP PARALLEL DO &
44    !$OMP PRIVATE ( ij, k1, k, j, i)
45    do ij = 1 , grid%num_tiles
47    ! 2.1 Psi, Chi --> Psi, Chi_u
48    !-------------------------
49    ! there is no need for Psi --> Psi transform
51    ! Chi --> Chi_u
52    !--------------------
53    if (cv_options /= 7) then
54       do k = kts, kte
55          do j = grid%j_start(ij), grid%j_end(ij)
56             do i = its, ite
57                vp%v2(i,j,k) = vp%v2(i,j,k) - be%reg_psi_chi(j,k)* vp%v1(i,j,k)
58             end do
59          end do
60       end do
61    end if
62   
63    ! 2.2 T --> T_u
64    !-------------------
65    if (cv_options /= 7) then ! - balance contri. from psi
66       do k1 = kts, kte
67          do k = kts, kte
68             do j = grid%j_start(ij), grid%j_end(ij)
69                do i = its, ite
70                   !vp%v3(i,j,k) = grid%xa%t(i,j,k) - be%reg_psi_t(j,k,k1)*vp%v1(i,j,k1)
71                   vp%v3(i,j,k) = vp%v3(i,j,k) - be%reg_psi_t(j,k,k1)*vp%v1(i,j,k1)
72                end do
73             end do
74          end do
75       end do
76    end if
78    if ( cv_options == 6 ) then ! - balance contri. from Chi_u 
79       do k1 = kts, kte
80          do k = kts, kte
81             do j = grid%j_start(ij), grid%j_end(ij)
82                do i = its, ite
83                   vp%v3(i,j,k) = vp%v3(i,j,k) - be%reg_chi_u_t(j,k,k1)*vp%v2(i,j,k1)
84                end do
85             end do
86          end do
87       end do
88    end if
90    ! 2.3 Ps --> Ps_u
91    !------------------------
92    !do j = grid%j_start(ij), grid%j_end(ij)
93    !   do i = its, ite
94    !      grid%xa%psfc(i,j) = vp%v5(i,j,1) 
95    !   end do
96    !end do
98    if (cv_options /= 7) then ! - balance contri. from psi
99       do k = kts,kte
100          do j = grid%j_start(ij), grid%j_end(ij)
101             do i = its, ite
102                !vp%v5(i,j,1) = grid%xa%psfc(i,j) - be%reg_psi_ps(j,k)*vp%v1(i,j,k)
103                vp%v5(i,j,1) = vp%v5(i,j,1) - be%reg_psi_ps(j,k)*vp%v1(i,j,k)
104             end do
105          end do
106       end do
107    end if
109    if ( cv_options == 6 ) then ! - balance contri. from Chi_u
110       do k = kts,kte
111          do j = grid%j_start(ij), grid%j_end(ij)
112             do i = its, ite
113                vp%v5(i,j,1) = vp%v5(i,j,1) - be%reg_chi_u_ps(j,k)*vp%v2(i,j,k)
114             end do
115          end do
116       end do
117    end if
119    ! 2.4 q --> pseudo rh=q/qs(background)
120    !----------------------------
121    ! if cv5 or cv7, no need for pseudo rh transform
123      !do k = kts, kte
124      !   do j = grid%j_start(ij), grid%j_end(ij)
125      !      do i = its, ite
126      !          vp%v4(i,j,k) =  grid%xa % q(i,j,k) / grid%xb%qs(i,j,k)
127      !      enddo
128      !   enddo
129      !enddo
131    if ( cv_options == 6 ) then
132       do k1 = kts, kte
133          do k = kts, kte
134             do j = grid%j_start(ij), grid%j_end(ij)
135                do i = its, ite
136                   vp%v4(i,j,k1) = vp%v4(i,j,k1) -  &
137                                   be%reg_psi_rh(j,k1,k)*vp%v1(i,j,k) - &
138                                 be%reg_chi_u_rh(j,k1,k)*vp%v2(i,j,k) - &
139                                   be%reg_t_u_rh(j,k1,k)*vp%v3(i,j,k)
140                end do
141             end do
142          end do
143       end do
145       do k = kts, kte
146          do j = grid%j_start(ij), grid%j_end(ij)
147             do i = its, ite
148                vp%v4(i,j,k) = vp%v4(i,j,k) - be%reg_ps_u_rh(j,k)*vp%v5(i,j,1)
149             end do
150          end do
151       end do
152    end if
154    end do
156    !---------------------------------------------------------------------------
157    !  [4] Add flow-dependent increments in model space (grid%xa):
158    !---------------------------------------------------------------------------
160 !  if (be % ne > 0 .and. alphacv_method == alphacv_method_xa) then
161 !     call da_add_flow_dependence_xa(grid, be % ne, ep, vp)
162 !  end if
163 !   if (be % ne > 0 .and. alphacv_method == alphacv_method_xa) then
164 !      if ( anal_type_hybrid_dual_res ) then
165 !         call da_add_flow_dependence_xa_dual_res(grid, be % ne, ep, vp)
166 !      else
167 !         call da_add_flow_dependence_xa(grid, be % ne, ep, vp)
168 !      endif
169 !   end if
171    if (trace_use) call da_trace_exit("da_transform_vptox_inv") 
173 end subroutine da_transform_vptox_inv