68 integer,
intent(in) :: KA, KS, KE
69 integer,
intent(in) :: IA, IS, IE
70 integer,
intent(in) :: JA, JS, JE
72 real(RP),
intent(in) :: DENS(ka,ia,ja)
73 real(RP),
intent(in) :: MOMZ(ka,ia,ja)
74 real(RP),
intent(in) :: MOMX(ka,ia,ja)
75 real(RP),
intent(in) :: MOMY(ka,ia,ja)
77 real(RP),
intent(out) :: W (ka,ia,ja)
78 real(RP),
intent(out) :: U (ka,ia,ja)
79 real(RP),
intent(out) :: V (ka,ia,ja)
92 w(k,i,j) = 0.5_rp * ( momz(k-1,i,j)+momz(k,i,j) ) / dens(k,i,j)
102 momws = momz(ks,i,j) &
103 + ( j13g(ks,i,j,
i_xyw) * ( momx(ks,i,j) + momx(ks,i-1,j) + momx(ks+1,i,j) + momx(ks+1,i-1,j) ) &
104 + j23g(ks,i,j,
i_xyw) * ( momy(ks,i,j) + momy(ks,i,j-1) + momy(ks+1,i,j) + momy(ks+1,i,j-1) ) ) &
105 * 0.25_rp / gsqrt(ks,i,j,
i_xyw)
108 w(ks,i,j) = momws * 0.5_rp &
109 - ( j13g(ks,i,j,
i_xyz) * ( momx(ks,i,j) + momx(ks,i-1,j) ) &
110 + j23g(ks,i,j,
i_xyz) * ( momy(ks,i,j) + momy(ks,i,j-1) ) ) &
111 * 0.5_rp / ( dens(ks,i,j) * gsqrt(ks,i,j,
i_xyz) )
118 w(ke,i,j) = 0.5_rp * ( momz(ke-1,i,j) ) / dens(ke,i,j)
127 u(k,i,j) = 0.5_rp * ( momx(k,i-1,j)+momx(k,i,j) ) / dens(k,i,j)
135 u(k,1,j) = momx(k,1,j) / dens(k,1,j)
144 v(k,i,j) = 0.5_rp * ( momy(k,i,j-1)+momy(k,i,j) ) / dens(k,i,j)
152 v(k,i,1) = momy(k,i,1) / dens(k,i,1)
159 w( 1:ks-1,i,j) = w(ks,i,j)
160 u( 1:ks-1,i,j) = u(ks,i,j)
161 v( 1:ks-1,i,j) = v(ks,i,j)
162 w(ke+1:ka, i,j) = w(ke,i,j)
163 u(ke+1:ka, i,j) = u(ke,i,j)
164 v(ke+1:ka, i,j) = v(ke,i,j)
168 call comm_vars8( w(:,:,:), 1 )
169 call comm_vars8( u(:,:,:), 2 )
170 call comm_vars8( v(:,:,:), 3 )
171 call comm_wait ( w(:,:,:), 1, .false. )
172 call comm_wait ( u(:,:,:), 2, .false. )
173 call comm_wait ( v(:,:,:), 3, .false. )
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_gsqrt
transformation metrics from Z to Xi, {G}^1/2
module Atmosphere Grid CartesianC metirc
subroutine, public atmos_diagnostic_cartesc_get_vel(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, MOMZ, MOMX, MOMY, W, U, V)
ATMOS_DIAGNOSTIC_CARTESC_get_vel W, U, V.
module atmosphere / grid / cartesC index
module atmosphere / diagnostic / CartesianC
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j23g
(2,3) element of Jacobian matrix * {G}^1/2
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j13g
(1,3) element of Jacobian matrix * {G}^1/2