59 use constants_mod
, only: rdgas, rvgas, cp_air, cp_vapor, hlv, hlf, kappa, grav
60 use tracer_manager_mod
, only: get_tracer_index
61 use field_manager_mod
, only: model_atmos
63 use gfdl_cloud_microphys_mod
, only: wqs1, wqs2, wqsat2_moist
65 use fv_mp_mod, only: mp_reduce_min, is_master
74 real,
parameter::
esl = 0.621971831
75 real,
parameter::
tice = 273.16
78 real,
parameter::
c_liq = 4.1855e+3
80 real,
parameter::
cv_vap = cp_vapor - rvgas
87 real,
parameter::
hlv0 = 2.5e6
88 real,
parameter::
hlf0 = 3.3358e5
91 real,
parameter::
t_ice = 273.16
101 real,
parameter::
zvir = rvgas/rdgas - 1.
113 subroutine fv_subgrid_z( isd, ied, jsd, jed, is, ie, js, je, km, nq, dt, &
114 tau, nwat, delp, pe, peln, pkz, ta, qa, ua, va, &
115 hydrostatic, w, delz, u_dt, v_dt, t_dt, k_bot )
118 integer,
intent(in):: is, ie, js, je, km, nq, nwat
119 integer,
intent(in):: isd, ied, jsd, jed
120 integer,
intent(in):: tau
121 real,
intent(in):: dt
122 real,
intent(in):: pe(is-1:ie+1,km+1,js-1:je+1)
123 real,
intent(in):: peln(is :ie, km+1,js :je)
124 real,
intent(in):: delp(isd:ied,jsd:jed,km)
125 real,
intent(in):: delz(is:,js:,1:)
126 real,
intent(in):: pkz(is:ie,js:je,km)
127 logical,
intent(in):: hydrostatic
128 integer,
intent(in),
optional:: k_bot
130 real,
intent(inout):: ua(isd:ied,jsd:jed,km)
131 real,
intent(inout):: va(isd:ied,jsd:jed,km)
132 real,
intent(inout):: w(isd:,jsd:,1:)
133 real,
intent(inout):: ta(isd:ied,jsd:jed,km)
134 real,
intent(inout):: qa(isd:ied,jsd:jed,km,nq)
135 real,
intent(inout):: u_dt(isd:ied,jsd:jed,km)
136 real,
intent(inout):: v_dt(isd:ied,jsd:jed,km)
137 real,
intent(inout):: t_dt(is:ie,js:je,km)
139 real,
dimension(is:ie,km):: u0, v0, w0, t0, hd, te, gz, tvm, pm, den
140 real q0(is:ie,km,nq), qcon(is:ie,km)
141 real,
dimension(is:ie):: gzh, lcp2, icp2, cvm, cpm, qs
143 real :: rkx, rdx, rzx, c_air
145 real ri_ref, ri, pt1, pt2, ratio, tv, cv, tmp, q_liq, q_sol
146 real tv1, tv2, g2, h0, mc, fra, rk, rz, rdt, tvd, tv_surf
147 real dh, dq, qsw, dqsdt, tcp3, t_max, t_min
148 integer i, j, k, kk, n, m, iq, km1, im, kbot, l
149 real,
parameter:: ustar2 = 1.e-4
151 integer :: sphum, liq_wat, rainwat, snowwat, graupel, ice_wat, cld_amt
153 cv_air = cp_air - rdgas
154 rk = cp_air/rdgas + 1.
162 if (
present(k_bot) )
then 163 if ( k_bot < 3 )
return 168 if ( pe(is,1,js) < 2. )
then 174 if ( k_bot < min(km,24) )
then 180 sphum = get_tracer_index(model_atmos,
'sphum')
181 if ( nwat == 0 )
then 187 liq_wat = get_tracer_index(model_atmos,
'liq_wat')
188 ice_wat = get_tracer_index(model_atmos,
'ice_wat')
189 rainwat = get_tracer_index(model_atmos,
'rainwat')
190 snowwat = get_tracer_index(model_atmos,
'snowwat')
191 graupel = get_tracer_index(model_atmos,
'graupel')
192 cld_amt = get_tracer_index(model_atmos,
'cld_amt')
221 q0(i,k,iq) = qa(i,j,k,iq)
230 tvm(i,k) = t0(i,k)*
virq(q0(i,k,:))
232 tvm(i,k) = t0(i,k)*(1.+xvir*q0(i,k,sphum))
236 pm(i,k) = delp(i,j,k)/(peln(i,k+1,j)-peln(i,k,j))
244 if( hydrostatic )
then 248 den(i,k) = pm(i,k)/tv
249 gz(i,k) = gzh(i) + tv*(1.-pe(i,k,j)/pm(i,k))
250 hd(i,k) = cp_air*tvm(i,k)+gz(i,k)+0.5*(u0(i,k)**2+v0(i,k)**2)
251 gzh(i) = gzh(i) + tv*(peln(i,k+1,j)-peln(i,k,j))
256 if ( nwat == 0 )
then 259 cpm(i) = cp_air*
vicpqd(q0(i,k,:))
260 cvm(i) = cv_air*
vicvqd(q0(i,k,:))
266 elseif ( nwat==1 )
then 269 cpm(i) = (1.-q0(i,k,sphum))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor
270 cvm(i) = (1.-q0(i,k,sphum))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap 272 cpm(i) = (1.-q0(i,k,sphum))*cp_air + q0(i,k,sphum)*cp_vapor
273 cvm(i) = (1.-q0(i,k,sphum))*cv_air + q0(i,k,sphum)*
cv_vap 276 elseif ( nwat==2 )
then 279 c_air = cp_air*
vicpqd(q0(i,k,:))
280 cpm(i) = c_air + (cp_vapor-c_air)*q0(i,k,sphum)/(1.0-q0(i,k,nwat))
281 c_air = cv_air*
vicvqd(q0(i,k,:))
282 cvm(i) = c_air + (
cv_vap-c_air)*q0(i,k,sphum)/(1.0-q0(i,k,nwat))
284 cpm(i) = (1.-q0(i,k,sphum))*cp_air + q0(i,k,sphum)*cp_vapor
285 cvm(i) = (1.-q0(i,k,sphum))*cv_air + q0(i,k,sphum)*
cv_vap 288 elseif ( nwat==3 )
then 290 q_liq = q0(i,k,liq_wat)
291 q_sol = q0(i,k,ice_wat)
293 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 294 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 296 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 297 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 300 elseif ( nwat==4 )
then 303 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
305 cpm(i) = (1.-(q0(i,k,sphum)+q_liq))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq 306 cvm(i) = (1.-(q0(i,k,sphum)+q_liq))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq 308 cpm(i) = (1.-(q0(i,k,sphum)+q_liq))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq 309 cvm(i) = (1.-(q0(i,k,sphum)+q_liq))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq 313 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
314 q_sol = q0(i,k,ice_wat)
316 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 317 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 319 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 320 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 326 elseif ( nwat==5 )
then 328 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
329 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat)
331 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 332 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 334 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 335 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 340 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
341 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat) + q0(i,k,graupel)
343 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 344 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 346 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 347 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 353 den(i,k) = -delp(i,j,k)/(grav*delz(i,j,k))
355 gz(i,k) = gzh(i) - g2*delz(i,j,k)
356 tmp = gz(i,k) + 0.5*(u0(i,k)**2+v0(i,k)**2+w0(i,k)**2)
357 hd(i,k) = cpm(i)*t0(i,k) + tmp
358 te(i,k) = cvm(i)*t0(i,k) + tmp
359 gzh(i) = gzh(i) - grav*delz(i,j,k)
367 if ( n==1) ratio = 0.25
368 if ( n==2) ratio = 0.5
369 if ( n==3) ratio = 0.999
371 ratio =
real(n)/
real(m)
385 elseif ( nwat==2 )
then 388 qcon(i,k) = q0(i,k,liq_wat)
391 elseif ( nwat==3 )
then 394 qcon(i,k) = q0(i,k,liq_wat) + q0(i,k,ice_wat)
397 elseif ( nwat==4 )
then 401 qcon(i,k) = q0(i,k,liq_wat) + q0(i,k,rainwat)
403 qcon(i,k) = q0(i,k,liq_wat) + q0(i,k,rainwat) + q0(i,k,ice_wat)
407 elseif ( nwat==5 )
then 410 qcon(i,k) = q0(i,k,liq_wat)+q0(i,k,ice_wat)+q0(i,k,snowwat)+q0(i,k,rainwat)
416 qcon(i,k) = q0(i,k,liq_wat)+q0(i,k,ice_wat)+q0(i,k,snowwat)+q0(i,k,rainwat)+q0(i,k,graupel)
427 tv1 = t0(i,km1)*
virq_qpz(q0(i,km1,:),qcon(i,km1))
428 tv2 = t0(i,k )*
virq_qpz(q0(i,k, :),qcon(i,k ))
430 tv1 = t0(i,km1)*(1.+xvir*q0(i,km1,sphum)-qcon(i,km1))
431 tv2 = t0(i,k )*(1.+xvir*q0(i,k ,sphum)-qcon(i,k))
433 pt1 = tv1 / pkz(i,j,km1)
434 pt2 = tv2 / pkz(i,j,k )
436 ri = (gz(i,km1)-gz(i,k))*(pt1-pt2)/( 0.5*(pt1+pt2)* &
437 ((u0(i,km1)-u0(i,k))**2+(v0(i,km1)-v0(i,k))**2+ustar2) )
438 if ( tv1>t_max .and. tv1>tv2 )
then 441 elseif ( tv2<t_min )
then 458 if ( ri < ri_ref )
then 459 mc = ratio*delp(i,j,km1)*delp(i,j,k)/(delp(i,j,km1)+delp(i,j,k))*(1.-max(0.0,ri/ri_ref))**2
461 h0 = mc*(q0(i,k,iq)-q0(i,km1,iq))
462 q0(i,km1,iq) = q0(i,km1,iq) + h0/delp(i,j,km1)
463 q0(i,k ,iq) = q0(i,k ,iq) - h0/delp(i,j,k )
468 elseif ( nwat==2 )
then 469 qcon(i,km1) = q0(i,km1,liq_wat)
470 elseif ( nwat==3 )
then 471 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat)
472 elseif ( nwat==4 )
then 474 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,rainwat)
476 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat) + &
479 elseif ( nwat==5 )
then 480 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat) + &
481 q0(i,km1,snowwat) + q0(i,km1,rainwat)
483 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat) + &
484 q0(i,km1,snowwat) + q0(i,km1,rainwat) + q0(i,km1,graupel)
487 h0 = mc*(u0(i,k)-u0(i,k-1))
488 u0(i,k-1) = u0(i,k-1) + h0/delp(i,j,k-1)
489 u0(i,k ) = u0(i,k ) - h0/delp(i,j,k )
491 h0 = mc*(v0(i,k)-v0(i,k-1))
492 v0(i,k-1) = v0(i,k-1) + h0/delp(i,j,k-1)
493 v0(i,k ) = v0(i,k ) - h0/delp(i,j,k )
495 if ( hydrostatic )
then 497 h0 = mc*(hd(i,k)-hd(i,k-1))
498 hd(i,k-1) = hd(i,k-1) + h0/delp(i,j,k-1)
499 hd(i,k ) = hd(i,k ) - h0/delp(i,j,k )
502 h0 = mc*(hd(i,k)-hd(i,k-1))
503 te(i,k-1) = te(i,k-1) + h0/delp(i,j,k-1)
504 te(i,k ) = te(i,k ) - h0/delp(i,j,k )
506 h0 = mc*(w0(i,k)-w0(i,k-1))
507 w0(i,k-1) = w0(i,k-1) + h0/delp(i,j,k-1)
508 w0(i,k ) = w0(i,k ) - h0/delp(i,j,k )
516 if ( hydrostatic )
then 520 rkx = cp_air/rdgas * (
vicpqd(q0(i,kk,:))/
virqd(q0(i,kk,:))) + 1
521 rdx = rdgas*
virqd(q0(i,kk,:))
523 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
524 / ( rkx - pe(i,kk,j)/pm(i,kk) )
525 gzh(i) = gzh(i) + t0(i,kk)*(peln(i,kk+1,j)-peln(i,kk,j))
526 rdx = rdx + rzx*q0(i,kk,sphum)/(1.0-sum(q0(i,kk,sphum+1:sphum+nwat-1)))
527 t0(i,kk) = t0(i,kk) / rdx
529 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
530 / ( rk - pe(i,kk,j)/pm(i,kk) )
531 gzh(i) = gzh(i) + t0(i,kk)*(peln(i,kk+1,j)-peln(i,kk,j))
532 t0(i,kk) = t0(i,kk) / ( rdgas + rz*q0(i,kk,sphum) )
538 rkx = cp_air/rdgas * (
vicpqd(q0(i,kk,:))/
virqd(q0(i,kk,:))) + 1
539 rdx = rdgas*
virqd(q0(i,kk,:))
541 rdx = rdx + rzx*q0(i,kk,sphum)/(1.0-sum(q0(i,kk,sphum+1:sphum+nwat-1)))
542 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
543 / ((rkx-pe(i,kk,j)/pm(i,kk))*rdx)
545 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
546 / ((rk-pe(i,kk,j)/pm(i,kk))*(rdgas+rz*q0(i,kk,sphum)))
552 if ( nwat == 0 )
then 555 cpm(i) = cp_air*
vicpqd(q0(i,kk,:))
556 cvm(i) = cv_air*
vicvqd(q0(i,kk,:))
562 elseif ( nwat == 1 )
then 565 cpm(i) = (1.-q0(i,kk,sphum))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor
566 cvm(i) = (1.-q0(i,kk,sphum))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap 568 cpm(i) = (1.-q0(i,kk,sphum))*cp_air + q0(i,kk,sphum)*cp_vapor
569 cvm(i) = (1.-q0(i,kk,sphum))*cv_air + q0(i,kk,sphum)*
cv_vap 572 elseif ( nwat == 2 )
then 575 q_liq = q0(i,kk,liq_wat)
576 c_air = cp_air*
vicpqd(q0(i,kk,:))
577 cpm(i) = c_air + (cp_vapor-c_air)*q0(i,kk,sphum)/(1.0-q_liq)
578 c_air = cv_air*
vicvqd(q0(i,kk,:))
579 cvm(i) = c_air + (
cv_vap-c_air)*q0(i,kk,sphum)/(1.0-q_liq)
581 cpm(i) = (1.-q0(i,kk,sphum))*cp_air + q0(i,kk,sphum)*cp_vapor
582 cvm(i) = (1.-q0(i,kk,sphum))*cv_air + q0(i,kk,sphum)*
cv_vap 585 elseif ( nwat == 3 )
then 587 q_liq = q0(i,kk,liq_wat)
588 q_sol = q0(i,kk,ice_wat)
590 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 591 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 593 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 594 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 597 elseif ( nwat == 4 )
then 600 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
602 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq 603 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq 605 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq 606 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq 609 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
610 q_sol = q0(i,kk,ice_wat)
612 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 613 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 615 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 616 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 622 elseif ( nwat == 5 )
then 624 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
625 q_sol = q0(i,kk,ice_wat) + q0(i,kk,snowwat)
627 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 628 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 630 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 631 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 636 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
637 q_sol = q0(i,kk,ice_wat) + q0(i,kk,snowwat) + q0(i,kk,graupel)
639 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 640 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 642 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 643 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 649 tv = gz(i,kk) + 0.5*(u0(i,kk)**2+v0(i,kk)**2+w0(i,kk)**2)
650 t0(i,kk) = (te(i,kk)- tv) / cvm(i)
651 hd(i,kk) = cpm(i)*t0(i,kk) + tv
662 t0(i,k) = ta(i,j,k) + (t0(i,k) - ta(i,j,k))*fra
663 u0(i,k) = ua(i,j,k) + (u0(i,k) - ua(i,j,k))*fra
664 v0(i,k) = va(i,j,k) + (v0(i,k) - va(i,j,k))*fra
668 if ( .not. hydrostatic )
then 671 w0(i,k) = w(i,j,k) + (w0(i,k) - w(i,j,k))*fra
679 q0(i,k,iq) = qa(i,j,k,iq) + (q0(i,k,iq) - qa(i,j,k,iq))*fra
687 u_dt(i,j,k) = rdt*(u0(i,k) - ua(i,j,k))
688 v_dt(i,j,k) = rdt*(v0(i,k) - va(i,j,k))
697 qa(i,j,k,iq) = q0(i,k,iq)
702 if ( .not. hydrostatic )
then 715 subroutine fv_subgrid_z( isd, ied, jsd, jed, is, ie, js, je, km, nq, dt, &
716 tau, nwat, delp, pe, peln, pkz, ta, qa, ua, va, &
717 hydrostatic, w, delz, u_dt, v_dt, t_dt, q_dt, k_bot )
720 integer,
intent(in):: is, ie, js, je, km, nq, nwat
721 integer,
intent(in):: isd, ied, jsd, jed
722 integer,
intent(in):: tau
723 real,
intent(in):: dt
724 real,
intent(in):: pe(is-1:ie+1,km+1,js-1:je+1)
725 real,
intent(in):: peln(is :ie, km+1,js :je)
726 real,
intent(in):: delp(isd:ied,jsd:jed,km)
727 real,
intent(in):: delz(is:,js:,1:)
728 real,
intent(in):: pkz(is:ie,js:je,km)
729 logical,
intent(in):: hydrostatic
730 integer,
intent(in),
optional:: k_bot
732 real,
intent(inout):: ua(isd:ied,jsd:jed,km)
733 real,
intent(inout):: va(isd:ied,jsd:jed,km)
734 real,
intent(inout):: w(isd:,jsd:,1:)
735 real,
intent(inout):: ta(isd:ied,jsd:jed,km)
736 real,
intent(inout):: qa(isd:ied,jsd:jed,km,nq)
737 real,
intent(inout):: u_dt(isd:ied,jsd:jed,km)
738 real,
intent(inout):: v_dt(isd:ied,jsd:jed,km)
739 real,
intent(inout):: t_dt(is:ie,js:je,km)
740 real,
intent(inout):: q_dt(is:ie,js:je,km,nq)
742 real,
dimension(is:ie,km):: u0, v0, w0, t0, hd, te, gz, tvm, pm, den
743 real q0(is:ie,km,nq), qcon(is:ie,km)
744 real,
dimension(is:ie):: gzh, lcp2, icp2, cvm, cpm, qs
746 real :: rkx, rdx, rzx, c_air
748 real ri_ref, ri, pt1, pt2, ratio, tv, cv, tmp, q_liq, q_sol
749 real tv1, tv2, g2, h0, mc, fra, rk, rz, rdt, tvd, tv_surf
750 real dh, dq, qsw, dqsdt, tcp3
751 integer i, j, k, kk, n, m, iq, km1, im, kbot
752 real,
parameter:: ustar2 = 1.e-4
754 integer :: sphum, liq_wat, rainwat, snowwat, graupel, ice_wat, cld_amt
756 cv_air = cp_air - rdgas
757 rk = cp_air/rdgas + 1.
765 if (
present(k_bot) )
then 766 if ( k_bot < 3 )
return 772 sphum = get_tracer_index(model_atmos,
'sphum')
773 if ( nwat == 0 )
then 779 liq_wat = get_tracer_index(model_atmos,
'liq_wat')
780 ice_wat = get_tracer_index(model_atmos,
'ice_wat')
781 rainwat = get_tracer_index(model_atmos,
'rainwat')
782 snowwat = get_tracer_index(model_atmos,
'snowwat')
783 graupel = get_tracer_index(model_atmos,
'graupel')
784 cld_amt = get_tracer_index(model_atmos,
'cld_amt')
809 q0(i,k,iq) = qa(i,j,k,iq)
818 tvm(i,k) = t0(i,k)*
virq(q0(i,k,:))
820 tvm(i,k) = t0(i,k)*(1.+xvir*q0(i,k,sphum))
824 pm(i,k) = delp(i,j,k)/(peln(i,k+1,j)-peln(i,k,j))
832 if( hydrostatic )
then 836 den(i,k) = pm(i,k)/tv
837 gz(i,k) = gzh(i) + tv*(1.-pe(i,k,j)/pm(i,k))
838 hd(i,k) = cp_air*tvm(i,k)+gz(i,k)+0.5*(u0(i,k)**2+v0(i,k)**2)
839 gzh(i) = gzh(i) + tv*(peln(i,k+1,j)-peln(i,k,j))
844 if ( nwat == 0 )
then 847 cpm(i) = cp_air*
vicpqd(q0(i,k,:))
848 cvm(i) = cv_air*
vicvqd(q0(i,k,:))
854 elseif ( nwat==1 )
then 857 cpm(i) = (1.-q0(i,k,sphum))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor
858 cvm(i) = (1.-q0(i,k,sphum))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap 860 cpm(i) = (1.-q0(i,k,sphum))*cp_air + q0(i,k,sphum)*cp_vapor
861 cvm(i) = (1.-q0(i,k,sphum))*cv_air + q0(i,k,sphum)*
cv_vap 864 elseif ( nwat==2 )
then 866 q_sol = q0(i,k,liq_wat)
868 c_air = cp_air*
vicpqd(q0(i,k,:))
869 cpm(i) = c_air + (cp_vapor-c_air)*q0(i,k,sphum)/(1.0-q0(i,k,nwat))
870 c_air = cv_air*
vicvqd(q0(i,k,:))
871 cvm(i) = c_air + (
cv_vap-c_air)*q0(i,k,sphum)/(1.0-q0(i,k,nwat))
873 cpm(i) = (1.-q0(i,k,sphum))*cp_air + q0(i,k,sphum)*cp_vapor
874 cvm(i) = (1.-q0(i,k,sphum))*cv_air + q0(i,k,sphum)*
cv_vap 877 elseif ( nwat==3 )
then 879 q_liq = q0(i,k,liq_wat)
880 q_sol = q0(i,k,ice_wat)
882 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 883 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 885 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 886 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 889 elseif ( nwat==4 )
then 892 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
894 cpm(i) = (1.-(q0(i,k,sphum)+q_liq))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq 895 cvm(i) = (1.-(q0(i,k,sphum)+q_liq))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq 897 cpm(i) = (1.-(q0(i,k,sphum)+q_liq))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq 898 cvm(i) = (1.-(q0(i,k,sphum)+q_liq))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq 901 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
902 q_sol = q0(i,k,ice_wat)
904 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 905 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 907 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 908 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 913 elseif ( nwat==5 )
then 915 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
916 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat)
918 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 919 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 921 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 922 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 927 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
928 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat) + q0(i,k,graupel)
930 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 931 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 933 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 934 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 940 den(i,k) = -delp(i,j,k)/(grav*delz(i,j,k))
942 gz(i,k) = gzh(i) - g2*delz(i,j,k)
943 tmp = gz(i,k) + 0.5*(u0(i,k)**2+v0(i,k)**2+w0(i,k)**2)
944 hd(i,k) = cpm(i)*t0(i,k) + tmp
945 te(i,k) = cvm(i)*t0(i,k) + tmp
946 gzh(i) = gzh(i) - grav*delz(i,j,k)
954 if ( n==1) ratio = 0.25
955 if ( n==2) ratio = 0.5
956 if ( n==3) ratio = 0.999
958 ratio =
real(n)/
real(m)
972 elseif ( nwat==2 )
then 975 qcon(i,k) = q0(i,k,liq_wat)
978 elseif ( nwat==3 )
then 981 qcon(i,k) = q0(i,k,liq_wat) + q0(i,k,ice_wat)
984 elseif ( nwat==4 )
then 988 qcon(i,k) = q0(i,k,liq_wat) + q0(i,k,rainwat)
990 qcon(i,k) = q0(i,k,liq_wat) + q0(i,k,rainwat) + q0(i,k,ice_wat)
994 elseif ( nwat==5 )
then 997 qcon(i,k) = q0(i,k,liq_wat)+q0(i,k,ice_wat)+q0(i,k,snowwat)+q0(i,k,rainwat)
1002 qcon(i,k) = q0(i,k,liq_wat)+q0(i,k,ice_wat)+q0(i,k,snowwat)+q0(i,k,rainwat)+q0(i,k,graupel)
1011 if(nwat>0)
call qsmith(im, km, im, 1, t0(is,km1), pm(is,km1), q0(is,km1,sphum), qs)
1013 if(nwat>0)
call qsmith(im, 1, 1, t0(is,km1), pm(is,km1), q0(is,km1,sphum), qs)
1020 tv1 = t0(i,km1)*
virq_qpz(q0(i,km1,:),qcon(i,km1))
1021 tv2 = t0(i,k )*
virq_qpz(q0(i,k ,:),qcon(i,k ))
1023 tv1 = t0(i,km1)*(1.+xvir*q0(i,km1,sphum)-qcon(i,km1))
1024 tv2 = t0(i,k )*(1.+xvir*q0(i,k ,sphum)-qcon(i,k))
1026 pt1 = tv1 / pkz(i,j,km1)
1027 pt2 = tv2 / pkz(i,j,k )
1028 ri = (gz(i,km1)-gz(i,k))*(pt1-pt2)/( 0.5*(pt1+pt2)* &
1029 ((u0(i,km1)-u0(i,k))**2+(v0(i,km1)-v0(i,k))**2+ustar2) )
1036 if ( nwat > 0 )
then 1037 h0 = hd(i,k)-hd(i,km1) +
hlv0*(q0(i,k,sphum)-qs(i))
1038 if ( h0 > 0. ) ri_ref = 5.*ri_ref
1041 if ( ri < ri_ref )
then 1042 mc = ratio*delp(i,j,km1)*delp(i,j,k)/(delp(i,j,km1)+delp(i,j,k))*(1.-max(0.0,ri/ri_ref))**2
1044 h0 = mc*(q0(i,k,iq)-q0(i,km1,iq))
1045 q0(i,km1,iq) = q0(i,km1,iq) + h0/delp(i,j,km1)
1046 q0(i,k ,iq) = q0(i,k ,iq) - h0/delp(i,j,k )
1051 elseif ( nwat==2 )
then 1052 qcon(i,km1) = q0(i,km1,liq_wat)
1053 elseif ( nwat==3 )
then 1054 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat)
1055 elseif ( nwat==4 )
then 1057 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,rainwat)
1059 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,rainwat) + q0(i,km1,ice_wat)
1061 elseif ( nwat==5 )
then 1062 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat) + &
1063 q0(i,km1,snowwat) + q0(i,km1,rainwat)
1065 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat) + &
1066 q0(i,km1,snowwat) + q0(i,km1,rainwat) + q0(i,km1,graupel)
1069 h0 = mc*(u0(i,k)-u0(i,k-1))
1070 u0(i,k-1) = u0(i,k-1) + h0/delp(i,j,k-1)
1071 u0(i,k ) = u0(i,k ) - h0/delp(i,j,k )
1073 h0 = mc*(v0(i,k)-v0(i,k-1))
1074 v0(i,k-1) = v0(i,k-1) + h0/delp(i,j,k-1)
1075 v0(i,k ) = v0(i,k ) - h0/delp(i,j,k )
1077 if ( hydrostatic )
then 1079 h0 = mc*(hd(i,k)-hd(i,k-1))
1080 hd(i,k-1) = hd(i,k-1) + h0/delp(i,j,k-1)
1081 hd(i,k ) = hd(i,k ) - h0/delp(i,j,k )
1084 h0 = mc*(hd(i,k)-hd(i,k-1))
1085 te(i,k-1) = te(i,k-1) + h0/delp(i,j,k-1)
1086 te(i,k ) = te(i,k ) - h0/delp(i,j,k )
1088 h0 = mc*(w0(i,k)-w0(i,k-1))
1089 w0(i,k-1) = w0(i,k-1) + h0/delp(i,j,k-1)
1090 w0(i,k ) = w0(i,k ) - h0/delp(i,j,k )
1098 if ( hydrostatic )
then 1102 rkx = cp_air/rdgas * (
vicpqd(q0(i,kk,:))/
virqd(q0(i,kk,:))) + 1
1103 rdx = rdgas*
virqd(q0(i,kk,:))
1105 rdx = rdx + rzx*q0(i,kk,sphum)/(1.0-sum(q0(i,kk,sphum+1:sphum+nwat-1)))
1106 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
1107 / ( rkx - pe(i,kk,j)/pm(i,kk) )
1108 gzh(i) = gzh(i) + t0(i,kk)*(peln(i,kk+1,j)-peln(i,kk,j))
1109 t0(i,kk) = t0(i,kk) / ( rdx + rzx*q0(i,kk,sphum) )
1111 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
1112 / ( rk - pe(i,kk,j)/pm(i,kk) )
1113 gzh(i) = gzh(i) + t0(i,kk)*(peln(i,kk+1,j)-peln(i,kk,j))
1114 t0(i,kk) = t0(i,kk) / ( rdgas + rz*q0(i,kk,sphum) )
1120 rkx = cp_air/rdgas * (
vicpqd(q0(i,kk,:))/
virqd(q0(i,kk,:))) + 1
1121 rdx = rdgas*
virqd(i)
1123 rdx = rdx + rzx*q0(i,kk,sphum)/(1.0-sum(q0(i,kk,sphum+1:sphum+nwat-1)))
1124 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
1125 / ((rkx-pe(i,kk,j)/pm(i,kk))*rdx)
1127 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
1128 / ((rk-pe(i,kk,j)/pm(i,kk))*(rdgas+rz*q0(i,kk,sphum)))
1134 if ( nwat == 0 )
then 1137 cpm(i) = cp_air*
vicpqd(q0(i,kk,:))
1138 cvm(i) = cv_air*
vicvqd(q0(i,kk,:))
1144 elseif ( nwat == 1 )
then 1147 cpm(i) = (1.-q0(i,kk,sphum))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor
1148 cvm(i) = (1.-q0(i,kk,sphum))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap 1150 cpm(i) = (1.-q0(i,kk,sphum))*cp_air + q0(i,kk,sphum)*cp_vapor
1151 cvm(i) = (1.-q0(i,kk,sphum))*cv_air + q0(i,kk,sphum)*
cv_vap 1154 elseif ( nwat == 2 )
then 1157 c_air = cp_air*
vicpqd(q0(i,kk,:))
1158 cpm(i) = c_air+(cp_vapor-c_air)*q0(i,kk,sphum)/(1.0-q0(i,kk,sphum+nwat-1))
1159 c_air = cv_air*
vicvqd(q0(i,kk,:))
1160 cvm(i) = c_air+(
cv_vap -c_air)*q0(i,kk,sphum)/(1.0-q0(i,kk,sphum+nwat-1))
1162 cpm(i) = (1.-q0(i,kk,sphum))*cp_air + q0(i,kk,sphum)*cp_vapor
1163 cvm(i) = (1.-q0(i,kk,sphum))*cv_air + q0(i,kk,sphum)*
cv_vap 1166 elseif ( nwat == 3 )
then 1168 q_liq = q0(i,kk,liq_wat)
1169 q_sol = q0(i,kk,ice_wat)
1171 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1172 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1174 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1175 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1178 elseif ( nwat == 4 )
then 1181 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
1183 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq 1184 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq 1186 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq 1187 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq 1190 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
1191 q_sol = q0(i,kk,ice_wat)
1193 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1194 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1196 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1197 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1202 elseif ( nwat == 5 )
then 1204 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
1205 q_sol = q0(i,kk,ice_wat) + q0(i,kk,snowwat)
1207 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1208 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1210 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1211 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1216 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
1217 q_sol = q0(i,kk,ice_wat) + q0(i,kk,snowwat) + q0(i,kk,graupel)
1219 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1220 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1222 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1223 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq+q_sol))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1230 tv = gz(i,kk) + 0.5*(u0(i,kk)**2+v0(i,kk)**2+w0(i,kk)**2)
1231 t0(i,kk) = (te(i,kk)- tv) / cvm(i)
1232 hd(i,kk) = cpm(i)*t0(i,kk) + tv
1240 if ( fra < 1. )
then 1243 t0(i,k) = ta(i,j,k) + (t0(i,k) - ta(i,j,k))*fra
1244 u0(i,k) = ua(i,j,k) + (u0(i,k) - ua(i,j,k))*fra
1245 v0(i,k) = va(i,j,k) + (v0(i,k) - va(i,j,k))*fra
1249 if ( .not. hydrostatic )
then 1252 w0(i,k) = w(i,j,k) + (w0(i,k) - w(i,j,k))*fra
1260 q0(i,k,iq) = qa(i,j,k,iq) + (q0(i,k,iq) - qa(i,j,k,iq))*fra
1270 if ( nwat > 5 )
then 1272 if ( hydrostatic )
then 1276 den(i,k) = pm(i,k)/(rdgas*t0(i,k)*
virq(q0(i,k,:)))
1278 den(i,k) = pm(i,k)/(rdgas*t0(i,k)*(1.+xvir*q0(i,k,sphum)))
1280 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
1281 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat) + q0(i,k,graupel)
1283 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1285 cpm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1287 lcp2(i) = hlv / cpm(i)
1288 icp2(i) = hlf / cpm(i)
1292 den(i,k) = -delp(i,j,k)/(grav*delz(i,j,k))
1293 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
1294 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat) + q0(i,k,graupel)
1296 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1298 cvm(i) = (1.-(q0(i,k,sphum)+q_liq+q_sol))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1307 qsw = wqs2(t0(i,k), den(i,k), dqsdt)
1308 dq = q0(i,k,sphum) - qsw
1310 tcp3 = lcp2(i) + icp2(i)*min(1., dim(
tice,t0(i,k))/40.)
1311 tmp = dq/(1.+tcp3*dqsdt)
1312 t0(i,k) = t0(i,k) + tmp*lcp2(i)
1313 q0(i,k, sphum) = q0(i,k, sphum) - tmp
1314 q0(i,k,liq_wat) = q0(i,k,liq_wat) + tmp
1316 if (cld_amt .gt. 0)
then 1317 qa(i,j,k,cld_amt) = max(0.5, min(1., qa(i,j,k,cld_amt)+25.*tmp/qsw))
1321 tmp =
tice-40. - t0(i,k)
1322 if( tmp>0.0 .and. q0(i,k,liq_wat)>0. )
then 1323 dh = min( q0(i,k,liq_wat), q0(i,k,liq_wat)*tmp*0.125, tmp/icp2(i) )
1324 q0(i,k,liq_wat) = q0(i,k,liq_wat) - dh
1325 q0(i,k,ice_wat) = q0(i,k,ice_wat) + dh
1326 t0(i,k) = t0(i,k) + dh*icp2(i)
1335 u_dt(i,j,k) = rdt*(u0(i,k) - ua(i,j,k))
1336 v_dt(i,j,k) = rdt*(v0(i,k) - va(i,j,k))
1344 if (iq .ne. cld_amt )
then 1346 qa(i,j,k,iq) = q0(i,k,iq)
1352 if ( .not. hydrostatic )
then 1368 integer,
parameter:: length=2621
1371 if( .not.
allocated(
table) )
then 1374 allocate (
table(length) )
1375 allocate (
des(length) )
1382 des(length) =
des(length-1)
1389 subroutine qsmith(im, km, imx, kmx, t, p, q, qs, dqdt)
1391 subroutine qsmith(im, km, k1, t, p, q, qs, dqdt)
1394 real,
intent(in),
dimension(im,km):: t, p
1396 real,
intent(in),
dimension(im,km,*):: q
1397 integer,
intent(in):: im, km, imx, kmx
1399 integer,
intent(in):: im, km, k1
1400 real,
intent(in),
dimension(im,km):: q
1402 real,
intent(out),
dimension(im,km):: qs
1403 real,
intent(out),
optional:: dqdt(im,km)
1422 ap1 = 10.*dim(t(i,k), tmin) + 1.
1423 ap1 = min(2621., ap1)
1425 es(i,k) =
table(it) + (ap1-it)*
des(it)
1429 qs(i,k) =
esl*es(i,k)*(1.+
zvir*q(i,k))/p(i,k)
1434 if (
present(dqdt) )
then 1442 ap1 = 10.*dim(t(i,k), tmin) + 1.
1443 ap1 = min(2621., ap1) - 0.5
1448 dqdt(i,k) = eps10*(
des(it)+(ap1-it)*(
des(it+1)-
des(it)))*(1.+
zvir*q(i,k))/p(i,k)
1458 integer,
intent(in):: n
1461 real esbasw, tbasw, esbasi, tbasi, Tmin, tem, aa, b, c, d, e, esh20
1472 tem = tmin+dt*
real(i-1)
1473 aa = -7.90298*(tbasw/tem-1)
1474 b = 5.02808*alog10(tbasw/tem)
1475 c = -1.3816e-07*(10**((1-tem/tbasw)*11.344)-1)
1476 d = 8.1328e-03*(10**((tbasw/tem-1)*(-3.49149))-1)
1478 table(i) = 0.1*10**(aa+b+c+d+e)
1485 integer,
intent(in):: n
1489 real esbasw, tbasw, esbasi, tbasi, Tmin, tem, aa, b, c, d, e, esh20
1503 tem = tmin+dt*
real(i-1)
1504 aa = -9.09718 *(tbasi/tem-1.0)
1505 b = -3.56654 *alog10(tbasi/tem)
1506 c = 0.876793*(1.0-tem/tbasi)
1508 table(i)=10**(aa+b+c+e)
1514 tem = 253.16+dt*
real(i-1)
1515 aa = -7.90298*(tbasw/tem-1)
1516 b = 5.02808*alog10(tbasw/tem)
1517 c = -1.3816e-07*(10**((1-tem/tbasw)*11.344)-1)
1518 d = 8.1328e-03*(10**((tbasw/tem-1)*(-3.49149))-1)
1520 esh20 = 10**(aa+b+c+d+e)
1524 table(i+1400) = esh20
1530 tem = 253.16+dt*
real(i-1)
1531 wice = 0.05*(273.16-tem)
1532 wh2o = 0.05*(tem-253.16)
1533 table(i+1400) = wice*table(i+1400)+wh2o*esupc(i)
1537 table(i) = table(i)*0.1
1542 subroutine neg_adj3(is, ie, js, je, ng, kbot, hydrostatic, peln, delz, pt, dp, &
1548 ql, qr, qi, qs, qg, qa, check_negative)
1551 integer,
intent(in):: is, ie, js, je, ng, kbot
1552 logical,
intent(in):: hydrostatic
1553 real,
intent(in):: dp(is-ng:ie+ng,js-ng:je+ng,kbot)
1554 real,
intent(in):: delz(is:,js:,1:)
1555 real,
intent(in):: peln(is:ie,kbot+1,js:je)
1556 logical,
intent(in),
OPTIONAL :: check_negative
1558 real,
intent(inout),
dimension(is-ng:ie+ng,js-ng:je+ng,kbot,*):: qvi
1560 real,
intent(inout),
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: qv
1562 real,
intent(inout),
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: &
1563 pt, ql, qr, qi, qs, qg
1564 real,
intent(inout),
OPTIONAL,
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: qa
1566 logical:: sat_adj = .false.
1567 real,
parameter :: t48 =
tice - 48.
1569 real,
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: qv
1571 real,
dimension(is:ie,kbot):: dpk, q2
1572 real,
dimension(is:ie,js:je):: pt2, qv2, ql2, qi2, qs2, qr2, qg2, dp2, p2, icpk, lcpk
1574 real:: dq, qsum, dq1, q_liq, q_sol, cpm, sink, qsw, dwsdt
1577 cv_air = cp_air - rdgas
1580 qv(:,:,:) = qvi(:,:,:,1)
1583 if (
present(check_negative) )
then 1584 if ( check_negative )
then 1585 call prt_negative(
'Temperature', pt, is, ie, js, je, ng, kbot, 165.)
1586 call prt_negative(
'sphum', qv, is, ie, js, je, ng, kbot, -1.e-8)
1587 call prt_negative(
'liq_wat', ql, is, ie, js, je, ng, kbot, -1.e-7)
1588 call prt_negative(
'rainwat', qr, is, ie, js, je, ng, kbot, -1.e-7)
1589 call prt_negative(
'ice_wat', qi, is, ie, js, je, ng, kbot, -1.e-7)
1590 call prt_negative(
'snowwat', qs, is, ie, js, je, ng, kbot, -1.e-7)
1591 call prt_negative(
'graupel', qg, is, ie, js, je, ng, kbot, -1.e-7)
1593 qvi(:,:,:,1) = qv(:,:,:)
1598 if ( hydrostatic )
then 1615 qv2(i,j) = qv(i,j,k)
1616 ql2(i,j) = ql(i,j,k)
1617 qi2(i,j) = qi(i,j,k)
1618 qs2(i,j) = qs(i,j,k)
1619 qr2(i,j) = qr(i,j,k)
1620 qg2(i,j) = qg(i,j,k)
1621 dp2(i,j) = dp(i,j,k)
1622 pt2(i,j) = pt(i,j,k)
1626 if ( hydrostatic )
then 1629 p2(i,j) = dp2(i,j)/(peln(i,k+1,j)-peln(i,k,j))
1630 q_liq = max(0., ql2(i,j) + qr2(i,j))
1631 q_sol = max(0., qi2(i,j) + qs2(i,j))
1633 cpm = (1.-(qv2(i,j)+q_liq+q_sol))*cp_air*
vicpqd_qpz(qvi(i,j,k,1:
num_gas),qv2(i,j)+q_liq+q_sol) + &
1634 qv2(i,j)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1636 cpm = (1.-(qv2(i,j)+q_liq+q_sol))*cp_air + qv2(i,j)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1638 lcpk(i,j) = hlv / cpm
1639 icpk(i,j) = hlf / cpm
1645 q_liq = max(0., ql2(i,j) + qr2(i,j))
1646 q_sol = max(0., qi2(i,j) + qs2(i,j))
1648 p2(i,j) = -dp2(i,j)/(grav*delz(i,j,k))*rdgas*pt2(i,j)*
virq(qvi(i,j,k,1:
num_gas))
1649 cpm = (1.-(qv2(i,j)+q_liq+q_sol))*cv_air*
vicvqd_qpz(qvi(i,j,k,1:
num_gas),qv2(i,j)+q_liq+q_sol) + &
1652 p2(i,j) = -dp2(i,j)/(grav*delz(i,j,k))*rdgas*pt2(i,j)*(1.+
zvir*qv2(i,j))
1653 cpm = (1.-(qv2(i,j)+q_liq+q_sol))*cv_air + qv2(i,j)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1667 qsum = qi2(i,j) + qs2(i,j)
1668 if ( qsum > 0. )
then 1669 if ( qi2(i,j) < 0. )
then 1672 elseif ( qs2(i,j) < 0. )
then 1680 qg2(i,j) = qg2(i,j) + qsum
1685 if ( qg2(i,j) < 0. )
then 1686 dq = min( qs2(i,j), -qg2(i,j) )
1687 qs2(i,j) = qs2(i,j) - dq
1688 qg2(i,j) = qg2(i,j) + dq
1689 if ( qg2(i,j) < 0. )
then 1691 dq = min( qi2(i,j), -qg2(i,j) )
1692 qi2(i,j) = qi2(i,j) - dq
1693 qg2(i,j) = qg2(i,j) + dq
1698 if ( qg2(i,j)<0. .and. qr2(i,j)>0. )
then 1699 dq = min( qr2(i,j), -qg2(i,j) )
1700 qg2(i,j) = qg2(i,j) + dq
1701 qr2(i,j) = qr2(i,j) - dq
1702 pt2(i,j) = pt2(i,j) + dq*icpk(i,j)
1705 if ( qg2(i,j)<0. .and. ql2(i,j)>0. )
then 1706 dq = min( ql2(i,j), -qg2(i,j) )
1707 qg2(i,j) = qg2(i,j) + dq
1708 ql2(i,j) = ql2(i,j) - dq
1709 pt2(i,j) = pt2(i,j) + dq*icpk(i,j)
1712 if ( qg2(i,j)<0. .and. qv2(i,j)>0. )
then 1713 dq = min( 0.999*qv2(i,j), -qg2(i,j) )
1714 qg2(i,j) = qg2(i,j) + dq
1715 qv2(i,j) = qv2(i,j) - dq
1716 pt2(i,j) = pt2(i,j) + dq*(icpk(i,j)+lcpk(i,j))
1722 qsum = ql2(i,j) + qr2(i,j)
1723 if ( qsum > 0. )
then 1724 if ( qr2(i,j) < 0. )
then 1727 elseif ( ql2(i,j) < 0. )
then 1735 dq = min( max(0.0, qg2(i,j)), -qr2(i,j) )
1736 qr2(i,j) = qr2(i,j) + dq
1737 qg2(i,j) = qg2(i,j) - dq
1738 pt2(i,j) = pt2(i,j) - dq*icpk(i,j)
1739 if ( qr(i,j,k) < 0. )
then 1741 dq = min( qi2(i,j)+qs2(i,j), -qr2(i,j) )
1742 qr2(i,j) = qr2(i,j) + dq
1743 dq1 = min( dq, qs2(i,j) )
1744 qs2(i,j) = qs2(i,j) - dq1
1745 qi2(i,j) = qi2(i,j) + dq1 - dq
1746 pt2(i,j) = pt2(i,j) - dq*icpk(i,j)
1749 if ( qr2(i,j)<0. .and. qv2(i,j)>0. )
then 1750 dq = min( 0.999*qv2(i,j), -qr2(i,j) )
1751 qv2(i,j) = qv2(i,j) - dq
1752 qr2(i,j) = qr2(i,j) + dq
1753 pt2(i,j) = pt2(i,j) + dq*lcpk(i,j)
1768 if ( qi2(i,j)>1.e-8 .and. pt2(i,j) >
tice )
then 1769 sink = min( qi2(i,j), (pt2(i,j)-
tice)/icpk(i,j) )
1770 ql2(i,j) = ql2(i,j) + sink
1771 qi2(i,j) = qi2(i,j) - sink
1772 pt2(i,j) = pt2(i,j) - sink*icpk(i,j)
1776 qsw = wqsat2_moist(pt2(i,j), qv2(i,j), p2(i,j), dwsdt)
1777 sink = min( ql2(i,j), (qsw-qv2(i,j))/(1.+lcpk(i,j)*dwsdt) )
1778 qv2(i,j) = qv2(i,j) + sink
1779 ql2(i,j) = ql2(i,j) - sink
1780 pt2(i,j) = pt2(i,j) - sink*lcpk(i,j)
1784 if( ql2(i,j)>1.e-8 .and. pt2(i,j) < t48 )
then 1786 sink = min( ql2(i,j), (t48-pt2(i,j))/icpk(i,j) )
1787 ql2(i,j) = ql2(i,j) - sink
1788 qi2(i,j) = qi2(i,j) + sink
1789 pt2(i,j) = pt2(i,j) + sink*icpk(i,j)
1800 qv(i,j,k) = qv2(i,j)
1801 ql(i,j,k) = ql2(i,j)
1802 qi(i,j,k) = qi2(i,j)
1803 qs(i,j,k) = qs2(i,j)
1804 qr(i,j,k) = qr2(i,j)
1805 qg(i,j,k) = qg2(i,j)
1806 pt(i,j,k) = pt2(i,j)
1810 qvi(:,:,k,1) = qv(:,:,k)
1821 dpk(i,k) = dp(i,j,k)
1825 call fillq(ie-is+1, kbot, q2, dpk)
1837 call fillq(ie-is+1, kbot, q2, dpk)
1853 if( qv(i,j,k) < 0. )
then 1854 qv(i,j,k+1) = qv(i,j,k+1) + qv(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
1866 if( qv(i,j,k) < 0. .and. qv(i,j,k-1) > 0. )
then 1867 dq = min(-qv(i,j,k)*dp(i,j,k), qv(i,j,k-1)*dp(i,j,k-1))
1868 qv(i,j,k-1) = qv(i,j,k-1) - dq/dp(i,j,k-1)
1869 qv(i,j,k ) = qv(i,j,k ) + dq/dp(i,j,k )
1871 if( qv(i,j,k) < 0. )
then 1872 qv(i,j,k+1) = qv(i,j,k+1) + qv(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
1883 if( qv(i,j,kbot) < 0. )
then 1885 if ( qv(i,j,kbot)>=0. )
goto 123
1886 if ( qv(i,j,k) > 0. )
then 1887 dq = min(-qv(i,j,kbot)*dp(i,j,kbot), qv(i,j,k)*dp(i,j,k))
1888 qv(i,j,k ) = qv(i,j,k ) - dq/dp(i,j,k)
1889 qv(i,j,kbot) = qv(i,j,kbot) + dq/dp(i,j,kbot)
1897 qvi(:,:,:,1) = qv(:,:,:)
1901 if (
present(qa))
then 1910 if( qa(i,j,k) < 0. )
then 1911 qa(i,j,k+1) = qa(i,j,k+1) + qa(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
1923 if( qa(i,j,kbot) < 0. .and. qa(i,j,kbot-1)>0.)
then 1924 dq = min(-qa(i,j,kbot)*dp(i,j,kbot), qa(i,j,kbot-1)*dp(i,j,kbot-1))
1925 qa(i,j,kbot-1) = qa(i,j,kbot-1) - dq/dp(i,j,kbot-1)
1926 qa(i,j,kbot ) = qa(i,j,kbot ) + dq/dp(i,j,kbot )
1929 qa(i,j,kbot) = max(0., qa(i,j,kbot))
1937 subroutine neg_adj2(is, ie, js, je, ng, kbot, hydrostatic, &
1938 peln, delz, pt, dp, qv, ql, qr, qi, qs, qa, check_negative)
1941 integer,
intent(in):: is, ie, js, je, ng, kbot
1942 logical,
intent(in):: hydrostatic
1943 real,
intent(in):: dp(is-ng:ie+ng,js-ng:je+ng,kbot)
1944 real,
intent(in):: delz(is-ng:,js-ng:,1:)
1945 real,
intent(in):: peln(is:ie,kbot+1,js:je)
1946 logical,
intent(in),
OPTIONAL :: check_negative
1947 real,
intent(inout),
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: &
1948 pt, qv, ql, qr, qi, qs
1949 real,
intent(inout),
OPTIONAL,
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: qa
1951 logical:: sat_adj = .false.
1952 real,
parameter :: t48 =
tice - 48.
1953 real,
dimension(is:ie,kbot):: dpk, q2
1954 real,
dimension(is:ie,js:je):: pt2, qv2, ql2, qi2, qs2, qr2, dp2, p2, icpk, lcpk
1956 real:: dq, qsum, dq1, q_liq, q_sol, oneocpm, sink, qsw, dwsdt, tx1
1959 cv_air = cp_air - rdgas
1961 if (
present(check_negative) )
then 1962 if ( check_negative )
then 1963 call prt_negative(
'Temperature', pt, is, ie, js, je, ng, kbot, 165.)
1964 call prt_negative(
'sphum', qv, is, ie, js, je, ng, kbot, -1.e-8)
1965 call prt_negative(
'liq_wat', ql, is, ie, js, je, ng, kbot, -1.e-7)
1966 call prt_negative(
'rainwat', qr, is, ie, js, je, ng, kbot, -1.e-7)
1967 call prt_negative(
'ice_wat', qi, is, ie, js, je, ng, kbot, -1.e-7)
1968 call prt_negative(
'snowwat', qs, is, ie, js, je, ng, kbot, -1.e-7)
1972 if ( hydrostatic )
then 1986 qv2(i,j) = qv(i,j,k)
1987 ql2(i,j) = ql(i,j,k)
1988 qi2(i,j) = qi(i,j,k)
1989 qs2(i,j) = qs(i,j,k)
1990 qr2(i,j) = qr(i,j,k)
1991 dp2(i,j) = dp(i,j,k)
1992 pt2(i,j) = pt(i,j,k)
1996 if ( hydrostatic )
then 1999 p2(i,j) = dp2(i,j)/(peln(i,k+1,j)-peln(i,k,j))
2000 q_liq = max(0., ql2(i,j) + qr2(i,j))
2001 q_sol = max(0., qi2(i,j) + qs2(i,j))
2002 oneocpm = 1.0 / ((1.-(qv2(i,j)+q_liq+q_sol))*cp_air + qv2(i,j)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice)
2003 lcpk(i,j) = hlv * oneocpm
2004 icpk(i,j) = hlf * oneocpm
2010 p2(i,j) = -dp2(i,j)/(grav*delz(i,j,k))*rdgas*pt2(i,j)*(1.+
zvir*qv2(i,j))
2011 q_liq = max(0., ql2(i,j) + qr2(i,j))
2012 q_sol = max(0., qi2(i,j) + qs2(i,j))
2013 oneocpm = 1.0 / ((1.-(qv2(i,j)+q_liq+q_sol))*cv_air + qv2(i,j)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice)
2015 icpk(i,j) = (
li0+
dc_ice*pt2(i,j)) * oneocpm
2026 qsum = qi2(i,j) + qs2(i,j)
2027 if ( qsum > 0. )
then 2028 if ( qi2(i,j) < 0. )
then 2031 elseif ( qs2(i,j) < 0. )
then 2040 if ( qs2(i,j) < 0. .and. qr2(i,j) > 0. )
then 2041 dq = min( qr2(i,j), -qs2(i,j) )
2042 qs2(i,j) = qs2(i,j) + dq
2043 qr2(i,j) = qr2(i,j) - dq
2044 pt2(i,j) = pt2(i,j) + dq*icpk(i,j)
2047 if ( qs2(i,j) < 0. .and. ql2(i,j) > 0. )
then 2048 dq = min( ql2(i,j), -qs2(i,j) )
2049 qs2(i,j) = qs2(i,j) + dq
2050 ql2(i,j) = ql2(i,j) - dq
2051 pt2(i,j) = pt2(i,j) + dq*icpk(i,j)
2054 if ( qs2(i,j) < 0. .and. qv2(i,j) > 0. )
then 2055 dq = min( 0.999*qv2(i,j), -qs2(i,j) )
2056 qs2(i,j) = qs2(i,j) + dq
2057 qv2(i,j) = qv2(i,j) - dq
2058 pt2(i,j) = pt2(i,j) + dq*(icpk(i,j)+lcpk(i,j))
2065 qsum = ql2(i,j) + qr2(i,j)
2066 if ( qsum > 0. )
then 2067 if ( qr2(i,j) < 0. )
then 2070 elseif ( ql2(i,j) < 0. )
then 2077 if ( qr(i,j,k) < 0. )
then 2079 dq = min( qi2(i,j)+qs2(i,j), -qr2(i,j) )
2080 qr2(i,j) = qr2(i,j) + dq
2081 dq1 = min( dq, qs2(i,j) )
2082 qs2(i,j) = qs2(i,j) - dq1
2083 qi2(i,j) = qi2(i,j) + dq1 - dq
2084 pt2(i,j) = pt2(i,j) - dq*icpk(i,j)
2087 if ( qr2(i,j) < 0. .and. qv2(i,j) > 0. )
then 2088 dq = min( 0.999*qv2(i,j), -qr2(i,j) )
2089 qv2(i,j) = qv2(i,j) - dq
2090 qr2(i,j) = qr2(i,j) + dq
2091 pt2(i,j) = pt2(i,j) + dq*lcpk(i,j)
2106 if ( qi2(i,j)>1.e-8 .and. pt2(i,j) >
tice )
then 2107 sink = min( qi2(i,j), (pt2(i,j)-
tice)/icpk(i,j) )
2108 ql2(i,j) = ql2(i,j) + sink
2109 qi2(i,j) = qi2(i,j) - sink
2110 pt2(i,j) = pt2(i,j) - sink*icpk(i,j)
2114 qsw = wqsat2_moist(pt2(i,j), qv2(i,j), p2(i,j), dwsdt)
2115 sink = min( ql2(i,j), (qsw-qv2(i,j))/(1.+lcpk(i,j)*dwsdt) )
2116 qv2(i,j) = qv2(i,j) + sink
2117 ql2(i,j) = ql2(i,j) - sink
2118 pt2(i,j) = pt2(i,j) - sink*lcpk(i,j)
2122 if( ql2(i,j)>1.e-8 .and. pt2(i,j) < t48 )
then 2124 sink = min( ql2(i,j), (t48-pt2(i,j))/icpk(i,j) )
2125 ql2(i,j) = ql2(i,j) - sink
2126 qi2(i,j) = qi2(i,j) + sink
2127 pt2(i,j) = pt2(i,j) + sink*icpk(i,j)
2138 qv(i,j,k) = qv2(i,j)
2139 ql(i,j,k) = ql2(i,j)
2140 qi(i,j,k) = qi2(i,j)
2141 qs(i,j,k) = qs2(i,j)
2142 qr(i,j,k) = qr2(i,j)
2143 pt(i,j,k) = pt2(i,j)
2155 dpk(i,k) = dp(i,j,k)
2159 call fillq(ie-is+1, kbot, q2, dpk)
2175 if( qv(i,j,k) < 0. )
then 2176 qv(i,j,k+1) = qv(i,j,k+1) + qv(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
2188 if( qv(i,j,k) < 0. .and. qv(i,j,k-1) > 0. )
then 2189 dq = min(-qv(i,j,k)*dp(i,j,k), qv(i,j,k-1)*dp(i,j,k-1))
2190 qv(i,j,k-1) = qv(i,j,k-1) - dq/dp(i,j,k-1)
2191 qv(i,j,k ) = qv(i,j,k ) + dq/dp(i,j,k )
2193 if( qv(i,j,k) < 0. )
then 2194 qv(i,j,k+1) = qv(i,j,k+1) + qv(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
2205 if( qv(i,j,kbot) < 0. )
then 2206 tx1 = 1.0 / dp(i,j,kbot)
2208 if ( qv(i,j,kbot)>= 0. )
goto 123
2209 if ( qv(i,j,k) > 0. )
then 2210 dq = min(-qv(i,j,kbot)*dp(i,j,kbot), qv(i,j,k)*dp(i,j,k))
2211 qv(i,j,k ) = qv(i,j,k ) - dq / dp(i,j,k)
2212 qv(i,j,kbot) = qv(i,j,kbot) + dq * tx1
2221 if (
present(qa))
then 2230 if( qa(i,j,k) < 0. )
then 2231 qa(i,j,k+1) = qa(i,j,k+1) + qa(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
2243 if( qa(i,j,kbot) < 0. .and. qa(i,j,kbot-1)>0.)
then 2244 dq = min(-qa(i,j,kbot)*dp(i,j,kbot), qa(i,j,kbot-1)*dp(i,j,kbot-1))
2245 qa(i,j,kbot-1) = qa(i,j,kbot-1) - dq/dp(i,j,kbot-1)
2246 qa(i,j,kbot ) = qa(i,j,kbot ) + dq/dp(i,j,kbot )
2249 qa(i,j,kbot) = max(0., qa(i,j,kbot))
2257 subroutine fillq(im, km, q, dp)
2259 integer,
intent(in):: im, km
2260 real,
intent(inout),
dimension(im,km):: q, dp
2262 real:: sum1, sum2, dq
2267 if ( q(i,k)>0. )
then 2268 sum1 = sum1 + q(i,k)*dp(i,k)
2271 if ( sum1<1.e-12 )
goto 500
2274 if ( q(i,k)<0.0 .and. sum1>0. )
then 2275 dq = min( sum1, -q(i,k)*dp(i,k) )
2278 q(i,k) = q(i,k) + dq/dp(i,k)
2282 if ( q(i,k)>0.0 .and. sum2>0. )
then 2283 dq = min( sum2, q(i,k)*dp(i,k) )
2285 q(i,k) = q(i,k) - dq/dp(i,k)
2290 end subroutine fillq 2292 subroutine prt_negative(qname, q, is, ie, js, je, n_g, km, threshold)
2293 character(len=*),
intent(in):: qname
2294 integer,
intent(in):: is, ie, js, je
2295 integer,
intent(in):: n_g, km
2296 real,
intent(in):: q(is-n_g:ie+n_g, js-n_g:je+n_g, km)
2297 real,
intent(in):: threshold
2305 qmin = min(qmin, q(i,j,k))
2309 call mp_reduce_min(qmin)
2310 if(is_master() .and. qmin<threshold)
write(6,*) qname,
' min (negative) = ', qmin
real, parameter dc_ice
= 2112.
pure real function, public vicpqd(q)
The module 'fv_mp_mod' is a single program multiple data (SPMD) parallel decompostion/communication m...
subroutine, public neg_adj3(is, ie, js, je, ng, kbot, hydrostatic, peln, delz, pt, dp, qv, ql, qr, qi, qs, qg, qa, check_negative)
real, parameter li0
= -2.431928e5
real, parameter c_ice
-15 C
The module 'multi_gases' peforms multi constitutents computations.
real, parameter cv_vap
1384.5
pure real function, public vicvqd_qpz(q, qpz)
pure real function, public virq_qpz(q, qpz)
The module 'fv_sg' performs FV sub-grid mixing.
pure real function, public virqd(q)
pure real function, public virq(q)
real, parameter zvir
= 0.607789855
subroutine fillq(im, km, q, dp)
subroutine, public qsmith(im, km, k1, t, p, q, qs, dqdt)
real, parameter c_liq
GFS.
real, parameter lv0
= 3.147782e6
subroutine prt_negative(qname, q, is, ie, js, je, n_g, km, threshold)
subroutine, public fv_subgrid_z(isd, ied, jsd, jed, is, ie, js, je, km, nq, dt, tau, nwat, delp, pe, peln, pkz, ta, qa, ua, va, hydrostatic, w, delz, u_dt, v_dt, t_dt, q_dt, k_bot)
pure real function, public vicpqd_qpz(q, qpz)
subroutine qs_table(n, table)
real, dimension(:), allocatable table
real, dimension(:), allocatable des
subroutine, public neg_adj2(is, ie, js, je, ng, kbot, hydrostatic, peln, delz, pt, dp, qv, ql, qr, qi, qs, qa, check_negative)
subroutine qs_table_m(n, table)
pure real function, public vicvqd(q)
real, parameter dc_vap
= -2368.