1 subroutine da_transform_vptox_inv(grid, vp, be, ep)
3 !-----------------------------------------------------------------------
4 ! Purpose: Inverse of balance (physical) transform of increment
6 ! Author: Zhiquan (Jake) Liu, NCAR/MMM, 2015-9
7 !-----------------------------------------------------------------------
11 type (domain), intent(inout) :: grid
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 !!! ??
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.
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 !--------------------------------------------------------------------------
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
53 if (cv_options /= 7) then
55 do j = grid%j_start(ij), grid%j_end(ij)
57 vp%v2(i,j,k) = vp%v2(i,j,k) - be%reg_psi_chi(j,k)* vp%v1(i,j,k)
65 if (cv_options /= 7) then ! - balance contri. from psi
68 do j = grid%j_start(ij), grid%j_end(ij)
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)
78 if ( cv_options == 6 ) then ! - balance contri. from Chi_u
81 do j = grid%j_start(ij), grid%j_end(ij)
83 vp%v3(i,j,k) = vp%v3(i,j,k) - be%reg_chi_u_t(j,k,k1)*vp%v2(i,j,k1)
91 !------------------------
92 !do j = grid%j_start(ij), grid%j_end(ij)
94 ! grid%xa%psfc(i,j) = vp%v5(i,j,1)
98 if (cv_options /= 7) then ! - balance contri. from psi
100 do j = grid%j_start(ij), grid%j_end(ij)
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)
109 if ( cv_options == 6 ) then ! - balance contri. from Chi_u
111 do j = grid%j_start(ij), grid%j_end(ij)
113 vp%v5(i,j,1) = vp%v5(i,j,1) - be%reg_chi_u_ps(j,k)*vp%v2(i,j,k)
119 ! 2.4 q --> pseudo rh=q/qs(background)
120 !----------------------------
121 ! if cv5 or cv7, no need for pseudo rh transform
124 ! do j = grid%j_start(ij), grid%j_end(ij)
126 ! vp%v4(i,j,k) = grid%xa % q(i,j,k) / grid%xb%qs(i,j,k)
131 if ( cv_options == 6 ) then
134 do j = grid%j_start(ij), grid%j_end(ij)
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)
146 do j = grid%j_start(ij), grid%j_end(ij)
148 vp%v4(i,j,k) = vp%v4(i,j,k) - be%reg_ps_u_rh(j,k)*vp%v5(i,j,1)
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)
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)
167 ! call da_add_flow_dependence_xa(grid, be % ne, ep, vp)
171 if (trace_use) call da_trace_exit("da_transform_vptox_inv")
173 end subroutine da_transform_vptox_inv