FV3DYCORE  Version 2.0.0
fv_cmp.F90
Go to the documentation of this file.
1 
2 !***********************************************************************
3 !* GNU Lesser General Public License
4 !*
5 !* This file is part of the GFDL Cloud Microphysics.
6 !*
7 !* The GFDL Cloud Microphysics is free software: you can
8 !8 redistribute it and/or modify it under the terms of the
9 !* GNU Lesser General Public License as published by the
10 !* Free Software Foundation, either version 3 of the License, or
11 !* (at your option) any later version.
12 !*
13 !* The GFDL Cloud Microphysics is distributed in the hope it will be
14 !* useful, but WITHOUT ANYWARRANTY; without even the implied warranty
15 !* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
16 !* See the GNU General Public License for more details.
17 !*
18 !* You should have received a copy of the GNU Lesser General Public
19 !* License along with the GFDL Cloud Microphysics.
20 !* If not, see <http://www.gnu.org/licenses/>.
21 !***********************************************************************
22 
26 ! Fast saturation adjustment is part of the gfdl cloud microphysics
27 ! =======================================================================
28 
29 module fv_cmp_mod
30 ! Modules Included:
31 ! <table>
32 ! <tr>
33 ! <th>Module Name</th>
34 ! <th>Functions Included</th>
35 ! </tr>
36 ! <tr>
37 ! <td>constants_mod</td>
38 ! <td>rvgas, rdgas, grav, hlv, hlf, cp_air</td>
39 ! </tr>
40 ! <tr>
41 ! <td>fv_arrays_mod</td>
42 ! <td> r_grid</td>
43 ! </tr>
44 ! <tr>
45 ! <tr>
46 ! <td>fv_mp_mod</td>
47 ! <td>is_master</td>
48 ! </tr>
49 ! <tr>
50 ! <td>gfdl_cloud_microphys_mod</td>
51 ! <td>ql_gen, qi_gen, qi0_max, ql_mlt, ql0_max, qi_lim, qs_mlt,
52 ! tau_r2g, tau_smlt, tau_i2s, tau_v2l, tau_l2v, tau_imlt, tau_l2r,
53 ! rad_rain, rad_snow, rad_graupel, dw_ocean, dw_land, tintqs</td>
54 ! </tr>
55 ! </table>
56 
57  use constants_mod, only: rvgas, rdgas, grav, hlv, hlf, cp_air
58  use fv_mp_mod, only: is_master
59  use fv_arrays_mod, only: r_grid
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
64 #ifdef MULTI_GASES
66 #endif
67 
68  implicit none
69 
70  private
71 
72  public fv_sat_adj, qs_init
73 
74  ! real, parameter :: cp_air = cp_air ! 1004.6, heat capacity of dry air at constant pressure, come from constants_mod
75  real, parameter :: cp_vap = 4.0 * rvgas
76  real, parameter :: cv_air = cp_air - rdgas
77  real, parameter :: cv_vap = 3.0 * rvgas
78 
79  ! http: / / www.engineeringtoolbox.com / ice - thermal - properties - d_576.html
80  ! c_ice = 2050.0 at 0 deg c
81  ! c_ice = 1972.0 at - 15 deg c
82  ! c_ice = 1818.0 at - 40 deg c
83  ! http: / / www.engineeringtoolbox.com / water - thermal - properties - d_162.html
84  ! c_liq = 4205.0 at 4 deg c
85  ! c_liq = 4185.5 at 15 deg c
86  ! c_liq = 4178.0 at 30 deg c
87 
88  ! real, parameter :: c_ice = 2106.0 ! ifs: heat capacity of ice at 0 deg c
89  ! real, parameter :: c_liq = 4218.0 ! ifs: heat capacity of liquid at 0 deg c
90  real, parameter :: c_ice = 1972.0
91  real, parameter :: c_liq = 4185.5
92 
93  real, parameter :: dc_vap = cp_vap - c_liq
94  real, parameter :: dc_ice = c_liq - c_ice
95 
96  real, parameter :: tice = 273.16
97  real, parameter :: t_wfr = tice - 40.
98 
99  real, parameter :: lv0 = hlv - dc_vap * tice
100  real, parameter :: li00 = hlf - dc_ice * tice
101 
102  ! real (kind = r_grid), parameter :: e00 = 610.71 ! gfdl: saturation vapor pressure at 0 deg c
103  real (kind = r_grid), parameter :: e00 = 611.21
104 
105  real (kind = r_grid), parameter :: d2ice = dc_vap + dc_ice
106  real (kind = r_grid), parameter :: li2 = lv0 + li00
107 
108  real, parameter :: lat2 = (hlv + hlf) ** 2
109 
110  real :: d0_vap
111  real :: lv00
112 
113  real, allocatable :: table (:), table2 (:), tablew (:), des2 (:), desw (:)
114 
115  logical :: mp_initialized = .false.
116 
117 contains
118 
122 subroutine fv_sat_adj (mdt, zvir, is, ie, js, je, ng, hydrostatic, consv_te, te0,&
123 #ifdef MULTI_GASES
124  km, qvi, &
125 #else
126  qv, &
127 #endif
128  ql, qi, qr, qs, qg, hs, dpln, delz, pt, dp, q_con, cappa, &
129  area, dtdt, out_dt, last_step, do_qa, qa)
130 
131  implicit none
132 
133  integer, intent (in) :: is, ie, js, je, ng
134 
135  logical, intent (in) :: hydrostatic, consv_te, out_dt, last_step, do_qa
136 
137  real, intent (in) :: zvir, mdt ! remapping time step
138 
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
141 
142 
143 #ifdef MULTI_GASES
144  integer, intent(in) :: km
145  real, intent (inout), dimension (is - ng:ie + ng, js - ng:je + ng,km,*) :: qvi
146 #else
147  real, intent (inout), dimension (is - ng:ie + ng, js - ng:je + ng) :: qv
148 #endif
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
152 
153  real, intent (out), dimension (is - ng:ie + ng, js - ng:je + ng) :: qa, te0
154 
155  real (kind = r_grid), intent (in), dimension (is - ng:ie + ng, js - ng:je + ng) :: area
156 
157 #ifdef MULTI_GASES
158  real, dimension (is - ng:ie + ng, js - ng:je + ng) :: qv
159 #endif
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
164 
165  real :: qsw, rh
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
171 
172  integer :: i, j
173 
174 #ifdef MULTI_GASES
175  qv(:,:) = qvi(:,:,1,1)
176 #endif
177  sdt = 0.5 * mdt ! half remapping time step
178  dt_bigg = mdt ! bigg mechinism time step
179 
180  tice0 = tice - 0.01 ! 273.15, standard freezing temperature
181 
182  ! -----------------------------------------------------------------------
183  ! define conversion scalar / factor
184  ! -----------------------------------------------------------------------
185 
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)
190 
191  fac_l2v = 1. - exp(- sdt / tau_l2v)
192  fac_l2v = min(sat_adj0, fac_l2v)
193 
194  fac_imlt = 1. - exp(- sdt / tau_imlt)
195  fac_smlt = 1. - exp(- mdt / tau_smlt)
196 
197  ! -----------------------------------------------------------------------
198  ! define heat capacity of dry air and water vapor based on hydrostatical property
199  ! -----------------------------------------------------------------------
200 
201  if (hydrostatic) then
202  c_air = cp_air
203  c_vap = cp_vap
204  else
205  c_air = cv_air
206  c_vap = cv_vap
207  endif
208  d0_vap = c_vap - c_liq
209  lv00 = hlv - d0_vap * tice
210  ! dc_vap = cp_vap - c_liq ! - 2339.5
211  ! d0_vap = cv_vap - c_liq ! - 2801.0
212 
213  do j = js, je ! start j loop
214 
215  do i = is, ie
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)
219 #ifdef MULTI_GASES
220  pt1(i) = pt(i, j) / virq_qpz(qvi(i,j,1,1:num_gas),qpz(i))
221 #else
222 #ifdef USE_COND
223  pt1(i) = pt(i, j) / ((1 + zvir * qv(i, j)) * (1 - qpz(i)))
224 #else
225  pt1(i) = pt(i, j) / (1 + zvir * qv(i, j))
226 #endif
227 #endif
228  t0(i) = pt1(i) ! true temperature
229  qpz(i) = qpz(i) + qv(i, j) ! total_wat conserved in this routine
230  enddo
231 
232  ! -----------------------------------------------------------------------
233  ! define air density based on hydrostatical property
234  ! -----------------------------------------------------------------------
235 
236  if (hydrostatic) then
237  do i = is, ie
238  den(i) = dp(i, j) / (dpln(i, j) * rdgas * pt(i, j))
239  enddo
240  else
241  do i = is, ie
242  den(i) = - dp(i, j) / (grav * delz(i, j)) ! moist_air density
243  enddo
244  endif
245 
246  ! -----------------------------------------------------------------------
247  ! define heat capacity and latend heat coefficient
248  ! -----------------------------------------------------------------------
249 
250  do i = is, ie
251 #ifdef MULTI_GASES
252  if (hydrostatic) then
253  c_air = cp_air * vicpqd_qpz(qvi(i,j,1,1:num_gas),qpz(i))
254  else
255  c_air = cv_air * vicvqd_qpz(qvi(i,j,1,1:num_gas),qpz(i))
256  endif
257 #endif
258  mc_air(i) = (1. - qpz(i)) * c_air ! constant
259  cvm(i) = mc_air(i) + qv(i, j) * c_vap + q_liq(i) * c_liq + q_sol(i) * c_ice
260  lhi(i) = li00 + dc_ice * pt1(i)
261  icp2(i) = lhi(i) / cvm(i)
262  enddo
263 
264  ! -----------------------------------------------------------------------
265  ! fix energy conservation
266  ! -----------------------------------------------------------------------
267 
268  if (consv_te) then
269  if (hydrostatic) then
270  do i = is, ie
271 #ifdef MULTI_GASES
272  c_air = cp_air * vicpqd_qpz(qvi(i,j,1,1:num_gas),qpz(i))
273 #endif
274  te0(i, j) = - c_air * t0(i)
275  enddo
276  else
277  do i = is, ie
278 #ifdef USE_COND
279  te0(i, j) = - cvm(i) * t0(i)
280 #else
281 #ifdef MULTI_GASES
282  c_air = cv_air * vicvqd_qpz(qvi(i,j,1,1:num_gas),qpz(i))
283 #endif
284  te0(i, j) = - c_air * t0(i)
285 #endif
286  enddo
287  endif
288  endif
289 
290  ! -----------------------------------------------------------------------
291  ! fix negative cloud ice with snow
292  ! -----------------------------------------------------------------------
293 
294  do i = is, ie
295  if (qi(i, j) < 0.) then
296  qs(i, j) = qs(i, j) + qi(i, j)
297  qi(i, j) = 0.
298  endif
299  enddo
300 
301  ! -----------------------------------------------------------------------
302  ! melting of cloud ice to cloud water and rain
303  ! -----------------------------------------------------------------------
304 
305  do i = is, ie
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)
309  ! sjl, may 17, 2017
310  ! tmp = min (sink (i), dim (ql_mlt, ql (i, j))) ! max ql amount
311  ! ql (i, j) = ql (i, j) + tmp
312  ! qr (i, j) = qr (i, j) + sink (i) - tmp
313  ! sjl, may 17, 2017
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)
319  endif
320  enddo
321 
322  ! -----------------------------------------------------------------------
323  ! update latend heat coefficient
324  ! -----------------------------------------------------------------------
325 
326  do i = is, ie
327  lhi(i) = li00 + dc_ice * pt1(i)
328  icp2(i) = lhi(i) / cvm(i)
329  enddo
330 
331  ! -----------------------------------------------------------------------
332  ! fix negative snow with graupel or graupel with available snow
333  ! -----------------------------------------------------------------------
334 
335  do i = is, ie
336  if (qs(i, j) < 0.) then
337  qg(i, j) = qg(i, j) + qs(i, j)
338  qs(i, j) = 0.
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
343  endif
344  enddo
345 
346  ! after this point cloud ice & snow are positive definite
347 
348  ! -----------------------------------------------------------------------
349  ! fix negative cloud water with rain or rain with available cloud water
350  ! -----------------------------------------------------------------------
351 
352  do i = is, ie
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
361  endif
362  enddo
363 
364  ! -----------------------------------------------------------------------
365  ! enforce complete freezing of cloud water to cloud ice below - 48 c
366  ! -----------------------------------------------------------------------
367 
368  do i = is, ie
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)
378  endif
379  enddo
380 
381  ! -----------------------------------------------------------------------
382  ! update latend heat coefficient
383  ! -----------------------------------------------------------------------
384 
385  do i = is, ie
386  lhl(i) = lv00 + d0_vap * pt1(i)
387  lhi(i) = li00 + dc_ice * pt1(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.)
391  enddo
392 
393  ! -----------------------------------------------------------------------
394  ! condensation / evaporation between water vapor and cloud water
395  ! -----------------------------------------------------------------------
396 
397  call wqs2_vect (is, ie, pt1, den, wqsat, dq2dt)
398 
399  adj_fac = sat_adj0
400  do i = is, ie
401  dq0 = (qv(i, j) - wqsat(i)) / (1. + tcp3(i) * dq2dt(i))
402  if (dq0 > 0.) then ! whole grid - box saturated
403  src(i) = min(adj_fac * dq0, max(ql_gen - ql(i, j), fac_v2l * dq0))
404  else ! evaporation of ql
405  ! sjl 20170703 added ql factor to prevent the situation of high ql and rh<1
406  ! factor = - min (1., fac_l2v * sqrt (max (0., ql (i, j)) / 1.e-5) * 10. * (1. - qv (i, j) / wqsat (i)))
407  ! factor = - fac_l2v
408  ! factor = - 1
409  factor = - min(1., fac_l2v * 10. * (1. - qv(i, j) / wqsat(i))) ! the rh dependent factor = 1 at 90%
410  src(i) = - min(ql(i, j), factor * dq0)
411  endif
412  qv(i, j) = qv(i, j) - src(i)
413 #ifdef MULTI_GASES
414  qvi(i,j,1,1) = qv(i, j)
415 #endif
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)
420  enddo
421 
422  ! -----------------------------------------------------------------------
423  ! update latend heat coefficient
424  ! -----------------------------------------------------------------------
425 
426  do i = is, ie
427  lhl(i) = lv00 + d0_vap * pt1(i)
428  lhi(i) = li00 + dc_ice * pt1(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.)
432  enddo
433 
434  if (last_step) then
435 
436  ! -----------------------------------------------------------------------
437  ! condensation / evaporation between water vapor and cloud water, last time step
438  ! enforce upper (no super_sat) & lower (critical rh) bounds
439  ! final iteration:
440  ! -----------------------------------------------------------------------
441 
442  call wqs2_vect (is, ie, pt1, den, wqsat, dq2dt)
443 
444  do i = is, ie
445  dq0 = (qv(i, j) - wqsat(i)) / (1. + tcp3(i) * dq2dt(i))
446  if (dq0 > 0.) then ! remove super - saturation, prevent super saturation over water
447  src(i) = dq0
448  else ! evaporation of ql
449  ! factor = - min (1., fac_l2v * sqrt (max (0., ql (i, j)) / 1.e-5) * 10. * (1. - qv (i, j) / wqsat (i))) ! the rh dependent factor = 1 at 90%
450  ! factor = - fac_l2v
451  ! factor = - 1
452  factor = - min(1., fac_l2v * 10. * (1. - qv(i, j) / wqsat(i))) ! the rh dependent factor = 1 at 90%
453  src(i) = - min(ql(i, j), factor * dq0)
454  endif
455  adj_fac = 1.
456  qv(i, j) = qv(i, j) - src(i)
457 #ifdef MULTI_GASES
458  qvi(i,j,1,1) = qv(i,j)
459 #endif
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)
464  enddo
465 
466  ! -----------------------------------------------------------------------
467  ! update latend heat coefficient
468  ! -----------------------------------------------------------------------
469 
470  do i = is, ie
471  lhl(i) = lv00 + d0_vap * pt1(i)
472  lhi(i) = li00 + dc_ice * pt1(i)
473  lcp2(i) = lhl(i) / cvm(i)
474  icp2(i) = lhi(i) / cvm(i)
475  enddo
476 
477  endif
478 
479  ! -----------------------------------------------------------------------
480  ! homogeneous freezing of cloud water to cloud ice
481  ! -----------------------------------------------------------------------
482 
483  do i = is, ie
484  dtmp = t_wfr - pt1(i) ! [ - 40, - 48]
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)
493  endif
494  enddo
495 
496  ! -----------------------------------------------------------------------
497  ! update latend heat coefficient
498  ! -----------------------------------------------------------------------
499 
500  do i = is, ie
501  lhi(i) = li00 + dc_ice * pt1(i)
502  icp2(i) = lhi(i) / cvm(i)
503  enddo
504 
505  ! -----------------------------------------------------------------------
506  ! bigg mechanism (heterogeneous freezing of cloud water to cloud ice)
507  ! -----------------------------------------------------------------------
508 
509  do i = is, ie
510  tc = tice0 - pt1(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)
520  endif
521  enddo
522 
523  ! -----------------------------------------------------------------------
524  ! update latend heat coefficient
525  ! -----------------------------------------------------------------------
526 
527  do i = is, ie
528  lhi(i) = li00 + dc_ice * pt1(i)
529  icp2(i) = lhi(i) / cvm(i)
530  enddo
531 
532  ! -----------------------------------------------------------------------
533  ! freezing of rain to graupel
534  ! -----------------------------------------------------------------------
535 
536  do i = is, ie
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) ! no limit on freezing below - 40 deg c
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)
547  endif
548  enddo
549 
550  ! -----------------------------------------------------------------------
551  ! update latend heat coefficient
552  ! -----------------------------------------------------------------------
553 
554  do i = is, ie
555  lhi(i) = li00 + dc_ice * pt1(i)
556  icp2(i) = lhi(i) / cvm(i)
557  enddo
558 
559  ! -----------------------------------------------------------------------
560  ! melting of snow to rain or cloud water
561  ! -----------------------------------------------------------------------
562 
563  do i = is, ie
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) ! no limter on melting above 10 deg c
567  sink(i) = min(tmp, fac_smlt * dtmp / icp2(i))
568  tmp = min(sink(i), dim(qs_mlt, ql(i, j))) ! max ql due to snow melt
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
572  ! qr (i, j) = qr (i, j) + sink (i)
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)
577  endif
578  enddo
579 
580  ! -----------------------------------------------------------------------
581  ! autoconversion from cloud water to rain
582  ! -----------------------------------------------------------------------
583 
584  do i = is, ie
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)
589  endif
590  enddo
591 
592  ! -----------------------------------------------------------------------
593  ! update latend heat coefficient
594  ! -----------------------------------------------------------------------
595 
596  do i = is, ie
597  lhi(i) = li00 + dc_ice * pt1(i)
598  lhl(i) = lv00 + d0_vap * pt1(i)
599  lcp2(i) = lhl(i) / cvm(i)
600  icp2(i) = lhi(i) / cvm(i)
601  tcp2(i) = lcp2(i) + icp2(i)
602  enddo
603 
604  ! -----------------------------------------------------------------------
605  ! sublimation / deposition between water vapor and cloud ice
606  ! -----------------------------------------------------------------------
607 
608  do i = is, ie
609  src(i) = 0.
610  if (pt1(i) < t_sub) then ! too cold to be accurate; freeze qv as a fix
611  src(i) = dim(qv(i, j), 1.e-6)
612  elseif (pt1(i) < tice0) then
613  qsi = iqs2(pt1(i), den(i), dqsdt)
614  dq = qv(i, j) - qsi
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)
619  else
620  pidep = 0.
621  endif
622  if (dq > 0.) then ! vapor - > ice
623  tmp = tice - pt1(i)
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))
626  else
627  pidep = pidep * min(1., dim(pt1(i), t_sub) * 0.2)
628  src(i) = max(pidep, sink(i), - qi(i, j))
629  endif
630  endif
631  qv(i, j) = qv(i, j) - src(i)
632 #ifdef MULTI_GASES
633  qvi(i,j,1,1) = qv(i,j)
634 #endif
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)
639  enddo
640 
641  ! -----------------------------------------------------------------------
642  ! virtual temp updated
643  ! -----------------------------------------------------------------------
644 
645  do i = is, ie
646 #ifdef USE_COND
647  q_con(i, j) = q_liq(i) + q_sol(i)
648 #ifdef MULTI_GASES
649  pt(i, j) = pt1(i) * virq_qpz(qvi(i,j,1,1:num_gas),q_con(i,j))
650 #else
651  tmp = 1. + zvir * qv(i, j)
652  pt(i, j) = pt1(i) * tmp * (1. - q_con(i, j))
653 #endif
654  tmp = rdgas * tmp
655  cappa(i, j) = tmp / (tmp + cvm(i))
656 #else
657 #ifdef MULTI_GASES
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))
660 #else
661  pt(i, j) = pt1(i) * (1. + zvir * qv(i, j))
662 #endif
663 #endif
664  enddo
665 
666  ! -----------------------------------------------------------------------
667  ! fix negative graupel with available cloud ice
668  ! -----------------------------------------------------------------------
669 
670  do i = is, ie
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
675  endif
676  enddo
677 
678  ! -----------------------------------------------------------------------
679  ! autoconversion from cloud ice to snow
680  ! -----------------------------------------------------------------------
681 
682  do i = is, ie
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)
688  endif
689  enddo
690 
691  if (out_dt) then
692  do i = is, ie
693  dtdt(i, j) = dtdt(i, j) + pt1(i) - t0(i)
694  enddo
695  endif
696 
697  ! -----------------------------------------------------------------------
698  ! fix energy conservation
699  ! -----------------------------------------------------------------------
700 
701  if (consv_te) then
702  do i = is, ie
703  if (hydrostatic) then
704 #ifdef MULTI_GASES
705  c_air = cp_air * vicpqd_qpz(qvi(i,j,1,1:num_gas),qpz(i))
706 #endif
707  te0(i, j) = dp(i, j) * (te0(i, j) + c_air * pt1(i))
708  else
709 #ifdef USE_COND
710  te0(i, j) = dp(i, j) * (te0(i, j) + cvm(i) * pt1(i))
711 #else
712 #ifdef MULTI_GASES
713  c_air = cv_air * vicvqd_qpz(qvi(i,j,1,1:num_gas),qpz(i))
714 #endif
715  te0(i, j) = dp(i, j) * (te0(i, j) + c_air * pt1(i))
716 #endif
717  endif
718  enddo
719  endif
720 
721  ! -----------------------------------------------------------------------
722  ! update latend heat coefficient
723  ! -----------------------------------------------------------------------
724 
725  do i = is, ie
726  lhi(i) = li00 + dc_ice * pt1(i)
727  lhl(i) = lv00 + d0_vap * 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)
731  enddo
732 
733  ! -----------------------------------------------------------------------
734  ! compute cloud fraction
735  ! -----------------------------------------------------------------------
736 
737  if (do_qa .and. last_step) then
738 
739  ! -----------------------------------------------------------------------
740  ! combine water species
741  ! -----------------------------------------------------------------------
742 
743  if (rad_snow) then
744  if (rad_graupel) then
745  do i = is, ie
746  q_sol(i) = qi(i, j) + qs(i, j) + qg(i, j)
747  enddo
748  else
749  do i = is, ie
750  q_sol(i) = qi(i, j) + qs(i, j)
751  enddo
752  endif
753  else
754  do i = is, ie
755  q_sol(i) = qi(i, j)
756  enddo
757  endif
758  if (rad_rain) then
759  do i = is, ie
760  q_liq(i) = ql(i, j) + qr(i, j)
761  enddo
762  else
763  do i = is, ie
764  q_liq(i) = ql(i, j)
765  enddo
766  endif
767  do i = is, ie
768  q_cond(i) = q_sol(i) + q_liq(i)
769  enddo
770 
771  ! -----------------------------------------------------------------------
772  ! use the "liquid - frozen water temperature" (tin) to compute saturated specific humidity
773  ! -----------------------------------------------------------------------
774 
775  do i = is, ie
776 
777  if(tintqs) then
778  tin = pt1(i)
779  else
780  tin = pt1(i) - (lcp2(i) * q_cond(i) + icp2(i) * q_sol(i)) ! minimum temperature
781  ! tin = pt1 (i) - ((lv00 + d0_vap * pt1 (i)) * q_cond (i) + &
782  ! (li00 + dc_ice * pt1 (i)) * q_sol (i)) / (mc_air (i) + qpz (i) * c_vap)
783  endif
784 
785  ! -----------------------------------------------------------------------
786  ! determine saturated specific humidity
787  ! -----------------------------------------------------------------------
788 
789  if (tin <= t_wfr) then
790  ! ice phase:
791  qstar(i) = iqs1(tin, den(i))
792  elseif (tin >= tice) then
793  ! liquid phase:
794  qstar(i) = wqs1(tin, den(i))
795  else
796  ! mixed phase:
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)
801  else
802  ! mostly liquid water clouds at initial cloud development stage
803  rqi = ((tice - tin) / (tice - t_wfr))
804  endif
805  qstar(i) = rqi * qsi + (1. - rqi) * qsw
806  endif
807 
808  ! higher than 10 m is considered "land" and will have higher subgrid variability
809  dw = dw_ocean + (dw_land - dw_ocean) * min(1., abs(hs(i, j)) / (10. * grav))
810  ! "scale - aware" subgrid variability: 100 - km as the base
811  hvar(i) = min(0.2, max(0.01, dw * sqrt(sqrt(area(i, j)) / 100.e3)))
812 
813  ! -----------------------------------------------------------------------
814  ! partial cloudiness by pdf:
815  ! assuming subgrid linear distribution in horizontal; this is effectively a smoother for the
816  ! binary cloud scheme; qa = 0.5 if qstar (i) == qpz
817  ! -----------------------------------------------------------------------
818 
819  rh = qpz(i) / qstar(i)
820 
821  ! -----------------------------------------------------------------------
822  ! icloud_f = 0: bug - fixed
823  ! icloud_f = 1: old fvgfs gfdl) mp implementation
824  ! icloud_f = 2: binary cloud scheme (0 / 1)
825  ! -----------------------------------------------------------------------
826 
827  if (rh > 0.75 .and. qpz(i) > 1.e-8) then
828  dq = hvar(i) * qpz(i)
829  q_plus = qpz(i) + dq
830  q_minus = qpz(i) - dq
831  if (icloud_f == 2) then
832  if (qpz(i) > qstar(i)) then
833  qa(i, j) = 1.
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))
837  else
838  qa(i, j) = 0.
839  endif
840  else
841  if (qstar(i) < q_minus) then
842  qa(i, j) = 1.
843  else
844  if (qstar(i) < q_plus) then
845  if (icloud_f == 0) then
846  qa(i, j) = (q_plus - qstar(i)) / (dq + dq)
847  else
848  qa(i, j) = (q_plus - qstar(i)) / (2. * dq * (1. - q_cond(i)))
849  endif
850  else
851  qa(i, j) = 0.
852  endif
853  ! impose minimum cloudiness if substantial q_cond (i) exist
854  if (q_cond(i) > 1.e-8) then
855  qa(i, j) = max(cld_min, qa(i, j))
856  endif
857  qa(i, j) = min(1., qa(i, j))
858  endif
859  endif
860  else
861  qa(i, j) = 0.
862  endif
863 
864  enddo
865 
866  endif
867 
868  enddo ! end j loop
869 
870 end subroutine fv_sat_adj
871 
872 ! =======================================================================
875 ! =======================================================================
876 real function wqs1 (ta, den)
877 
878  implicit none
879 
880  ! pure water phase; universal dry / moist formular using air density
881  ! input "den" can be either dry or moist air density
882 
883  real, intent (in) :: ta, den
884 
885  real :: es, ap1, tmin
886 
887  integer :: it
888 
889  tmin = tice - 160.
890  ap1 = 10. * dim(ta, tmin) + 1.
891  ap1 = min(2621., ap1)
892  it = ap1
893  es = tablew(it) + (ap1 - it) * desw(it)
894  wqs1 = es / (rvgas * ta * den)
895 
896 end function wqs1
897 
898 ! =======================================================================
901 ! =======================================================================
902 real function iqs1 (ta, den)
903 
904  implicit none
905 
906  ! water - ice phase; universal dry / moist formular using air density
907  ! input "den" can be either dry or moist air density
908 
909  real, intent (in) :: ta, den
910 
911  real :: es, ap1, tmin
912 
913  integer :: it
914 
915  tmin = tice - 160.
916  ap1 = 10. * dim(ta, tmin) + 1.
917  ap1 = min(2621., ap1)
918  it = ap1
919  es = table2(it) + (ap1 - it) * des2(it)
920  iqs1 = es / (rvgas * ta * den)
921 
922 end function iqs1
923 
924 ! =======================================================================
927 ! =======================================================================
928 real function wqs2 (ta, den, dqdt)
929 
930  implicit none
931 
932  ! pure water phase; universal dry / moist formular using air density
933  ! input "den" can be either dry or moist air density
934 
935  real, intent (in) :: ta, den
936 
937  real, intent (out) :: dqdt
938 
939  real :: es, ap1, tmin
940 
941  integer :: it
942 
943  tmin = tice - 160.
944  ap1 = 10. * dim(ta, tmin) + 1.
945  ap1 = min(2621., ap1)
946  it = ap1
947  es = tablew(it) + (ap1 - it) * desw(it)
948  wqs2 = es / (rvgas * ta * den)
949  it = ap1 - 0.5
950  ! finite diff, del_t = 0.1:
951  dqdt = 10. * (desw(it) + (ap1 - it) * (desw(it + 1) - desw(it))) / (rvgas * ta * den)
952 
953 end function wqs2
954 
955 ! =======================================================================
959 ! =======================================================================
960 subroutine wqs2_vect (is, ie, ta, den, wqsat, dqdt)
961 
962  implicit none
963 
964  ! pure water phase; universal dry / moist formular using air density
965  ! input "den" can be either dry or moist air density
966 
967  integer, intent (in) :: is, ie
968 
969  real, intent (in), dimension (is:ie) :: ta, den
970 
971  real, intent (out), dimension (is:ie) :: wqsat, dqdt
972 
973  real :: es, ap1, tmin
974 
975  integer :: i, it
976 
977  tmin = tice - 160.
978 
979  do i = is, ie
980  ap1 = 10. * dim(ta(i), tmin) + 1.
981  ap1 = min(2621., ap1)
982  it = ap1
983  es = tablew(it) + (ap1 - it) * desw(it)
984  wqsat(i) = es / (rvgas * ta(i) * den(i))
985  it = ap1 - 0.5
986  ! finite diff, del_t = 0.1:
987  dqdt(i) = 10. * (desw(it) + (ap1 - it) * (desw(it + 1) - desw(it))) / (rvgas * ta(i) * den(i))
988  enddo
989 
990 end subroutine wqs2_vect
991 
992 ! =======================================================================
995 ! =======================================================================
996 real function iqs2 (ta, den, dqdt)
997 
998  implicit none
999 
1000  ! water - ice phase; universal dry / moist formular using air density
1001  ! input "den" can be either dry or moist air density
1002 
1003  real, intent (in) :: ta, den
1004 
1005  real, intent (out) :: dqdt
1006 
1007  real :: es, ap1, tmin
1008 
1009  integer :: it
1010 
1011  tmin = tice - 160.
1012  ap1 = 10. * dim(ta, tmin) + 1.
1013  ap1 = min(2621., ap1)
1014  it = ap1
1015  es = table2(it) + (ap1 - it) * des2(it)
1016  iqs2 = es / (rvgas * ta * den)
1017  it = ap1 - 0.5
1018  ! finite diff, del_t = 0.1:
1019  dqdt = 10. * (des2(it) + (ap1 - it) * (des2(it + 1) - des2(it))) / (rvgas * ta * den)
1020 
1021 end function iqs2
1022 
1023 ! =======================================================================
1024 ! initialization
1025 ! prepare saturation water vapor pressure tables
1026 ! =======================================================================
1028 subroutine qs_init (kmp)
1030  implicit none
1031 
1032  integer, intent (in) :: kmp
1033 
1034  integer, parameter :: length = 2621
1035 
1036  integer :: i
1037 
1038  if (mp_initialized) return
1039 
1040  if (is_master()) write (*, *) 'top layer for gfdl_mp = ', kmp
1041 
1042  ! generate es table (dt = 0.1 deg c)
1043 
1044  allocate (table(length))
1045  allocate (table2(length))
1046  allocate (tablew(length))
1047  allocate (des2(length))
1048  allocate (desw(length))
1049 
1050  call qs_table (length)
1051  call qs_table2 (length)
1052  call qs_tablew (length)
1053 
1054  do i = 1, length - 1
1055  des2(i) = max(0., table2(i + 1) - table2(i))
1056  desw(i) = max(0., tablew(i + 1) - tablew(i))
1057  enddo
1058  des2(length) = des2(length - 1)
1059  desw(length) = desw(length - 1)
1060 
1061  mp_initialized = .true.
1062 
1063 end subroutine qs_init
1064 
1065 ! =======================================================================
1066 ! saturation water vapor pressure table i
1067 ! 3 - phase table
1068 ! =======================================================================
1069 
1070 subroutine qs_table (n)
1072  implicit none
1073 
1074  integer, intent (in) :: n
1075 
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)
1080 
1081  integer :: i
1082 
1083  tmin = tice - 160.
1084 
1085  ! -----------------------------------------------------------------------
1086  ! compute es over ice between - 160 deg c and 0 deg c.
1087  ! -----------------------------------------------------------------------
1088 
1089  do i = 1, 1600
1090  tem = tmin + delt * real(i - 1)
1091  fac0 = (tem - tice) / (tem * tice)
1092  fac1 = fac0 * li2
1093  fac2 = (d2ice * log(tem / tice) + fac1) / rvgas
1094  table(i) = e00 * exp(fac2)
1095  enddo
1096 
1097  ! -----------------------------------------------------------------------
1098  ! compute es over water between - 20 deg c and 102 deg c.
1099  ! -----------------------------------------------------------------------
1100 
1101  do i = 1, 1221
1102  tem = 253.16 + delt * real(i - 1)
1103  fac0 = (tem - tice) / (tem * tice)
1104  fac1 = fac0 * lv0
1105  fac2 = (dc_vap * log(tem / tice) + fac1) / rvgas
1106  esh20 = e00 * exp(fac2)
1107  if (i <= 200) then
1108  esupc(i) = esh20
1109  else
1110  table(i + 1400) = esh20
1111  endif
1112  enddo
1113 
1114  ! -----------------------------------------------------------------------
1115  ! derive blended es over ice and supercooled water between - 20 deg c and 0 deg c
1116  ! -----------------------------------------------------------------------
1117 
1118  do i = 1, 200
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)
1123  enddo
1124 
1125 end subroutine qs_table
1126 
1127 ! =======================================================================
1128 ! saturation water vapor pressure table ii
1129 ! 1 - phase table
1130 ! =======================================================================
1131 
1132 subroutine qs_tablew (n)
1134  implicit none
1135 
1136  integer, intent (in) :: n
1137 
1138  real (kind = r_grid) :: delt = 0.1
1139  real (kind = r_grid) :: tmin, tem, fac0, fac1, fac2
1140 
1141  integer :: i
1142 
1143  tmin = tice - 160.
1144 
1145  ! -----------------------------------------------------------------------
1146  ! compute es over water
1147  ! -----------------------------------------------------------------------
1148 
1149  do i = 1, n
1150  tem = tmin + delt * real(i - 1)
1151  fac0 = (tem - tice) / (tem * tice)
1152  fac1 = fac0 * lv0
1153  fac2 = (dc_vap * log(tem / tice) + fac1) / rvgas
1154  tablew(i) = e00 * exp(fac2)
1155  enddo
1156 
1157 end subroutine qs_tablew
1158 
1159 ! =======================================================================
1160 ! saturation water vapor pressure table iii
1161 ! 2 - phase table
1162 ! =======================================================================
1163 
1164 subroutine qs_table2 (n)
1166  implicit none
1167 
1168  integer, intent (in) :: n
1169 
1170  real (kind = r_grid) :: delt = 0.1
1171  real (kind = r_grid) :: tmin, tem0, tem1, fac0, fac1, fac2
1172 
1173  integer :: i, i0, i1
1174 
1175  tmin = tice - 160.
1176 
1177  do i = 1, n
1178  tem0 = tmin + delt * real(i - 1)
1179  fac0 = (tem0 - tice) / (tem0 * tice)
1180  if (i <= 1600) then
1181  ! -----------------------------------------------------------------------
1182  ! compute es over ice between - 160 deg c and 0 deg c.
1183  ! -----------------------------------------------------------------------
1184  fac1 = fac0 * li2
1185  fac2 = (d2ice * log(tem0 / tice) + fac1) / rvgas
1186  else
1187  ! -----------------------------------------------------------------------
1188  ! compute es over water between 0 deg c and 102 deg c.
1189  ! -----------------------------------------------------------------------
1190  fac1 = fac0 * lv0
1191  fac2 = (dc_vap * log(tem0 / tice) + fac1) / rvgas
1192  endif
1193  table2(i) = e00 * exp(fac2)
1194  enddo
1195 
1196  ! -----------------------------------------------------------------------
1197  ! smoother around 0 deg c
1198  ! -----------------------------------------------------------------------
1199 
1200  i0 = 1600
1201  i1 = 1601
1202  tem0 = 0.25 * (table2(i0 - 1) + 2. * table(i0) + table2(i0 + 1))
1203  tem1 = 0.25 * (table2(i1 - 1) + 2. * table(i1) + table2(i1 + 1))
1204  table2(i0) = tem0
1205  table2(i1) = tem1
1206 
1207 end subroutine qs_table2
1208 
1209 end module fv_cmp_mod
real, dimension(:), allocatable des2
Definition: fv_cmp.F90:113
real, parameter cp_vap
1846.0, heat capacity of water vapor at constant pressure
Definition: fv_cmp.F90:75
real, dimension(:), allocatable table
Definition: fv_cmp.F90:113
real, parameter dc_vap
Definition: fv_cmp.F90:93
The module &#39;fv_mp_mod&#39; is a single program multiple data (SPMD) parallel decompostion/communication m...
Definition: fv_mp_mod.F90:24
real, dimension(:), allocatable tablew
Definition: fv_cmp.F90:113
real, parameter tice
freezing temperature
Definition: fv_cmp.F90:96
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 &#39;fv_sat_adj&#39; performs the fast processes in the GFDL microphysics. ...
Definition: fv_cmp.F90:130
integer, public num_gas
Definition: multi_gases.F90:46
real, dimension(:), allocatable table2
Definition: fv_cmp.F90:113
real, parameter cv_air
717.55, heat capacity of dry air at constant volume
Definition: fv_cmp.F90:76
The module &#39;multi_gases&#39; peforms multi constitutents computations.
Definition: multi_gases.F90:25
real, parameter cv_vap
1384.5, heat capacity of water vapor at constant volume
Definition: fv_cmp.F90:77
real, parameter c_ice
gfdl: heat capacity of ice at - 15 deg c
Definition: fv_cmp.F90:90
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
Definition: fv_cmp.F90:110
pure real function, public virq_qpz(q, qpz)
integer, parameter, public r_grid
Definition: fv_arrays.F90:34
subroutine qs_table2(n)
Definition: fv_cmp.F90:1165
real function iqs1(ta, den)
the function &#39;wqs1&#39; computes the saturated specific humidity for table iii
Definition: fv_cmp.F90:903
subroutine qs_table(n)
Definition: fv_cmp.F90:1071
real(kind=r_grid), parameter d2ice
Definition: fv_cmp.F90:105
real, dimension(:), allocatable desw
Definition: fv_cmp.F90:113
The module &#39;fv_arrays&#39; contains the &#39;fv_atmos_type&#39; and associated datatypes.
Definition: fv_arrays.F90:24
real lv00
the same as lv0, except that cp_vap can be cp_vap or cv_vap
Definition: fv_cmp.F90:111
real, parameter c_liq
gfdl: heat capacity of liquid at 15 deg c
Definition: fv_cmp.F90:91
real function iqs2(ta, den, dqdt)
The function &#39;iqs2&#39; computes the gradient of saturated specific humidity for table iii...
Definition: fv_cmp.F90:997
real, parameter dc_ice
2213.5, isobaric heating / colling
Definition: fv_cmp.F90:94
subroutine qs_tablew(n)
Definition: fv_cmp.F90:1133
real, parameter lat2
used in bigg mechanism
Definition: fv_cmp.F90:108
real, parameter li00
Definition: fv_cmp.F90:100
real, parameter t_wfr
homogeneous freezing temperature
Definition: fv_cmp.F90:97
real, parameter lv0
3.13905782e6, evaporation latent heat coefficient at 0 deg k
Definition: fv_cmp.F90:99
pure real function, public vicpqd_qpz(q, qpz)
real(kind=r_grid), parameter li2
2.86799816e6, sublimation latent heat coefficient at 0 deg k
Definition: fv_cmp.F90:106
logical mp_initialized
Definition: fv_cmp.F90:115
real function wqs2(ta, den, dqdt)
The function &#39;wqs2&#39;computes the gradient of saturated specific humidity for table ii...
Definition: fv_cmp.F90:929
real function wqs1(ta, den)
the function &#39;wqs1&#39; computes the saturated specific humidity for table ii
Definition: fv_cmp.F90:877
subroutine, public qs_init(kmp)
The subroutine &#39;qs_init&#39; initializes lookup tables for the saturation mixing ratio.
Definition: fv_cmp.F90:1029
real(kind=r_grid), parameter e00
ifs: saturation vapor pressure at 0 deg c
Definition: fv_cmp.F90:103
subroutine wqs2_vect(is, ie, ta, den, wqsat, dqdt)
The function wqs2_vect computes the gradient of saturated specific humidity for table ii...
Definition: fv_cmp.F90:961
The module &#39;fv_cmp&#39; implements the fast procesesses in the GFDL microphysics >
Definition: fv_cmp.F90:29