1 SUBROUTINE da_transform_bal_adj( vp, be, grid )
5 TYPE (vp_type), INTENT(INOUT) :: vp ! CV on grid structure.out
6 TYPE (be_type), INTENT(IN) :: be ! Background errors.
7 type (domain) , intent(inout) :: grid ! Dimensions and xpose buffers.
9 INTEGER :: i, j, k, kk, ij, iunit ! Loop counters.
10 !-------------------------------------------------------------------
12 !-------------------------------------------------------------------
14 ! linear balance btw psi and t-b, Psfc_b and chi_b
17 ! [3.4] adj of Transform psi and chi to u and v:
20 #include "HALO_PSICHI_UV_ADJ.inc"
24 !$OMP PRIVATE ( i, j, k, ij )
25 DO ij = 1 , grid%num_tiles
27 DO j = grid%j_start(ij), grid%j_end(ij)
37 call da_psichi_to_uv_adj( grid%xa % u, grid%xa % v, grid%xb % coefx, &
38 grid%xb % coefy, vp % v1, vp % v2 )
40 ! [3.3] adj Calculate Psfc_b from psi
42 !--convert from delt.ps to delt.ln(ps)
44 !$OMP PRIVATE ( i, j, ij )
45 DO ij = 1 , grid%num_tiles
46 DO j = grid%j_start(ij), grid%j_end(ij)
48 grid%xa%psfc(i,j) = grid%xa%psfc(i,j) * grid%xb%psfc(i,j)
55 !$OMP PRIVATE ( i, j, k, ij )
56 DO ij = 1 , grid%num_tiles
58 DO j = grid%j_start(ij), grid%j_end(ij)
60 vp % v1(i,j,k)= vp % v1(i,j,k)+ &
61 be % wgvz(i,j,k) * grid%xa % psfc(i,j)
69 !$OMP PRIVATE ( i, j, ij )
70 DO ij = 1 , grid%num_tiles
71 DO j = grid%j_start(ij), grid%j_end(ij)
73 vp % v5(i,j,1)=grid%xa % psfc(i,j)
79 ! [3.2] adj Calculate chi_b from psi
82 !$OMP PRIVATE ( i, j, k, ij )
83 DO ij = 1 , grid%num_tiles
85 DO j = grid%j_start(ij), grid%j_end(ij)
87 vp % v1(i,j,k)= vp % v1(i,j,k)+ &
88 be % bvz(i,j,k) * vp % v2(i,j,k)
95 ! [3.1] Calculate t_b from psi
98 !$OMP PRIVATE ( i, j, k, kk, ij )
99 DO ij = 1 , grid%num_tiles
102 DO j = grid%j_start(ij), grid%j_end(ij)
104 vp % v1(i,j,kk)= vp % v1(i,j,kk)+ &
105 be % agvz(i,j,k,kk) * grid%xa % t(i,j,k)
111 !$OMP END PARALLEL DO
114 !$OMP PRIVATE ( i, j, k, ij )
115 DO ij = 1 , grid%num_tiles
117 DO j = grid%j_start(ij), grid%j_end(ij)
119 vp % v3(i,j,k)= grid%xa % t(i,j,k)
124 !$OMP END PARALLEL DO
126 ! [3.5] treat humidity
128 IF ( cv_options == 3 ) THEN
130 !$OMP PRIVATE ( i, j, k, ij )
131 DO ij = 1 , grid%num_tiles
133 DO j = grid%j_start(ij), grid%j_end(ij)
135 vp % v4(i,j,k) = grid%xa % q(i,j,k) * grid%xb % qs(i,j,k)
140 !$OMP END PARALLEL DO
142 ELSE IF ( cv_options_hum == 1 ) THEN
144 vp % v4(its:ite,jts:jte,kts:kte) = vp % v4(its:ite,jts:jte,kts:kte)+&
145 grid%xa % q(its:ite,jts:jte,kts:kte)
147 ELSE IF ( cv_options_hum == 2 ) THEN
149 CALL DA_TPRH_To_Q_Adj( grid )
151 vp % v4(its:ite,jts:jte,kts:kte) = vp % v4(its:ite,jts:jte,kts:kte)+&
152 grid%xa % rh(its:ite,jts:jte,kts:kte)
156 END SUBROUTINE da_transform_bal_adj