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(isd:,jsd:,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 302 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
304 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 305 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 307 cpm(i) = (1.-(q0(i,k,sphum)+q_liq))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq 308 cvm(i) = (1.-(q0(i,k,sphum)+q_liq))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq 311 elseif ( nwat==5 )
then 313 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
314 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat)
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 325 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
326 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat) + q0(i,k,graupel)
328 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 329 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 331 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 332 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 338 den(i,k) = -delp(i,j,k)/(grav*delz(i,j,k))
340 gz(i,k) = gzh(i) - g2*delz(i,j,k)
341 tmp = gz(i,k) + 0.5*(u0(i,k)**2+v0(i,k)**2+w0(i,k)**2)
342 hd(i,k) = cpm(i)*t0(i,k) + tmp
343 te(i,k) = cvm(i)*t0(i,k) + tmp
344 gzh(i) = gzh(i) - grav*delz(i,j,k)
352 if ( n==1) ratio = 0.25
353 if ( n==2) ratio = 0.5
354 if ( n==3) ratio = 0.999
356 ratio =
real(n)/
real(m)
370 elseif ( nwat==2 )
then 373 qcon(i,k) = q0(i,k,liq_wat)
376 elseif ( nwat==3 )
then 379 qcon(i,k) = q0(i,k,liq_wat) + q0(i,k,ice_wat)
382 elseif ( nwat==4 )
then 385 qcon(i,k) = q0(i,k,liq_wat) + q0(i,k,rainwat)
388 elseif ( nwat==5 )
then 391 qcon(i,k) = q0(i,k,liq_wat)+q0(i,k,ice_wat)+q0(i,k,snowwat)+q0(i,k,rainwat)
397 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)
408 tv1 = t0(i,km1)*
virq_qpz(q0(i,km1,:),qcon(i,km1))
409 tv2 = t0(i,k )*
virq_qpz(q0(i,k, :),qcon(i,k ))
411 tv1 = t0(i,km1)*(1.+xvir*q0(i,km1,sphum)-qcon(i,km1))
412 tv2 = t0(i,k )*(1.+xvir*q0(i,k ,sphum)-qcon(i,k))
414 pt1 = tv1 / pkz(i,j,km1)
415 pt2 = tv2 / pkz(i,j,k )
417 ri = (gz(i,km1)-gz(i,k))*(pt1-pt2)/( 0.5*(pt1+pt2)* &
418 ((u0(i,km1)-u0(i,k))**2+(v0(i,km1)-v0(i,k))**2+ustar2) )
419 if ( tv1>t_max .and. tv1>tv2 )
then 422 elseif ( tv2<t_min )
then 439 if ( ri < ri_ref )
then 440 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
442 h0 = mc*(q0(i,k,iq)-q0(i,km1,iq))
443 q0(i,km1,iq) = q0(i,km1,iq) + h0/delp(i,j,km1)
444 q0(i,k ,iq) = q0(i,k ,iq) - h0/delp(i,j,k )
449 elseif ( nwat==2 )
then 450 qcon(i,km1) = q0(i,km1,liq_wat)
451 elseif ( nwat==3 )
then 452 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat)
453 elseif ( nwat==4 )
then 454 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,rainwat)
455 elseif ( nwat==5 )
then 456 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat) + &
457 q0(i,km1,snowwat) + q0(i,km1,rainwat)
459 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat) + &
460 q0(i,km1,snowwat) + q0(i,km1,rainwat) + q0(i,km1,graupel)
463 h0 = mc*(u0(i,k)-u0(i,k-1))
464 u0(i,k-1) = u0(i,k-1) + h0/delp(i,j,k-1)
465 u0(i,k ) = u0(i,k ) - h0/delp(i,j,k )
467 h0 = mc*(v0(i,k)-v0(i,k-1))
468 v0(i,k-1) = v0(i,k-1) + h0/delp(i,j,k-1)
469 v0(i,k ) = v0(i,k ) - h0/delp(i,j,k )
471 if ( hydrostatic )
then 473 h0 = mc*(hd(i,k)-hd(i,k-1))
474 hd(i,k-1) = hd(i,k-1) + h0/delp(i,j,k-1)
475 hd(i,k ) = hd(i,k ) - h0/delp(i,j,k )
478 h0 = mc*(hd(i,k)-hd(i,k-1))
479 te(i,k-1) = te(i,k-1) + h0/delp(i,j,k-1)
480 te(i,k ) = te(i,k ) - h0/delp(i,j,k )
482 h0 = mc*(w0(i,k)-w0(i,k-1))
483 w0(i,k-1) = w0(i,k-1) + h0/delp(i,j,k-1)
484 w0(i,k ) = w0(i,k ) - h0/delp(i,j,k )
492 if ( hydrostatic )
then 496 rkx = cp_air/rdgas * (
vicpqd(q0(i,kk,:))/
virqd(q0(i,kk,:))) + 1
497 rdx = rdgas*
virqd(q0(i,kk,:))
499 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
500 / ( rkx - pe(i,kk,j)/pm(i,kk) )
501 gzh(i) = gzh(i) + t0(i,kk)*(peln(i,kk+1,j)-peln(i,kk,j))
502 rdx = rdx + rzx*q0(i,kk,sphum)/(1.0-sum(q0(i,kk,sphum+1:sphum+nwat-1)))
503 t0(i,kk) = t0(i,kk) / rdx
505 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
506 / ( rk - pe(i,kk,j)/pm(i,kk) )
507 gzh(i) = gzh(i) + t0(i,kk)*(peln(i,kk+1,j)-peln(i,kk,j))
508 t0(i,kk) = t0(i,kk) / ( rdgas + rz*q0(i,kk,sphum) )
514 rkx = cp_air/rdgas * (
vicpqd(q0(i,kk,:))/
virqd(q0(i,kk,:))) + 1
515 rdx = rdgas*
virqd(q0(i,kk,:))
517 rdx = rdx + rzx*q0(i,kk,sphum)/(1.0-sum(q0(i,kk,sphum+1:sphum+nwat-1)))
518 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
519 / ((rkx-pe(i,kk,j)/pm(i,kk))*rdx)
521 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
522 / ((rk-pe(i,kk,j)/pm(i,kk))*(rdgas+rz*q0(i,kk,sphum)))
528 if ( nwat == 0 )
then 531 cpm(i) = cp_air*
vicpqd(q0(i,kk,:))
532 cvm(i) = cv_air*
vicvqd(q0(i,kk,:))
538 elseif ( nwat == 1 )
then 541 cpm(i) = (1.-q0(i,kk,sphum))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor
542 cvm(i) = (1.-q0(i,kk,sphum))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap 544 cpm(i) = (1.-q0(i,kk,sphum))*cp_air + q0(i,kk,sphum)*cp_vapor
545 cvm(i) = (1.-q0(i,kk,sphum))*cv_air + q0(i,kk,sphum)*
cv_vap 548 elseif ( nwat == 2 )
then 551 q_liq = q0(i,kk,liq_wat)
552 c_air = cp_air*
vicpqd(q0(i,kk,:))
553 cpm(i) = c_air + (cp_vapor-c_air)*q0(i,kk,sphum)/(1.0-q_liq)
554 c_air = cv_air*
vicvqd(q0(i,kk,:))
555 cvm(i) = c_air + (
cv_vap-c_air)*q0(i,kk,sphum)/(1.0-q_liq)
557 cpm(i) = (1.-q0(i,kk,sphum))*cp_air + q0(i,kk,sphum)*cp_vapor
558 cvm(i) = (1.-q0(i,kk,sphum))*cv_air + q0(i,kk,sphum)*
cv_vap 561 elseif ( nwat == 3 )
then 563 q_liq = q0(i,kk,liq_wat)
564 q_sol = q0(i,kk,ice_wat)
566 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 567 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 569 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 570 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 573 elseif ( nwat == 4 )
then 575 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
577 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 578 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 580 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq 581 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq 584 elseif ( nwat == 5 )
then 586 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
587 q_sol = q0(i,kk,ice_wat) + q0(i,kk,snowwat)
589 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 590 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 592 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 593 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 598 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
599 q_sol = q0(i,kk,ice_wat) + q0(i,kk,snowwat) + q0(i,kk,graupel)
601 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 602 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 604 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 605 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 611 tv = gz(i,kk) + 0.5*(u0(i,kk)**2+v0(i,kk)**2+w0(i,kk)**2)
612 t0(i,kk) = (te(i,kk)- tv) / cvm(i)
613 hd(i,kk) = cpm(i)*t0(i,kk) + tv
624 t0(i,k) = ta(i,j,k) + (t0(i,k) - ta(i,j,k))*fra
625 u0(i,k) = ua(i,j,k) + (u0(i,k) - ua(i,j,k))*fra
626 v0(i,k) = va(i,j,k) + (v0(i,k) - va(i,j,k))*fra
630 if ( .not. hydrostatic )
then 633 w0(i,k) = w(i,j,k) + (w0(i,k) - w(i,j,k))*fra
641 q0(i,k,iq) = qa(i,j,k,iq) + (q0(i,k,iq) - qa(i,j,k,iq))*fra
649 u_dt(i,j,k) = rdt*(u0(i,k) - ua(i,j,k))
650 v_dt(i,j,k) = rdt*(v0(i,k) - va(i,j,k))
659 qa(i,j,k,iq) = q0(i,k,iq)
664 if ( .not. hydrostatic )
then 677 subroutine fv_subgrid_z( isd, ied, jsd, jed, is, ie, js, je, km, nq, dt, &
678 tau, nwat, delp, pe, peln, pkz, ta, qa, ua, va, &
679 hydrostatic, w, delz, u_dt, v_dt, t_dt, q_dt, k_bot )
682 integer,
intent(in):: is, ie, js, je, km, nq, nwat
683 integer,
intent(in):: isd, ied, jsd, jed
684 integer,
intent(in):: tau
685 real,
intent(in):: dt
686 real,
intent(in):: pe(is-1:ie+1,km+1,js-1:je+1)
687 real,
intent(in):: peln(is :ie, km+1,js :je)
688 real,
intent(in):: delp(isd:ied,jsd:jed,km)
689 real,
intent(in):: delz(isd:,jsd:,1:)
690 real,
intent(in):: pkz(is:ie,js:je,km)
691 logical,
intent(in):: hydrostatic
692 integer,
intent(in),
optional:: k_bot
694 real,
intent(inout):: ua(isd:ied,jsd:jed,km)
695 real,
intent(inout):: va(isd:ied,jsd:jed,km)
696 real,
intent(inout):: w(isd:,jsd:,1:)
697 real,
intent(inout):: ta(isd:ied,jsd:jed,km)
698 real,
intent(inout):: qa(isd:ied,jsd:jed,km,nq)
699 real,
intent(inout):: u_dt(isd:ied,jsd:jed,km)
700 real,
intent(inout):: v_dt(isd:ied,jsd:jed,km)
701 real,
intent(inout):: t_dt(is:ie,js:je,km)
702 real,
intent(inout):: q_dt(is:ie,js:je,km,nq)
704 real,
dimension(is:ie,km):: u0, v0, w0, t0, hd, te, gz, tvm, pm, den
705 real q0(is:ie,km,nq), qcon(is:ie,km)
706 real,
dimension(is:ie):: gzh, lcp2, icp2, cvm, cpm, qs
708 real :: rkx, rdx, rzx, c_air
710 real ri_ref, ri, pt1, pt2, ratio, tv, cv, tmp, q_liq, q_sol
711 real tv1, tv2, g2, h0, mc, fra, rk, rz, rdt, tvd, tv_surf
712 real dh, dq, qsw, dqsdt, tcp3
713 integer i, j, k, kk, n, m, iq, km1, im, kbot
714 real,
parameter:: ustar2 = 1.e-4
716 integer :: sphum, liq_wat, rainwat, snowwat, graupel, ice_wat, cld_amt
718 cv_air = cp_air - rdgas
719 rk = cp_air/rdgas + 1.
727 if (
present(k_bot) )
then 728 if ( k_bot < 3 )
return 734 sphum = get_tracer_index(model_atmos,
'sphum')
735 if ( nwat == 0 )
then 741 liq_wat = get_tracer_index(model_atmos,
'liq_wat')
742 ice_wat = get_tracer_index(model_atmos,
'ice_wat')
743 rainwat = get_tracer_index(model_atmos,
'rainwat')
744 snowwat = get_tracer_index(model_atmos,
'snowwat')
745 graupel = get_tracer_index(model_atmos,
'graupel')
746 cld_amt = get_tracer_index(model_atmos,
'cld_amt')
771 q0(i,k,iq) = qa(i,j,k,iq)
780 tvm(i,k) = t0(i,k)*
virq(q0(i,k,:))
782 tvm(i,k) = t0(i,k)*(1.+xvir*q0(i,k,sphum))
786 pm(i,k) = delp(i,j,k)/(peln(i,k+1,j)-peln(i,k,j))
794 if( hydrostatic )
then 798 den(i,k) = pm(i,k)/tv
799 gz(i,k) = gzh(i) + tv*(1.-pe(i,k,j)/pm(i,k))
800 hd(i,k) = cp_air*tvm(i,k)+gz(i,k)+0.5*(u0(i,k)**2+v0(i,k)**2)
801 gzh(i) = gzh(i) + tv*(peln(i,k+1,j)-peln(i,k,j))
806 if ( nwat == 0 )
then 809 cpm(i) = cp_air*
vicpqd(q0(i,k,:))
810 cvm(i) = cv_air*
vicvqd(q0(i,k,:))
816 elseif ( nwat==1 )
then 819 cpm(i) = (1.-q0(i,k,sphum))*cp_air*
vicpqd(q0(i,k,:)) + q0(i,k,sphum)*cp_vapor
820 cvm(i) = (1.-q0(i,k,sphum))*cv_air*
vicvqd(q0(i,k,:)) + q0(i,k,sphum)*
cv_vap 822 cpm(i) = (1.-q0(i,k,sphum))*cp_air + q0(i,k,sphum)*cp_vapor
823 cvm(i) = (1.-q0(i,k,sphum))*cv_air + q0(i,k,sphum)*
cv_vap 826 elseif ( nwat==2 )
then 828 q_sol = q0(i,k,liq_wat)
830 c_air = cp_air*
vicpqd(q0(i,k,:))
831 cpm(i) = c_air + (cp_vapor-c_air)*q0(i,k,sphum)/(1.0-q0(i,k,nwat))
832 c_air = cv_air*
vicvqd(q0(i,k,:))
833 cvm(i) = c_air + (
cv_vap-c_air)*q0(i,k,sphum)/(1.0-q0(i,k,nwat))
835 cpm(i) = (1.-q0(i,k,sphum))*cp_air + q0(i,k,sphum)*cp_vapor
836 cvm(i) = (1.-q0(i,k,sphum))*cv_air + q0(i,k,sphum)*
cv_vap 839 elseif ( nwat==3 )
then 841 q_liq = q0(i,k,liq_wat)
842 q_sol = q0(i,k,ice_wat)
844 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 845 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 847 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 848 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 851 elseif ( nwat==4 )
then 853 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
855 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 856 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 858 cpm(i) = (1.-(q0(i,k,sphum)+q_liq))*cp_air + q0(i,k,sphum)*cp_vapor + q_liq*
c_liq 859 cvm(i) = (1.-(q0(i,k,sphum)+q_liq))*cv_air + q0(i,k,sphum)*
cv_vap + q_liq*
c_liq 862 elseif ( nwat==5 )
then 864 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
865 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat)
867 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 868 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 870 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 871 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 876 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
877 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat) + q0(i,k,graupel)
879 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 880 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 882 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 883 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 den(i,k) = -delp(i,j,k)/(grav*delz(i,j,k))
891 gz(i,k) = gzh(i) - g2*delz(i,j,k)
892 tmp = gz(i,k) + 0.5*(u0(i,k)**2+v0(i,k)**2+w0(i,k)**2)
893 hd(i,k) = cpm(i)*t0(i,k) + tmp
894 te(i,k) = cvm(i)*t0(i,k) + tmp
895 gzh(i) = gzh(i) - grav*delz(i,j,k)
903 if ( n==1) ratio = 0.25
904 if ( n==2) ratio = 0.5
905 if ( n==3) ratio = 0.999
907 ratio =
real(n)/
real(m)
921 elseif ( nwat==2 )
then 924 qcon(i,k) = q0(i,k,liq_wat)
927 elseif ( nwat==3 )
then 930 qcon(i,k) = q0(i,k,liq_wat) + q0(i,k,ice_wat)
933 elseif ( nwat==4 )
then 936 qcon(i,k) = q0(i,k,liq_wat) + q0(i,k,rainwat)
939 elseif ( nwat==5 )
then 942 qcon(i,k) = q0(i,k,liq_wat)+q0(i,k,ice_wat)+q0(i,k,snowwat)+q0(i,k,rainwat)
947 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)
956 if(nwat>0)
call qsmith(im, km, im, 1, t0(is,km1), pm(is,km1), q0(is,km1,sphum), qs)
958 if(nwat>0)
call qsmith(im, 1, 1, t0(is,km1), pm(is,km1), q0(is,km1,sphum), qs)
965 tv1 = t0(i,km1)*
virq_qpz(q0(i,km1,:),qcon(i,km1))
966 tv2 = t0(i,k )*
virq_qpz(q0(i,k ,:),qcon(i,k ))
968 tv1 = t0(i,km1)*(1.+xvir*q0(i,km1,sphum)-qcon(i,km1))
969 tv2 = t0(i,k )*(1.+xvir*q0(i,k ,sphum)-qcon(i,k))
971 pt1 = tv1 / pkz(i,j,km1)
972 pt2 = tv2 / pkz(i,j,k )
973 ri = (gz(i,km1)-gz(i,k))*(pt1-pt2)/( 0.5*(pt1+pt2)* &
974 ((u0(i,km1)-u0(i,k))**2+(v0(i,km1)-v0(i,k))**2+ustar2) )
982 h0 = hd(i,k)-hd(i,km1) +
hlv0*(q0(i,k,sphum)-qs(i))
983 if ( h0 > 0. ) ri_ref = 5.*ri_ref
986 if ( ri < ri_ref )
then 987 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
989 h0 = mc*(q0(i,k,iq)-q0(i,km1,iq))
990 q0(i,km1,iq) = q0(i,km1,iq) + h0/delp(i,j,km1)
991 q0(i,k ,iq) = q0(i,k ,iq) - h0/delp(i,j,k )
996 elseif ( nwat==2 )
then 997 qcon(i,km1) = q0(i,km1,liq_wat)
998 elseif ( nwat==3 )
then 999 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat)
1000 elseif ( nwat==4 )
then 1001 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,rainwat)
1002 elseif ( nwat==5 )
then 1003 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat) + &
1004 q0(i,km1,snowwat) + q0(i,km1,rainwat)
1006 qcon(i,km1) = q0(i,km1,liq_wat) + q0(i,km1,ice_wat) + &
1007 q0(i,km1,snowwat) + q0(i,km1,rainwat) + q0(i,km1,graupel)
1010 h0 = mc*(u0(i,k)-u0(i,k-1))
1011 u0(i,k-1) = u0(i,k-1) + h0/delp(i,j,k-1)
1012 u0(i,k ) = u0(i,k ) - h0/delp(i,j,k )
1014 h0 = mc*(v0(i,k)-v0(i,k-1))
1015 v0(i,k-1) = v0(i,k-1) + h0/delp(i,j,k-1)
1016 v0(i,k ) = v0(i,k ) - h0/delp(i,j,k )
1018 if ( hydrostatic )
then 1020 h0 = mc*(hd(i,k)-hd(i,k-1))
1021 hd(i,k-1) = hd(i,k-1) + h0/delp(i,j,k-1)
1022 hd(i,k ) = hd(i,k ) - h0/delp(i,j,k )
1025 h0 = mc*(hd(i,k)-hd(i,k-1))
1026 te(i,k-1) = te(i,k-1) + h0/delp(i,j,k-1)
1027 te(i,k ) = te(i,k ) - h0/delp(i,j,k )
1029 h0 = mc*(w0(i,k)-w0(i,k-1))
1030 w0(i,k-1) = w0(i,k-1) + h0/delp(i,j,k-1)
1031 w0(i,k ) = w0(i,k ) - h0/delp(i,j,k )
1039 if ( hydrostatic )
then 1043 rkx = cp_air/rdgas * (
vicpqd(q0(i,kk,:))/
virqd(q0(i,kk,:))) + 1
1044 rdx = rdgas*
virqd(q0(i,kk,:))
1046 rdx = rdx + rzx*q0(i,kk,sphum)/(1.0-sum(q0(i,kk,sphum+1:sphum+nwat-1)))
1047 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
1048 / ( rkx - pe(i,kk,j)/pm(i,kk) )
1049 gzh(i) = gzh(i) + t0(i,kk)*(peln(i,kk+1,j)-peln(i,kk,j))
1050 t0(i,kk) = t0(i,kk) / ( rdx + rzx*q0(i,kk,sphum) )
1052 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
1053 / ( rk - pe(i,kk,j)/pm(i,kk) )
1054 gzh(i) = gzh(i) + t0(i,kk)*(peln(i,kk+1,j)-peln(i,kk,j))
1055 t0(i,kk) = t0(i,kk) / ( rdgas + rz*q0(i,kk,sphum) )
1061 rkx = cp_air/rdgas * (
vicpqd(q0(i,kk,:))/
virqd(q0(i,kk,:))) + 1
1062 rdx = rdgas*
virqd(i)
1064 rdx = rdx + rzx*q0(i,kk,sphum)/(1.0-sum(q0(i,kk,sphum+1:sphum+nwat-1)))
1065 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
1066 / ((rkx-pe(i,kk,j)/pm(i,kk))*rdx)
1068 t0(i,kk) = (hd(i,kk)-gzh(i)-0.5*(u0(i,kk)**2+v0(i,kk)**2)) &
1069 / ((rk-pe(i,kk,j)/pm(i,kk))*(rdgas+rz*q0(i,kk,sphum)))
1075 if ( nwat == 0 )
then 1078 cpm(i) = cp_air*
vicpqd(q0(i,kk,:))
1079 cvm(i) = cv_air*
vicvqd(q0(i,kk,:))
1085 elseif ( nwat == 1 )
then 1088 cpm(i) = (1.-q0(i,kk,sphum))*cp_air*
vicpqd(q0(i,kk,:)) + q0(i,kk,sphum)*cp_vapor
1089 cvm(i) = (1.-q0(i,kk,sphum))*cv_air*
vicvqd(q0(i,kk,:)) + q0(i,kk,sphum)*
cv_vap 1091 cpm(i) = (1.-q0(i,kk,sphum))*cp_air + q0(i,kk,sphum)*cp_vapor
1092 cvm(i) = (1.-q0(i,kk,sphum))*cv_air + q0(i,kk,sphum)*
cv_vap 1095 elseif ( nwat == 2 )
then 1098 c_air = cp_air*
vicpqd(q0(i,kk,:))
1099 cpm(i) = c_air+(cp_vapor-c_air)*q0(i,kk,sphum)/(1.0-q0(i,kk,sphum+nwat-1))
1100 c_air = cv_air*
vicvqd(q0(i,kk,:))
1101 cvm(i) = c_air+(
cv_vap -c_air)*q0(i,kk,sphum)/(1.0-q0(i,kk,sphum+nwat-1))
1103 cpm(i) = (1.-q0(i,kk,sphum))*cp_air + q0(i,kk,sphum)*cp_vapor
1104 cvm(i) = (1.-q0(i,kk,sphum))*cv_air + q0(i,kk,sphum)*
cv_vap 1107 elseif ( nwat == 3 )
then 1109 q_liq = q0(i,kk,liq_wat)
1110 q_sol = q0(i,kk,ice_wat)
1112 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 1113 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 1115 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 1116 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 1119 elseif ( nwat == 4 )
then 1121 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
1123 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 1124 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 1126 cpm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cp_air + q0(i,kk,sphum)*cp_vapor + q_liq*
c_liq 1127 cvm(i) = (1.-(q0(i,kk,sphum)+q_liq))*cv_air + q0(i,kk,sphum)*
cv_vap + q_liq*
c_liq 1130 elseif ( nwat == 5 )
then 1132 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
1133 q_sol = q0(i,kk,ice_wat) + q0(i,kk,snowwat)
1135 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 1136 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 1138 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 1139 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 1144 q_liq = q0(i,kk,liq_wat) + q0(i,kk,rainwat)
1145 q_sol = q0(i,kk,ice_wat) + q0(i,kk,snowwat) + q0(i,kk,graupel)
1147 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 1148 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 1150 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 1151 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 1158 tv = gz(i,kk) + 0.5*(u0(i,kk)**2+v0(i,kk)**2+w0(i,kk)**2)
1159 t0(i,kk) = (te(i,kk)- tv) / cvm(i)
1160 hd(i,kk) = cpm(i)*t0(i,kk) + tv
1168 if ( fra < 1. )
then 1171 t0(i,k) = ta(i,j,k) + (t0(i,k) - ta(i,j,k))*fra
1172 u0(i,k) = ua(i,j,k) + (u0(i,k) - ua(i,j,k))*fra
1173 v0(i,k) = va(i,j,k) + (v0(i,k) - va(i,j,k))*fra
1177 if ( .not. hydrostatic )
then 1180 w0(i,k) = w(i,j,k) + (w0(i,k) - w(i,j,k))*fra
1188 q0(i,k,iq) = qa(i,j,k,iq) + (q0(i,k,iq) - qa(i,j,k,iq))*fra
1198 if ( nwat > 5 )
then 1200 if ( hydrostatic )
then 1204 den(i,k) = pm(i,k)/(rdgas*t0(i,k)*
virq(q0(i,k,:)))
1206 den(i,k) = pm(i,k)/(rdgas*t0(i,k)*(1.+xvir*q0(i,k,sphum)))
1208 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
1209 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat) + q0(i,k,graupel)
1211 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 1213 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 1215 lcp2(i) = hlv / cpm(i)
1216 icp2(i) = hlf / cpm(i)
1220 den(i,k) = -delp(i,j,k)/(grav*delz(i,j,k))
1221 q_liq = q0(i,k,liq_wat) + q0(i,k,rainwat)
1222 q_sol = q0(i,k,ice_wat) + q0(i,k,snowwat) + q0(i,k,graupel)
1224 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 1226 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 1235 qsw = wqs2(t0(i,k), den(i,k), dqsdt)
1236 dq = q0(i,k,sphum) - qsw
1238 tcp3 = lcp2(i) + icp2(i)*min(1., dim(
tice,t0(i,k))/40.)
1239 tmp = dq/(1.+tcp3*dqsdt)
1240 t0(i,k) = t0(i,k) + tmp*lcp2(i)
1241 q0(i,k, sphum) = q0(i,k, sphum) - tmp
1242 q0(i,k,liq_wat) = q0(i,k,liq_wat) + tmp
1244 if (cld_amt .gt. 0)
then 1245 qa(i,j,k,cld_amt) = max(0.5, min(1., qa(i,j,k,cld_amt)+25.*tmp/qsw))
1249 tmp =
tice-40. - t0(i,k)
1250 if( tmp>0.0 .and. q0(i,k,liq_wat)>0. )
then 1251 dh = min( q0(i,k,liq_wat), q0(i,k,liq_wat)*tmp*0.125, tmp/icp2(i) )
1252 q0(i,k,liq_wat) = q0(i,k,liq_wat) - dh
1253 q0(i,k,ice_wat) = q0(i,k,ice_wat) + dh
1254 t0(i,k) = t0(i,k) + dh*icp2(i)
1263 u_dt(i,j,k) = rdt*(u0(i,k) - ua(i,j,k))
1264 v_dt(i,j,k) = rdt*(v0(i,k) - va(i,j,k))
1272 if (iq .ne. cld_amt )
then 1274 qa(i,j,k,iq) = q0(i,k,iq)
1280 if ( .not. hydrostatic )
then 1296 integer,
parameter:: length=2621
1299 if( .not.
allocated(
table) )
then 1302 allocate (
table(length) )
1303 allocate (
des(length) )
1310 des(length) =
des(length-1)
1317 subroutine qsmith(im, km, imx, kmx, t, p, q, qs, dqdt)
1319 subroutine qsmith(im, km, k1, t, p, q, qs, dqdt)
1322 real,
intent(in),
dimension(im,km):: t, p
1324 real,
intent(in),
dimension(im,km,*):: q
1325 integer,
intent(in):: im, km, imx, kmx
1327 integer,
intent(in):: im, km, k1
1328 real,
intent(in),
dimension(im,km):: q
1330 real,
intent(out),
dimension(im,km):: qs
1331 real,
intent(out),
optional:: dqdt(im,km)
1350 ap1 = 10.*dim(t(i,k), tmin) + 1.
1351 ap1 = min(2621., ap1)
1353 es(i,k) =
table(it) + (ap1-it)*
des(it)
1357 qs(i,k) =
esl*es(i,k)*(1.+
zvir*q(i,k))/p(i,k)
1362 if (
present(dqdt) )
then 1370 ap1 = 10.*dim(t(i,k), tmin) + 1.
1371 ap1 = min(2621., ap1) - 0.5
1376 dqdt(i,k) = eps10*(
des(it)+(ap1-it)*(
des(it+1)-
des(it)))*(1.+
zvir*q(i,k))/p(i,k)
1386 integer,
intent(in):: n
1389 real esbasw, tbasw, esbasi, tbasi, Tmin, tem, aa, b, c, d, e, esh20
1400 tem = tmin+dt*
real(i-1)
1401 aa = -7.90298*(tbasw/tem-1)
1402 b = 5.02808*alog10(tbasw/tem)
1403 c = -1.3816e-07*(10**((1-tem/tbasw)*11.344)-1)
1404 d = 8.1328e-03*(10**((tbasw/tem-1)*(-3.49149))-1)
1406 table(i) = 0.1*10**(aa+b+c+d+e)
1413 integer,
intent(in):: n
1417 real esbasw, tbasw, esbasi, tbasi, Tmin, tem, aa, b, c, d, e, esh20
1431 tem = tmin+dt*
real(i-1)
1432 aa = -9.09718 *(tbasi/tem-1.0)
1433 b = -3.56654 *alog10(tbasi/tem)
1434 c = 0.876793*(1.0-tem/tbasi)
1436 table(i)=10**(aa+b+c+e)
1442 tem = 253.16+dt*
real(i-1)
1443 aa = -7.90298*(tbasw/tem-1)
1444 b = 5.02808*alog10(tbasw/tem)
1445 c = -1.3816e-07*(10**((1-tem/tbasw)*11.344)-1)
1446 d = 8.1328e-03*(10**((tbasw/tem-1)*(-3.49149))-1)
1448 esh20 = 10**(aa+b+c+d+e)
1452 table(i+1400) = esh20
1458 tem = 253.16+dt*
real(i-1)
1459 wice = 0.05*(273.16-tem)
1460 wh2o = 0.05*(tem-253.16)
1461 table(i+1400) = wice*table(i+1400)+wh2o*esupc(i)
1465 table(i) = table(i)*0.1
1470 subroutine neg_adj3(is, ie, js, je, ng, kbot, hydrostatic, peln, delz, pt, dp, &
1476 ql, qr, qi, qs, qg, qa, check_negative)
1479 integer,
intent(in):: is, ie, js, je, ng, kbot
1480 logical,
intent(in):: hydrostatic
1481 real,
intent(in):: dp(is-ng:ie+ng,js-ng:je+ng,kbot)
1482 real,
intent(in):: delz(is-ng:,js-ng:,1:)
1483 real,
intent(in):: peln(is:ie,kbot+1,js:je)
1484 logical,
intent(in),
OPTIONAL :: check_negative
1486 real,
intent(inout),
dimension(is-ng:ie+ng,js-ng:je+ng,kbot,*):: qvi
1488 real,
intent(inout),
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: qv
1490 real,
intent(inout),
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: &
1491 pt, ql, qr, qi, qs, qg
1492 real,
intent(inout),
OPTIONAL,
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: qa
1494 logical:: sat_adj = .false.
1495 real,
parameter :: t48 =
tice - 48.
1497 real,
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: qv
1499 real,
dimension(is:ie,kbot):: dpk, q2
1500 real,
dimension(is:ie,js:je):: pt2, qv2, ql2, qi2, qs2, qr2, qg2, dp2, p2, icpk, lcpk
1502 real:: dq, qsum, dq1, q_liq, q_sol, cpm, sink, qsw, dwsdt
1505 cv_air = cp_air - rdgas
1508 qv(:,:,:) = qvi(:,:,:,1)
1511 if (
present(check_negative) )
then 1512 if ( check_negative )
then 1513 call prt_negative(
'Temperature', pt, is, ie, js, je, ng, kbot, 165.)
1514 call prt_negative(
'sphum', qv, is, ie, js, je, ng, kbot, -1.e-8)
1515 call prt_negative(
'liq_wat', ql, is, ie, js, je, ng, kbot, -1.e-7)
1516 call prt_negative(
'rainwat', qr, is, ie, js, je, ng, kbot, -1.e-7)
1517 call prt_negative(
'ice_wat', qi, is, ie, js, je, ng, kbot, -1.e-7)
1518 call prt_negative(
'snowwat', qs, is, ie, js, je, ng, kbot, -1.e-7)
1519 call prt_negative(
'graupel', qg, is, ie, js, je, ng, kbot, -1.e-7)
1521 qvi(:,:,:,1) = qv(:,:,:)
1526 if ( hydrostatic )
then 1543 qv2(i,j) = qv(i,j,k)
1544 ql2(i,j) = ql(i,j,k)
1545 qi2(i,j) = qi(i,j,k)
1546 qs2(i,j) = qs(i,j,k)
1547 qr2(i,j) = qr(i,j,k)
1548 qg2(i,j) = qg(i,j,k)
1549 dp2(i,j) = dp(i,j,k)
1550 pt2(i,j) = pt(i,j,k)
1554 if ( hydrostatic )
then 1557 p2(i,j) = dp2(i,j)/(peln(i,k+1,j)-peln(i,k,j))
1558 q_liq = max(0., ql2(i,j) + qr2(i,j))
1559 q_sol = max(0., qi2(i,j) + qs2(i,j))
1561 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) + &
1562 qv2(i,j)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1564 cpm = (1.-(qv2(i,j)+q_liq+q_sol))*cp_air + qv2(i,j)*cp_vapor + q_liq*
c_liq + q_sol*
c_ice 1566 lcpk(i,j) = hlv / cpm
1567 icpk(i,j) = hlf / cpm
1573 q_liq = max(0., ql2(i,j) + qr2(i,j))
1574 q_sol = max(0., qi2(i,j) + qs2(i,j))
1576 p2(i,j) = -dp2(i,j)/(grav*delz(i,j,k))*rdgas*pt2(i,j)*
virq(qvi(i,j,k,1:
num_gas))
1577 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) + &
1580 p2(i,j) = -dp2(i,j)/(grav*delz(i,j,k))*rdgas*pt2(i,j)*(1.+
zvir*qv2(i,j))
1581 cpm = (1.-(qv2(i,j)+q_liq+q_sol))*cv_air + qv2(i,j)*
cv_vap + q_liq*
c_liq + q_sol*
c_ice 1595 qsum = qi2(i,j) + qs2(i,j)
1596 if ( qsum > 0. )
then 1597 if ( qi2(i,j) < 0. )
then 1600 elseif ( qs2(i,j) < 0. )
then 1608 qg2(i,j) = qg2(i,j) + qsum
1613 if ( qg2(i,j) < 0. )
then 1614 dq = min( qs2(i,j), -qg2(i,j) )
1615 qs2(i,j) = qs2(i,j) - dq
1616 qg2(i,j) = qg2(i,j) + dq
1617 if ( qg2(i,j) < 0. )
then 1619 dq = min( qi2(i,j), -qg2(i,j) )
1620 qi2(i,j) = qi2(i,j) - dq
1621 qg2(i,j) = qg2(i,j) + dq
1626 if ( qg2(i,j)<0. .and. qr2(i,j)>0. )
then 1627 dq = min( qr2(i,j), -qg2(i,j) )
1628 qg2(i,j) = qg2(i,j) + dq
1629 qr2(i,j) = qr2(i,j) - dq
1630 pt2(i,j) = pt2(i,j) + dq*icpk(i,j)
1633 if ( qg2(i,j)<0. .and. ql2(i,j)>0. )
then 1634 dq = min( ql2(i,j), -qg2(i,j) )
1635 qg2(i,j) = qg2(i,j) + dq
1636 ql2(i,j) = ql2(i,j) - dq
1637 pt2(i,j) = pt2(i,j) + dq*icpk(i,j)
1640 if ( qg2(i,j)<0. .and. qv2(i,j)>0. )
then 1641 dq = min( 0.999*qv2(i,j), -qg2(i,j) )
1642 qg2(i,j) = qg2(i,j) + dq
1643 qv2(i,j) = qv2(i,j) - dq
1644 pt2(i,j) = pt2(i,j) + dq*(icpk(i,j)+lcpk(i,j))
1650 qsum = ql2(i,j) + qr2(i,j)
1651 if ( qsum > 0. )
then 1652 if ( qr2(i,j) < 0. )
then 1655 elseif ( ql2(i,j) < 0. )
then 1663 dq = min( max(0.0, qg2(i,j)), -qr2(i,j) )
1664 qr2(i,j) = qr2(i,j) + dq
1665 qg2(i,j) = qg2(i,j) - dq
1666 pt2(i,j) = pt2(i,j) - dq*icpk(i,j)
1667 if ( qr(i,j,k) < 0. )
then 1669 dq = min( qi2(i,j)+qs2(i,j), -qr2(i,j) )
1670 qr2(i,j) = qr2(i,j) + dq
1671 dq1 = min( dq, qs2(i,j) )
1672 qs2(i,j) = qs2(i,j) - dq1
1673 qi2(i,j) = qi2(i,j) + dq1 - dq
1674 pt2(i,j) = pt2(i,j) - dq*icpk(i,j)
1677 if ( qr2(i,j)<0. .and. qv2(i,j)>0. )
then 1678 dq = min( 0.999*qv2(i,j), -qr2(i,j) )
1679 qv2(i,j) = qv2(i,j) - dq
1680 qr2(i,j) = qr2(i,j) + dq
1681 pt2(i,j) = pt2(i,j) + dq*lcpk(i,j)
1696 if ( qi2(i,j)>1.e-8 .and. pt2(i,j) >
tice )
then 1697 sink = min( qi2(i,j), (pt2(i,j)-
tice)/icpk(i,j) )
1698 ql2(i,j) = ql2(i,j) + sink
1699 qi2(i,j) = qi2(i,j) - sink
1700 pt2(i,j) = pt2(i,j) - sink*icpk(i,j)
1704 qsw = wqsat2_moist(pt2(i,j), qv2(i,j), p2(i,j), dwsdt)
1705 sink = min( ql2(i,j), (qsw-qv2(i,j))/(1.+lcpk(i,j)*dwsdt) )
1706 qv2(i,j) = qv2(i,j) + sink
1707 ql2(i,j) = ql2(i,j) - sink
1708 pt2(i,j) = pt2(i,j) - sink*lcpk(i,j)
1712 if( ql2(i,j)>1.e-8 .and. pt2(i,j) < t48 )
then 1714 sink = min( ql2(i,j), (t48-pt2(i,j))/icpk(i,j) )
1715 ql2(i,j) = ql2(i,j) - sink
1716 qi2(i,j) = qi2(i,j) + sink
1717 pt2(i,j) = pt2(i,j) + sink*icpk(i,j)
1728 qv(i,j,k) = qv2(i,j)
1729 ql(i,j,k) = ql2(i,j)
1730 qi(i,j,k) = qi2(i,j)
1731 qs(i,j,k) = qs2(i,j)
1732 qr(i,j,k) = qr2(i,j)
1733 qg(i,j,k) = qg2(i,j)
1734 pt(i,j,k) = pt2(i,j)
1738 qvi(:,:,k,1) = qv(:,:,k)
1749 dpk(i,k) = dp(i,j,k)
1753 call fillq(ie-is+1, kbot, q2, dpk)
1765 call fillq(ie-is+1, kbot, q2, dpk)
1781 if( qv(i,j,k) < 0. )
then 1782 qv(i,j,k+1) = qv(i,j,k+1) + qv(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
1794 if( qv(i,j,k) < 0. .and. qv(i,j,k-1) > 0. )
then 1795 dq = min(-qv(i,j,k)*dp(i,j,k), qv(i,j,k-1)*dp(i,j,k-1))
1796 qv(i,j,k-1) = qv(i,j,k-1) - dq/dp(i,j,k-1)
1797 qv(i,j,k ) = qv(i,j,k ) + dq/dp(i,j,k )
1799 if( qv(i,j,k) < 0. )
then 1800 qv(i,j,k+1) = qv(i,j,k+1) + qv(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
1811 if( qv(i,j,kbot) < 0. )
then 1813 if ( qv(i,j,kbot)>=0. )
goto 123
1814 if ( qv(i,j,k) > 0. )
then 1815 dq = min(-qv(i,j,kbot)*dp(i,j,kbot), qv(i,j,k)*dp(i,j,k))
1816 qv(i,j,k ) = qv(i,j,k ) - dq/dp(i,j,k)
1817 qv(i,j,kbot) = qv(i,j,kbot) + dq/dp(i,j,kbot)
1825 qvi(:,:,:,1) = qv(:,:,:)
1829 if (
present(qa))
then 1838 if( qa(i,j,k) < 0. )
then 1839 qa(i,j,k+1) = qa(i,j,k+1) + qa(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
1851 if( qa(i,j,kbot) < 0. .and. qa(i,j,kbot-1)>0.)
then 1852 dq = min(-qa(i,j,kbot)*dp(i,j,kbot), qa(i,j,kbot-1)*dp(i,j,kbot-1))
1853 qa(i,j,kbot-1) = qa(i,j,kbot-1) - dq/dp(i,j,kbot-1)
1854 qa(i,j,kbot ) = qa(i,j,kbot ) + dq/dp(i,j,kbot )
1857 qa(i,j,kbot) = max(0., qa(i,j,kbot))
1865 subroutine neg_adj2(is, ie, js, je, ng, kbot, hydrostatic, &
1866 peln, delz, pt, dp, qv, ql, qr, qi, qs, qa, check_negative)
1869 integer,
intent(in):: is, ie, js, je, ng, kbot
1870 logical,
intent(in):: hydrostatic
1871 real,
intent(in):: dp(is-ng:ie+ng,js-ng:je+ng,kbot)
1872 real,
intent(in):: delz(is-ng:,js-ng:,1:)
1873 real,
intent(in):: peln(is:ie,kbot+1,js:je)
1874 logical,
intent(in),
OPTIONAL :: check_negative
1875 real,
intent(inout),
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: &
1876 pt, qv, ql, qr, qi, qs
1877 real,
intent(inout),
OPTIONAL,
dimension(is-ng:ie+ng,js-ng:je+ng,kbot):: qa
1879 logical:: sat_adj = .false.
1880 real,
parameter :: t48 =
tice - 48.
1881 real,
dimension(is:ie,kbot):: dpk, q2
1882 real,
dimension(is:ie,js:je):: pt2, qv2, ql2, qi2, qs2, qr2, dp2, p2, icpk, lcpk
1884 real:: dq, qsum, dq1, q_liq, q_sol, oneocpm, sink, qsw, dwsdt, tx1
1887 cv_air = cp_air - rdgas
1889 if (
present(check_negative) )
then 1890 if ( check_negative )
then 1891 call prt_negative(
'Temperature', pt, is, ie, js, je, ng, kbot, 165.)
1892 call prt_negative(
'sphum', qv, is, ie, js, je, ng, kbot, -1.e-8)
1893 call prt_negative(
'liq_wat', ql, is, ie, js, je, ng, kbot, -1.e-7)
1894 call prt_negative(
'rainwat', qr, is, ie, js, je, ng, kbot, -1.e-7)
1895 call prt_negative(
'ice_wat', qi, is, ie, js, je, ng, kbot, -1.e-7)
1896 call prt_negative(
'snowwat', qs, is, ie, js, je, ng, kbot, -1.e-7)
1900 if ( hydrostatic )
then 1914 qv2(i,j) = qv(i,j,k)
1915 ql2(i,j) = ql(i,j,k)
1916 qi2(i,j) = qi(i,j,k)
1917 qs2(i,j) = qs(i,j,k)
1918 qr2(i,j) = qr(i,j,k)
1919 dp2(i,j) = dp(i,j,k)
1920 pt2(i,j) = pt(i,j,k)
1924 if ( hydrostatic )
then 1927 p2(i,j) = dp2(i,j)/(peln(i,k+1,j)-peln(i,k,j))
1928 q_liq = max(0., ql2(i,j) + qr2(i,j))
1929 q_sol = max(0., qi2(i,j) + qs2(i,j))
1930 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)
1931 lcpk(i,j) = hlv * oneocpm
1932 icpk(i,j) = hlf * oneocpm
1938 p2(i,j) = -dp2(i,j)/(grav*delz(i,j,k))*rdgas*pt2(i,j)*(1.+
zvir*qv2(i,j))
1939 q_liq = max(0., ql2(i,j) + qr2(i,j))
1940 q_sol = max(0., qi2(i,j) + qs2(i,j))
1941 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)
1943 icpk(i,j) = (
li0+
dc_ice*pt2(i,j)) * oneocpm
1954 qsum = qi2(i,j) + qs2(i,j)
1955 if ( qsum > 0. )
then 1956 if ( qi2(i,j) < 0. )
then 1959 elseif ( qs2(i,j) < 0. )
then 1968 if ( qs2(i,j) < 0. .and. qr2(i,j) > 0. )
then 1969 dq = min( qr2(i,j), -qs2(i,j) )
1970 qs2(i,j) = qs2(i,j) + dq
1971 qr2(i,j) = qr2(i,j) - dq
1972 pt2(i,j) = pt2(i,j) + dq*icpk(i,j)
1975 if ( qs2(i,j) < 0. .and. ql2(i,j) > 0. )
then 1976 dq = min( ql2(i,j), -qs2(i,j) )
1977 qs2(i,j) = qs2(i,j) + dq
1978 ql2(i,j) = ql2(i,j) - dq
1979 pt2(i,j) = pt2(i,j) + dq*icpk(i,j)
1982 if ( qs2(i,j) < 0. .and. qv2(i,j) > 0. )
then 1983 dq = min( 0.999*qv2(i,j), -qs2(i,j) )
1984 qs2(i,j) = qs2(i,j) + dq
1985 qv2(i,j) = qv2(i,j) - dq
1986 pt2(i,j) = pt2(i,j) + dq*(icpk(i,j)+lcpk(i,j))
1993 qsum = ql2(i,j) + qr2(i,j)
1994 if ( qsum > 0. )
then 1995 if ( qr2(i,j) < 0. )
then 1998 elseif ( ql2(i,j) < 0. )
then 2005 if ( qr(i,j,k) < 0. )
then 2007 dq = min( qi2(i,j)+qs2(i,j), -qr2(i,j) )
2008 qr2(i,j) = qr2(i,j) + dq
2009 dq1 = min( dq, qs2(i,j) )
2010 qs2(i,j) = qs2(i,j) - dq1
2011 qi2(i,j) = qi2(i,j) + dq1 - dq
2012 pt2(i,j) = pt2(i,j) - dq*icpk(i,j)
2015 if ( qr2(i,j) < 0. .and. qv2(i,j) > 0. )
then 2016 dq = min( 0.999*qv2(i,j), -qr2(i,j) )
2017 qv2(i,j) = qv2(i,j) - dq
2018 qr2(i,j) = qr2(i,j) + dq
2019 pt2(i,j) = pt2(i,j) + dq*lcpk(i,j)
2034 if ( qi2(i,j)>1.e-8 .and. pt2(i,j) >
tice )
then 2035 sink = min( qi2(i,j), (pt2(i,j)-
tice)/icpk(i,j) )
2036 ql2(i,j) = ql2(i,j) + sink
2037 qi2(i,j) = qi2(i,j) - sink
2038 pt2(i,j) = pt2(i,j) - sink*icpk(i,j)
2042 qsw = wqsat2_moist(pt2(i,j), qv2(i,j), p2(i,j), dwsdt)
2043 sink = min( ql2(i,j), (qsw-qv2(i,j))/(1.+lcpk(i,j)*dwsdt) )
2044 qv2(i,j) = qv2(i,j) + sink
2045 ql2(i,j) = ql2(i,j) - sink
2046 pt2(i,j) = pt2(i,j) - sink*lcpk(i,j)
2050 if( ql2(i,j)>1.e-8 .and. pt2(i,j) < t48 )
then 2052 sink = min( ql2(i,j), (t48-pt2(i,j))/icpk(i,j) )
2053 ql2(i,j) = ql2(i,j) - sink
2054 qi2(i,j) = qi2(i,j) + sink
2055 pt2(i,j) = pt2(i,j) + sink*icpk(i,j)
2066 qv(i,j,k) = qv2(i,j)
2067 ql(i,j,k) = ql2(i,j)
2068 qi(i,j,k) = qi2(i,j)
2069 qs(i,j,k) = qs2(i,j)
2070 qr(i,j,k) = qr2(i,j)
2071 pt(i,j,k) = pt2(i,j)
2083 dpk(i,k) = dp(i,j,k)
2087 call fillq(ie-is+1, kbot, q2, dpk)
2103 if( qv(i,j,k) < 0. )
then 2104 qv(i,j,k+1) = qv(i,j,k+1) + qv(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
2116 if( qv(i,j,k) < 0. .and. qv(i,j,k-1) > 0. )
then 2117 dq = min(-qv(i,j,k)*dp(i,j,k), qv(i,j,k-1)*dp(i,j,k-1))
2118 qv(i,j,k-1) = qv(i,j,k-1) - dq/dp(i,j,k-1)
2119 qv(i,j,k ) = qv(i,j,k ) + dq/dp(i,j,k )
2121 if( qv(i,j,k) < 0. )
then 2122 qv(i,j,k+1) = qv(i,j,k+1) + qv(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
2133 if( qv(i,j,kbot) < 0. )
then 2134 tx1 = 1.0 / dp(i,j,kbot)
2136 if ( qv(i,j,kbot)>= 0. )
goto 123
2137 if ( qv(i,j,k) > 0. )
then 2138 dq = min(-qv(i,j,kbot)*dp(i,j,kbot), qv(i,j,k)*dp(i,j,k))
2139 qv(i,j,k ) = qv(i,j,k ) - dq / dp(i,j,k)
2140 qv(i,j,kbot) = qv(i,j,kbot) + dq * tx1
2149 if (
present(qa))
then 2158 if( qa(i,j,k) < 0. )
then 2159 qa(i,j,k+1) = qa(i,j,k+1) + qa(i,j,k)*dp(i,j,k)/dp(i,j,k+1)
2171 if( qa(i,j,kbot) < 0. .and. qa(i,j,kbot-1)>0.)
then 2172 dq = min(-qa(i,j,kbot)*dp(i,j,kbot), qa(i,j,kbot-1)*dp(i,j,kbot-1))
2173 qa(i,j,kbot-1) = qa(i,j,kbot-1) - dq/dp(i,j,kbot-1)
2174 qa(i,j,kbot ) = qa(i,j,kbot ) + dq/dp(i,j,kbot )
2177 qa(i,j,kbot) = max(0., qa(i,j,kbot))
2185 subroutine fillq(im, km, q, dp)
2187 integer,
intent(in):: im, km
2188 real,
intent(inout),
dimension(im,km):: q, dp
2190 real:: sum1, sum2, dq
2195 if ( q(i,k)>0. )
then 2196 sum1 = sum1 + q(i,k)*dp(i,k)
2199 if ( sum1<1.e-12 )
goto 500
2202 if ( q(i,k)<0.0 .and. sum1>0. )
then 2203 dq = min( sum1, -q(i,k)*dp(i,k) )
2206 q(i,k) = q(i,k) + dq/dp(i,k)
2210 if ( q(i,k)>0.0 .and. sum2>0. )
then 2211 dq = min( sum2, q(i,k)*dp(i,k) )
2213 q(i,k) = q(i,k) - dq/dp(i,k)
2218 end subroutine fillq 2220 subroutine prt_negative(qname, q, is, ie, js, je, n_g, km, threshold)
2221 character(len=*),
intent(in):: qname
2222 integer,
intent(in):: is, ie, js, je
2223 integer,
intent(in):: n_g, km
2224 real,
intent(in):: q(is-n_g:ie+n_g, js-n_g:je+n_g, km)
2225 real,
intent(in):: threshold
2233 qmin = min(qmin, q(i,j,k))
2237 call mp_reduce_min(qmin)
2238 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.