48 use constants_mod
, only: rdgas, cp_air, grav
59 real,
parameter::
r3 = 1./3.
63 subroutine riem_solver3(ms, dt, is, ie, js, je, km, ng, &
64 isd, ied, jsd, jed, akap, cappa, cp, &
68 ptop, zs, q_con, w, delz, pt, &
69 delp, zh, pe, ppe, pk3, pk, peln, &
70 ws, scale_m, p_fac, a_imp, &
71 use_logp, last_call, fp_out)
78 integer,
intent(in):: ms, is, ie, js, je, km, ng
79 integer,
intent(in):: isd, ied, jsd, jed
81 real,
intent(in):: akap, cp, ptop, p_fac, a_imp, scale_m
82 real,
intent(in):: zs(isd:ied,jsd:jed)
83 logical,
intent(in):: last_call, use_logp, fp_out
84 real,
intent(in):: ws(is:ie,js:je)
85 real,
intent(in),
dimension(isd:,jsd:,1:):: q_con, cappa
87 real,
intent(in),
dimension(isd:ied,jsd:jed,km):: kapad
89 real,
intent(in),
dimension(isd:ied,jsd:jed,km):: delp, pt
90 real,
intent(inout),
dimension(isd:ied,jsd:jed,km+1):: zh
91 real,
intent(inout),
dimension(isd:ied,jsd:jed,km):: w
92 real,
intent(inout):: pe(is-1:ie+1,km+1,js-1:je+1)
93 real,
intent(out):: peln(is:ie,km+1,js:je)
94 real,
intent(out),
dimension(isd:ied,jsd:jed,km+1):: ppe
95 real,
intent(out):: delz(is:ie,js:je,km)
96 real,
intent(out):: pk(is:ie,js:je,km+1)
97 real,
intent(out):: pk3(isd:ied,jsd:jed,km+1)
99 real,
dimension(is:ie,km):: dm, dz2, pm2, w2, gm2, cp2
100 real,
dimension(is:ie,km+1)::pem, pe2, peln2, peg, pelng
102 real,
dimension(is:ie,km):: kapad2
104 real gama, rgrav, ptk, peln1
110 ptk = exp(akap*peln1)
125 dm(i,k) = delp(i,j,k)
127 cp2(i,k) = cappa(i,j,k)
130 kapad2(i,k) = kapad(i,j,k)
146 pem(i,k) = pem(i,k-1) + dm(i,k-1)
147 peln2(i,k) = log(pem(i,k))
151 peg(i,k) = peg(i,k-1) + dm(i,k-1)*(1.-q_con(i,j,k-1))
152 pelng(i,k) = log(peg(i,k))
155 pk3(i,j,k) = exp(akap*peln2(i,k))
162 pm2(i,k) = (peg(i,k+1)-peg(i,k))/(pelng(i,k+1)-pelng(i,k))
165 gm2(i,k) = 1. / (1.-cp2(i,k))
169 pm2(i,k) = dm(i,k)/(peln2(i,k+1)-peln2(i,k))
171 dm(i,k) = dm(i,k) * rgrav
172 dz2(i,k) = zh(i,j,k+1) - zh(i,j,k)
178 if ( a_imp < -0.999 )
then 179 call sim3p0_solver(dt, is, ie, km, rdgas, gama, akap, &
184 pem, w2, dz2, pt(is:ie,j,1:km), ws(is,j), p_fac, scale_m )
185 elseif ( a_imp < -0.5 )
then 186 call sim3_solver(dt, is, ie, km, rdgas, gama, akap, &
191 pem, w2, dz2, pt(is:ie,j,1:km), ws(is,j), abs(a_imp), p_fac, scale_m)
192 elseif ( a_imp <= 0.5 )
then 193 call rim_2d(ms, dt, is, ie, km, rdgas, gama, gm2, &
198 dm, pm2, w2, dz2, pt(is:ie,j,1:km), ws(is,j), .false.)
199 elseif ( a_imp > 0.999 )
then 200 call sim1_solver(dt, is, ie, km, rdgas, gama, gm2, cp2, akap, &
205 pm2, pem, w2, dz2, pt(is:ie,j,1:km), ws(is,j), p_fac)
207 call sim_solver(dt, is, ie, km, rdgas, gama, gm2, cp2, akap, &
212 pm2, pem, w2, dz2, pt(is:ie,j,1:km), ws(is,j), &
213 a_imp, p_fac, scale_m)
220 delz(i,j,k) = dz2(i,k)
224 if ( last_call )
then 227 peln(i,k,j) = peln2(i,k)
228 pk(i,j,k) = pk3(i,j,k)
237 ppe(i,j,k) = pe2(i,k) + pem(i,k)
243 ppe(i,j,k) = pe2(i,k)
251 pk3(i,j,k) = peln2(i,k)
257 zh(i,j,km+1) = zs(i,j)
261 zh(i,j,k) = zh(i,j,k+1) - dz2(i,k)
subroutine, public sim_solver(dt, is, ie, km, rgas, gama, gm2, cp2, kappa, pe2, dm2, pm2, pem, w2, dz2, pt2, ws, alpha, p_fac, scale_m)
subroutine, public sim3_solver(dt, is, ie, km, rgas, gama, kappa, pe2, dm, pem, w2, dz2, pt2, ws, alpha, p_fac, scale_m)
subroutine, public fv_tp_2d(q, crx, cry, npx, npy, hord, fx, fy, xfx, yfx, gridstruct, bd, ra_x, ra_y, lim_fac, mfx, mfy, mass, nord, damp_c)
The subroutine 'fv_tp_2d' contains the FV advection scheme .
subroutine, public nh_bc(ptop, grav, kappa, cp, delp, delzBC, pt, phis, pkc, gz, pk3, BC_step, BC_split, npx, npy, npz, bounded_domain, pkc_pertn, computepk3, fullhalo, bd)
The module 'nh_utils' peforms non-hydrostatic computations.
The module 'tp_core' is a collection of routines to support FV transport.
subroutine, public update_dz_c(is, ie, js, je, km, ng, dt, dp0, zs, area, ut, vt, gz, ws, npx, npy, sw_corner, se_corner, ne_corner, nw_corner, bd, grid_type)
subroutine, public riem_solver3(ms, dt, is, ie, js, je, km, ng, isd, ied, jsd, jed, akap, cappa, cp, ptop, zs, q_con, w, delz, pt, delp, zh, pe, ppe, pk3, pk, peln, ws, scale_m, p_fac, a_imp, use_logp, last_call, fp_out)
subroutine, public sim1_solver(dt, is, ie, km, rgas, gama, gm2, cp2, kappa, pe, dm2, pm2, pem, w2, dz2, pt2, ws, p_fac)
subroutine, public sim3p0_solver(dt, is, ie, km, rgas, gama, kappa, pe2, dm, pem, w2, dz2, pt2, ws, p_fac, scale_m)
subroutine, public rim_2d(ms, bdt, is, ie, km, rgas, gama, gm2, pe2, dm2, pm2, w2, dz2, pt2, ws, c_core)
The module 'nh_core' peforms non-hydrostatic computations include moisture effect in pt...
subroutine, public update_dz_d(ndif, damp, hord, is, ie, js, je, km, ng, npx, npy, area, rarea, dp0, zs, zh, crx, cry, xfx, yfx, ws, rdt, gridstruct, bd, lim_fac)
subroutine, public riem_solver_c(ms, dt, is, ie, js, je, km, ng, akap, cappa, cp, ptop, hs, w3, pt, q_con, delp, gz, pef, ws, p_fac, a_imp, scale_m)