57 use constants_mod
, only: rvgas, rdgas, grav, hlv, hlf, cp_air
60 use gfdl_cloud_microphys_mod
, only: ql_gen, qi_gen, qi0_max, ql_mlt, ql0_max, qi_lim, qs_mlt
61 use gfdl_cloud_microphys_mod
, only: icloud_f, sat_adj0, t_sub, cld_min
62 use gfdl_cloud_microphys_mod
, only: tau_r2g, tau_smlt, tau_i2s, tau_v2l, tau_l2v, tau_imlt, tau_l2r
63 use gfdl_cloud_microphys_mod
, only: rad_rain, rad_snow, rad_graupel, dw_ocean, dw_land, tintqs
75 real,
parameter ::
cp_vap = 4.0 * rvgas
76 real,
parameter ::
cv_air = cp_air - rdgas
77 real,
parameter ::
cv_vap = 3.0 * rvgas
90 real,
parameter ::
c_ice = 1972.0
91 real,
parameter ::
c_liq = 4185.5
96 real,
parameter ::
tice = 273.16
103 real (kind = r_grid),
parameter ::
e00 = 611.21
108 real,
parameter ::
lat2 = (hlv + hlf) ** 2
122 subroutine fv_sat_adj (mdt, zvir, is, ie, js, je, ng, hydrostatic, consv_te, te0,&
128 ql, qi, qr, qs, qg, hs, dpln, delz, pt, dp, q_con, cappa, &
129 area, dtdt, out_dt, last_step, do_qa, qa)
133 integer,
intent (in) :: is, ie, js, je, ng
135 logical,
intent (in) :: hydrostatic, consv_te, out_dt, last_step, do_qa
137 real,
intent (in) :: zvir, mdt
139 real,
intent (in),
dimension (is - ng:ie + ng, js - ng:je + ng) :: dp, hs
140 real,
intent (in),
dimension (is:ie, js:je) :: dpln, delz
144 integer,
intent(in) :: km
145 real,
intent (inout),
dimension (is - ng:ie + ng, js - ng:je + ng,km,*) :: qvi
147 real,
intent (inout),
dimension (is - ng:ie + ng, js - ng:je + ng) :: qv
149 real,
intent (inout),
dimension (is - ng:ie + ng, js - ng:je + ng) :: pt, ql, qi, qr, qs, qg
150 real,
intent (inout),
dimension (is - ng:, js - ng:) :: q_con, cappa
151 real,
intent (inout),
dimension (is:ie, js:je) :: dtdt
153 real,
intent (out),
dimension (is - ng:ie + ng, js - ng:je + ng) :: qa, te0
155 real (kind = r_grid),
intent (in),
dimension (is - ng:ie + ng, js - ng:je + ng) :: area
158 real,
dimension (is - ng:ie + ng, js - ng:je + ng) :: qv
160 real,
dimension (is:ie) :: wqsat, dq2dt, qpz, cvm, t0, pt1, qstar
161 real,
dimension (is:ie) :: icp2, lcp2, tcp2, tcp3
162 real,
dimension (is:ie) :: den, q_liq, q_sol, q_cond, src, sink, hvar
163 real,
dimension (is:ie) :: mc_air, lhl, lhi
166 real :: tc, qsi, dqsdt, dq, dq0, pidep, qi_crt, tmp, dtmp
167 real :: tin, rqi, q_plus, q_minus
168 real :: sdt, dt_bigg, adj_fac
169 real :: fac_smlt, fac_r2g, fac_i2s, fac_imlt, fac_l2r, fac_v2l, fac_l2v
170 real :: factor, qim, tice0, c_air, c_vap, dw
175 qv(:,:) = qvi(:,:,1,1)
186 fac_i2s = 1. - exp(- mdt / tau_i2s)
187 fac_v2l = 1. - exp(- sdt / tau_v2l)
188 fac_r2g = 1. - exp(- mdt / tau_r2g)
189 fac_l2r = 1. - exp(- mdt / tau_l2r)
191 fac_l2v = 1. - exp(- sdt / tau_l2v)
192 fac_l2v = min(sat_adj0, fac_l2v)
194 fac_imlt = 1. - exp(- sdt / tau_imlt)
195 fac_smlt = 1. - exp(- mdt / tau_smlt)
201 if (hydrostatic)
then 216 q_liq(i) = ql(i, j) + qr(i, j)
217 q_sol(i) = qi(i, j) + qs(i, j) + qg(i, j)
218 qpz(i) = q_liq(i) + q_sol(i)
223 pt1(i) = pt(i, j) / ((1 + zvir * qv(i, j)) * (1 - qpz(i)))
225 pt1(i) = pt(i, j) / (1 + zvir * qv(i, j))
229 qpz(i) = qpz(i) + qv(i, j)
236 if (hydrostatic)
then 238 den(i) = dp(i, j) / (dpln(i, j) * rdgas * pt(i, j))
242 den(i) = - dp(i, j) / (grav * delz(i, j))
252 if (hydrostatic)
then 258 mc_air(i) = (1. - qpz(i)) * c_air
259 cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) *
c_liq + q_sol(i) *
c_ice 261 icp2(i) = lhi(i) / cvm(i)
269 if (hydrostatic)
then 274 te0(i, j) = - c_air * t0(i)
279 te0(i, j) = - cvm(i) * t0(i)
284 te0(i, j) = - c_air * t0(i)
295 if (qi(i, j) < 0.)
then 296 qs(i, j) = qs(i, j) + qi(i, j)
306 if (qi(i, j) > 1.e-8 .and. pt1(i) >
tice)
then 307 sink(i) = min(qi(i, j), fac_imlt * (pt1(i) -
tice) / icp2(i))
308 qi(i, j) = qi(i, j) - sink(i)
314 ql(i, j) = ql(i, j) + sink(i)
315 q_liq(i) = q_liq(i) + sink(i)
316 q_sol(i) = q_sol(i) - sink(i)
317 cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) *
c_liq + q_sol(i) *
c_ice 318 pt1(i) = pt1(i) - sink(i) * lhi(i) / cvm(i)
328 icp2(i) = lhi(i) / cvm(i)
336 if (qs(i, j) < 0.)
then 337 qg(i, j) = qg(i, j) + qs(i, j)
339 elseif (qg(i, j) < 0.)
then 340 tmp = min(- qg(i, j), max(0., qs(i, j)))
341 qg(i, j) = qg(i, j) + tmp
342 qs(i, j) = qs(i, j) - tmp
353 if (ql(i, j) < 0.)
then 354 tmp = min(- ql(i, j), max(0., qr(i, j)))
355 ql(i, j) = ql(i, j) + tmp
356 qr(i, j) = qr(i, j) - tmp
357 elseif (qr(i, j) < 0.)
then 358 tmp = min(- qr(i, j), max(0., ql(i, j)))
359 ql(i, j) = ql(i, j) - tmp
360 qr(i, j) = qr(i, j) + tmp
369 dtmp =
tice - 48. - pt1(i)
370 if (ql(i, j) > 0. .and. dtmp > 0.)
then 371 sink(i) = min(ql(i, j), dtmp / icp2(i))
372 ql(i, j) = ql(i, j) - sink(i)
373 qi(i, j) = qi(i, j) + sink(i)
374 q_liq(i) = q_liq(i) - sink(i)
375 q_sol(i) = q_sol(i) + sink(i)
376 cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) *
c_liq + q_sol(i) *
c_ice 377 pt1(i) = pt1(i) + sink(i) * lhi(i) / cvm(i)
388 lcp2(i) = lhl(i) / cvm(i)
389 icp2(i) = lhi(i) / cvm(i)
390 tcp3(i) = lcp2(i) + icp2(i) * min(1., dim(
tice, pt1(i)) / 48.)
397 call wqs2_vect (is, ie, pt1, den, wqsat, dq2dt)
401 dq0 = (qv(i, j) - wqsat(i)) / (1. + tcp3(i) * dq2dt(i))
403 src(i) = min(adj_fac * dq0, max(ql_gen - ql(i, j), fac_v2l * dq0))
409 factor = - min(1., fac_l2v * 10. * (1. - qv(i, j) / wqsat(i)))
410 src(i) = - min(ql(i, j), factor * dq0)
412 qv(i, j) = qv(i, j) - src(i)
414 qvi(i,j,1,1) = qv(i, j)
416 ql(i, j) = ql(i, j) + src(i)
417 q_liq(i) = q_liq(i) + src(i)
418 cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) *
c_liq + q_sol(i) *
c_ice 419 pt1(i) = pt1(i) + src(i) * lhl(i) / cvm(i)
429 lcp2(i) = lhl(i) / cvm(i)
430 icp2(i) = lhi(i) / cvm(i)
431 tcp3(i) = lcp2(i) + icp2(i) * min(1., dim(
tice, pt1(i)) / 48.)
442 call wqs2_vect (is, ie, pt1, den, wqsat, dq2dt)
445 dq0 = (qv(i, j) - wqsat(i)) / (1. + tcp3(i) * dq2dt(i))
452 factor = - min(1., fac_l2v * 10. * (1. - qv(i, j) / wqsat(i)))
453 src(i) = - min(ql(i, j), factor * dq0)
456 qv(i, j) = qv(i, j) - src(i)
458 qvi(i,j,1,1) = qv(i,j)
460 ql(i, j) = ql(i, j) + src(i)
461 q_liq(i) = q_liq(i) + src(i)
462 cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) *
c_liq + q_sol(i) *
c_ice 463 pt1(i) = pt1(i) + src(i) * lhl(i) / cvm(i)
473 lcp2(i) = lhl(i) / cvm(i)
474 icp2(i) = lhi(i) / cvm(i)
484 dtmp =
t_wfr - pt1(i)
485 if (ql(i, j) > 0. .and. dtmp > 0.)
then 486 sink(i) = min(ql(i, j), ql(i, j) * dtmp * 0.125, dtmp / icp2(i))
487 ql(i, j) = ql(i, j) - sink(i)
488 qi(i, j) = qi(i, j) + sink(i)
489 q_liq(i) = q_liq(i) - sink(i)
490 q_sol(i) = q_sol(i) + sink(i)
491 cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) *
c_liq + q_sol(i) *
c_ice 492 pt1(i) = pt1(i) + sink(i) * lhi(i) / cvm(i)
502 icp2(i) = lhi(i) / cvm(i)
511 if (ql(i, j) > 0.0 .and. tc > 0.)
then 512 sink(i) = 3.3333e-10 * dt_bigg * (exp(0.66 * tc) - 1.) * den(i) * ql(i, j) ** 2
513 sink(i) = min(ql(i, j), tc / icp2(i), sink(i))
514 ql(i, j) = ql(i, j) - sink(i)
515 qi(i, j) = qi(i, j) + sink(i)
516 q_liq(i) = q_liq(i) - sink(i)
517 q_sol(i) = q_sol(i) + sink(i)
518 cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) *
c_liq + q_sol(i) *
c_ice 519 pt1(i) = pt1(i) + sink(i) * lhi(i) / cvm(i)
529 icp2(i) = lhi(i) / cvm(i)
537 dtmp = (
tice - 0.1) - pt1(i)
538 if (qr(i, j) > 1.e-7 .and. dtmp > 0.)
then 539 tmp = min(1., (dtmp * 0.025) ** 2) * qr(i, j)
540 sink(i) = min(tmp, fac_r2g * dtmp / icp2(i))
541 qr(i, j) = qr(i, j) - sink(i)
542 qg(i, j) = qg(i, j) + sink(i)
543 q_liq(i) = q_liq(i) - sink(i)
544 q_sol(i) = q_sol(i) + sink(i)
545 cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) *
c_liq + q_sol(i) *
c_ice 546 pt1(i) = pt1(i) + sink(i) * lhi(i) / cvm(i)
556 icp2(i) = lhi(i) / cvm(i)
564 dtmp = pt1(i) - (
tice + 0.1)
565 if (qs(i, j) > 1.e-7 .and. dtmp > 0.)
then 566 tmp = min(1., (dtmp * 0.1) ** 2) * qs(i, j)
567 sink(i) = min(tmp, fac_smlt * dtmp / icp2(i))
568 tmp = min(sink(i), dim(qs_mlt, ql(i, j)))
569 qs(i, j) = qs(i, j) - sink(i)
570 ql(i, j) = ql(i, j) + tmp
571 qr(i, j) = qr(i, j) + sink(i) - tmp
573 q_liq(i) = q_liq(i) + sink(i)
574 q_sol(i) = q_sol(i) - sink(i)
575 cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) *
c_liq + q_sol(i) *
c_ice 576 pt1(i) = pt1(i) - sink(i) * lhi(i) / cvm(i)
585 if (ql(i, j) > ql0_max)
then 586 sink(i) = fac_l2r * (ql(i, j) - ql0_max)
587 qr(i, j) = qr(i, j) + sink(i)
588 ql(i, j) = ql(i, j) - sink(i)
599 lcp2(i) = lhl(i) / cvm(i)
600 icp2(i) = lhi(i) / cvm(i)
601 tcp2(i) = lcp2(i) + icp2(i)
610 if (pt1(i) < t_sub)
then 611 src(i) = dim(qv(i, j), 1.e-6)
612 elseif (pt1(i) < tice0)
then 613 qsi =
iqs2(pt1(i), den(i), dqsdt)
615 sink(i) = adj_fac * dq / (1. + tcp2(i) * dqsdt)
616 if (qi(i, j) > 1.e-8)
then 617 pidep = sdt * dq * 349138.78 * exp(0.875 * log(qi(i, j) * den(i))) &
618 / (qsi * den(i) *
lat2 / (0.0243 * rvgas * pt1(i) ** 2) + 4.42478e4)
624 qi_crt = qi_gen * min(qi_lim, 0.1 * tmp) / den(i)
625 src(i) = min(sink(i), max(qi_crt - qi(i, j), pidep), tmp / tcp2(i))
627 pidep = pidep * min(1., dim(pt1(i), t_sub) * 0.2)
628 src(i) = max(pidep, sink(i), - qi(i, j))
631 qv(i, j) = qv(i, j) - src(i)
633 qvi(i,j,1,1) = qv(i,j)
635 qi(i, j) = qi(i, j) + src(i)
636 q_sol(i) = q_sol(i) + src(i)
637 cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) *
c_liq + q_sol(i) *
c_ice 638 pt1(i) = pt1(i) + src(i) * (lhl(i) + lhi(i)) / cvm(i)
647 q_con(i, j) = q_liq(i) + q_sol(i)
651 tmp = 1. + zvir * qv(i, j)
652 pt(i, j) = pt1(i) * tmp * (1. - q_con(i, j))
655 cappa(i, j) = tmp / (tmp + cvm(i))
658 q_con(i, j) = q_liq(i) + q_sol(i)
659 pt(i, j) = pt1(i) *
virq_qpz(qvi(i,j,1,1:
num_gas),q_con(i,j)) * (1. - q_con(i,j))
661 pt(i, j) = pt1(i) * (1. + zvir * qv(i, j))
671 if (qg(i, j) < 0.)
then 672 tmp = min(- qg(i, j), max(0., qi(i, j)))
673 qg(i, j) = qg(i, j) + tmp
674 qi(i, j) = qi(i, j) - tmp
683 qim = qi0_max / den(i)
684 if (qi(i, j) > qim)
then 685 sink(i) = fac_i2s * (qi(i, j) - qim)
686 qi(i, j) = qi(i, j) - sink(i)
687 qs(i, j) = qs(i, j) + sink(i)
693 dtdt(i, j) = dtdt(i, j) + pt1(i) - t0(i)
703 if (hydrostatic)
then 707 te0(i, j) = dp(i, j) * (te0(i, j) + c_air * pt1(i))
710 te0(i, j) = dp(i, j) * (te0(i, j) + cvm(i) * pt1(i))
715 te0(i, j) = dp(i, j) * (te0(i, j) + c_air * pt1(i))
728 cvm(i) = mc_air(i) + (qv(i, j) + q_liq(i) + q_sol(i)) * c_vap
729 lcp2(i) = lhl(i) / cvm(i)
730 icp2(i) = lhi(i) / cvm(i)
737 if (do_qa .and. last_step)
then 744 if (rad_graupel)
then 746 q_sol(i) = qi(i, j) + qs(i, j) + qg(i, j)
750 q_sol(i) = qi(i, j) + qs(i, j)
760 q_liq(i) = ql(i, j) + qr(i, j)
768 q_cond(i) = q_sol(i) + q_liq(i)
780 tin = pt1(i) - (lcp2(i) * q_cond(i) + icp2(i) * q_sol(i))
789 if (tin <=
t_wfr)
then 791 qstar(i) =
iqs1(tin, den(i))
792 elseif (tin >=
tice)
then 794 qstar(i) =
wqs1(tin, den(i))
797 qsi =
iqs1(tin, den(i))
798 qsw =
wqs1(tin, den(i))
799 if (q_cond(i) > 1.e-6)
then 800 rqi = q_sol(i) / q_cond(i)
805 qstar(i) = rqi * qsi + (1. - rqi) * qsw
809 dw = dw_ocean + (dw_land - dw_ocean) * min(1., abs(hs(i, j)) / (10. * grav))
811 hvar(i) = min(0.2, max(0.01, dw * sqrt(sqrt(area(i, j)) / 100.e3)))
819 rh = qpz(i) / qstar(i)
827 if (rh > 0.75 .and. qpz(i) > 1.e-8)
then 828 dq = hvar(i) * qpz(i)
830 q_minus = qpz(i) - dq
831 if (icloud_f == 2)
then 832 if (qpz(i) > qstar(i))
then 834 elseif (qstar(i) < q_plus .and. q_cond(i) > 1.e-8)
then 835 qa(i, j) = ((q_plus - qstar(i)) / dq) ** 2
836 qa(i, j) = min(1., qa(i, j))
841 if (qstar(i) < q_minus)
then 844 if (qstar(i) < q_plus)
then 845 if (icloud_f == 0)
then 846 qa(i, j) = (q_plus - qstar(i)) / (dq + dq)
848 qa(i, j) = (q_plus - qstar(i)) / (2. * dq * (1. - q_cond(i)))
854 if (q_cond(i) > 1.e-8)
then 855 qa(i, j) = max(cld_min, qa(i, j))
857 qa(i, j) = min(1., qa(i, j))
876 real function wqs1 (ta, den)
883 real,
intent (in) :: ta, den
885 real :: es, ap1, tmin
890 ap1 = 10. * dim(ta, tmin) + 1.
891 ap1 = min(2621., ap1)
894 wqs1 = es / (rvgas * ta * den)
902 real function iqs1 (ta, den)
909 real,
intent (in) :: ta, den
911 real :: es, ap1, tmin
916 ap1 = 10. * dim(ta, tmin) + 1.
917 ap1 = min(2621., ap1)
920 iqs1 = es / (rvgas * ta * den)
928 real function wqs2 (ta, den, dqdt)
935 real,
intent (in) :: ta, den
937 real,
intent (out) :: dqdt
939 real :: es, ap1, tmin
944 ap1 = 10. * dim(ta, tmin) + 1.
945 ap1 = min(2621., ap1)
948 wqs2 = es / (rvgas * ta * den)
951 dqdt = 10. * (
desw(it) + (ap1 - it) * (
desw(it + 1) -
desw(it))) / (rvgas * ta * den)
960 subroutine wqs2_vect (is, ie, ta, den, wqsat, dqdt)
967 integer,
intent (in) :: is, ie
969 real,
intent (in),
dimension (is:ie) :: ta, den
971 real,
intent (out),
dimension (is:ie) :: wqsat, dqdt
973 real :: es, ap1, tmin
980 ap1 = 10. * dim(ta(i), tmin) + 1.
981 ap1 = min(2621., ap1)
984 wqsat(i) = es / (rvgas * ta(i) * den(i))
987 dqdt(i) = 10. * (
desw(it) + (ap1 - it) * (
desw(it + 1) -
desw(it))) / (rvgas * ta(i) * den(i))
996 real function iqs2 (ta, den, dqdt)
1003 real,
intent (in) :: ta, den
1005 real,
intent (out) :: dqdt
1007 real :: es, ap1, tmin
1012 ap1 = 10. * dim(ta, tmin) + 1.
1013 ap1 = min(2621., ap1)
1016 iqs2 = es / (rvgas * ta * den)
1019 dqdt = 10. * (
des2(it) + (ap1 - it) * (
des2(it + 1) -
des2(it))) / (rvgas * ta * den)
1032 integer,
intent (in) :: kmp
1034 integer,
parameter :: length = 2621
1040 if (is_master())
write (*, *)
'top layer for gfdl_mp = ', kmp
1044 allocate (
table(length))
1045 allocate (
table2(length))
1046 allocate (
tablew(length))
1047 allocate (
des2(length))
1048 allocate (
desw(length))
1054 do i = 1, length - 1
1074 integer,
intent (in) :: n
1076 real (kind = r_grid) :: delt = 0.1
1077 real (kind = r_grid) :: tmin, tem, esh20
1078 real (kind = r_grid) :: wice, wh2o, fac0, fac1, fac2
1079 real (kind = r_grid) :: esupc (200)
1090 tem = tmin + delt * real(i - 1)
1093 fac2 = (
d2ice * log(tem /
tice) + fac1) / rvgas
1102 tem = 253.16 + delt * real(i - 1)
1105 fac2 = (
dc_vap * log(tem /
tice) + fac1) / rvgas
1106 esh20 =
e00 * exp(fac2)
1110 table(i + 1400) = esh20
1119 tem = 253.16 + delt * real(i - 1)
1120 wice = 0.05 * (
tice - tem)
1121 wh2o = 0.05 * (tem - 253.16)
1122 table(i + 1400) = wice *
table(i + 1400) + wh2o * esupc(i)
1136 integer,
intent (in) :: n
1138 real (kind = r_grid) :: delt = 0.1
1139 real (kind = r_grid) :: tmin, tem, fac0, fac1, fac2
1150 tem = tmin + delt * real(i - 1)
1153 fac2 = (
dc_vap * log(tem /
tice) + fac1) / rvgas
1168 integer,
intent (in) :: n
1170 real (kind = r_grid) :: delt = 0.1
1171 real (kind = r_grid) :: tmin, tem0, tem1, fac0, fac1, fac2
1173 integer :: i, i0, i1
1178 tem0 = tmin + delt * real(i - 1)
1179 fac0 = (tem0 -
tice) / (tem0 *
tice)
1185 fac2 = (
d2ice * log(tem0 /
tice) + fac1) / rvgas
1191 fac2 = (
dc_vap * log(tem0 /
tice) + fac1) / rvgas
real, dimension(:), allocatable des2
real, parameter cp_vap
1846.0, heat capacity of water vapor at constant pressure
real, dimension(:), allocatable table
The module 'fv_mp_mod' is a single program multiple data (SPMD) parallel decompostion/communication m...
real, dimension(:), allocatable tablew
real, parameter tice
freezing temperature
subroutine, public fv_sat_adj(mdt, zvir, is, ie, js, je, ng, hydrostatic, consv_te, te0, qv, ql, qi, qr, qs, qg, hs, dpln, delz, pt, dp, q_con, cappa, area, dtdt, out_dt, last_step, do_qa, qa)
The subroutine 'fv_sat_adj' performs the fast processes in the GFDL microphysics. ...
real, dimension(:), allocatable table2
real, parameter cv_air
717.55, heat capacity of dry air at constant volume
The module 'multi_gases' peforms multi constitutents computations.
real, parameter cv_vap
1384.5, heat capacity of water vapor at constant volume
real, parameter c_ice
gfdl: heat capacity of ice at - 15 deg c
pure real function, public vicvqd_qpz(q, qpz)
real d0_vap
the same as dc_vap, except that cp_vap can be cp_vap or cv_vap
pure real function, public virq_qpz(q, qpz)
integer, parameter, public r_grid
real function iqs1(ta, den)
the function 'wqs1' computes the saturated specific humidity for table iii
real(kind=r_grid), parameter d2ice
real, dimension(:), allocatable desw
The module 'fv_arrays' contains the 'fv_atmos_type' and associated datatypes.
real lv00
the same as lv0, except that cp_vap can be cp_vap or cv_vap
real, parameter c_liq
gfdl: heat capacity of liquid at 15 deg c
real function iqs2(ta, den, dqdt)
The function 'iqs2' computes the gradient of saturated specific humidity for table iii...
real, parameter dc_ice
2213.5, isobaric heating / colling
real, parameter lat2
used in bigg mechanism
real, parameter t_wfr
homogeneous freezing temperature
real, parameter lv0
3.13905782e6, evaporation latent heat coefficient at 0 deg k
pure real function, public vicpqd_qpz(q, qpz)
real(kind=r_grid), parameter li2
2.86799816e6, sublimation latent heat coefficient at 0 deg k
real function wqs2(ta, den, dqdt)
The function 'wqs2'computes the gradient of saturated specific humidity for table ii...
real function wqs1(ta, den)
the function 'wqs1' computes the saturated specific humidity for table ii
subroutine, public qs_init(kmp)
The subroutine 'qs_init' initializes lookup tables for the saturation mixing ratio.
real(kind=r_grid), parameter e00
ifs: saturation vapor pressure at 0 deg c
subroutine wqs2_vect(is, ie, ta, den, wqsat, dqdt)
The function wqs2_vect computes the gradient of saturated specific humidity for table ii...
The module 'fv_cmp' implements the fast procesesses in the GFDL microphysics >