Line data Source code
1 : ! Prognostic cloud water data and methods (cldwat)
2 : ! Original authors as marked in subroutines
3 : ! CCPP-ized: Haipeng Lin, January 2025
4 : module prognostic_cloud_water
5 : use ccpp_kinds, only: kind_phys
6 :
7 : implicit none
8 : private
9 : save
10 :
11 : ! public CCPP-compliant subroutines
12 : public :: prognostic_cloud_water_init
13 : public :: prognostic_cloud_water_run
14 :
15 : ! tuning parameters used by prognostic cloud water (RK stratiform)
16 : real(kind_phys), public :: icritc ! threshold for autoconversion of cold ice [kg kg-1]
17 : ! REMOVECAM: icritc does not have to be public after CAM is retired (namelist option; CCPP framework can provide it)
18 : real(kind_phys) :: icritw ! threshold for autoconversion of warm ice [kg kg-1]
19 : real(kind_phys) :: conke ! tunable constant for evaporation of precipitation [kg-0.5 m s-0.5]
20 : real(kind_phys) :: r3lcrit ! critical radius where liq conversion begins [m]
21 :
22 : logical :: do_psrhmin ! set rh in stratosphere poleward of 50 degrees [flag]
23 : real(kind_phys) :: psrhmin ! rh set in stratosphere poleward of 50 degrees [%]
24 :
25 : ! module private variables
26 : real(kind_phys) :: rhonot ! air density at surface [g cm-3]
27 : real(kind_phys) :: rhos ! assumed snow density [g cm-3]
28 : real(kind_phys) :: rhow ! water density [g cm-3]
29 : real(kind_phys) :: rhoi ! ice density [g cm-3]
30 : real(kind_phys) :: esi ! Collection efficiency for ice by snow [1]
31 : real(kind_phys) :: esw ! Collection efficiency for water by snow [1]
32 : real(kind_phys) :: t0 ! Approx. freezing temperature [K]
33 : real(kind_phys) :: cldmin ! Assumed minimum cloud amount [1]
34 : real(kind_phys) :: small ! Small number compared to unity [1]
35 : real(kind_phys) :: c ! Constant for graupel-like snow [cm^(1-d)/s]
36 :
37 : real(kind_phys) :: d ! Constant for graupel-like snow [1]
38 : real(kind_phys) :: thrpd ! 3+d
39 : real(kind_phys) :: gam3pd ! gamma(3+d)
40 : real(kind_phys) :: gam4pd ! gamma(4+d)
41 :
42 : real(kind_phys) :: nos ! Particles snow concentration [cm^-4]
43 : real(kind_phys) :: prhonos ! pi * rhos * nos
44 :
45 : ! cloud microphysics constants
46 : real(kind_phys) :: mcon01
47 : real(kind_phys) :: mcon02
48 : real(kind_phys) :: mcon03
49 : real(kind_phys) :: mcon04
50 : real(kind_phys) :: mcon05
51 : real(kind_phys) :: mcon06
52 : real(kind_phys) :: mcon07
53 : real(kind_phys) :: mcon08
54 :
55 : ! Parameters for findmcnew
56 : real(kind_phys) :: capnw ! Warm continental cloud particle number concentration [cm-3]
57 : real(kind_phys) :: capnc ! Cold/oceanic cloud particle number concentration [cm-3]
58 : real(kind_phys) :: capnsi ! Sea ice cloud particle number concentration [cm-3]
59 : real(kind_phys) :: kconst ! Terminal velocity constant (Stokes regime) [1]
60 : real(kind_phys) :: effc ! Autoconv collection efficiency [1]
61 : real(kind_phys) :: alpha ! Ratio of 3rd moment radius to 2nd [1]
62 : real(kind_phys) :: capc ! Autoconversion constant [1]
63 : real(kind_phys) :: convfw ! Fall velocity calculation constant [1]
64 : real(kind_phys) :: cracw ! Rain accreting water constant [1]
65 : real(kind_phys) :: critpr ! Critical precip rate for collection efficiency [mm day-1]
66 : real(kind_phys) :: ciautb ! Ice autoconversion coefficient [s-1]
67 :
68 : contains
69 :
70 : ! Initialize prognostic cloud water module and calculate microphysical constants
71 : !> \section arg_table_prognostic_cloud_water_init Argument Table
72 : !! \htmlinclude arg_table_prognostic_cloud_water_init.html
73 1024 : subroutine prognostic_cloud_water_init( &
74 : amIRoot, iulog, &
75 : tmelt, rhodair, pi, &
76 : icritc_in, icritw_in, &
77 : conke_in, r3lcrit_in, &
78 : do_psrhmin_in, psrhmin_in, &
79 1024 : errmsg, errflg)
80 :
81 : ! Input arguments
82 : logical, intent(in) :: amIRoot ! are we on the MPI root task?
83 : integer, intent(in) :: iulog ! log output unit
84 :
85 : real(kind_phys), intent(in) :: tmelt ! freezing_point_of_water [K]
86 : real(kind_phys), intent(in) :: rhodair ! density_of_dry_air_at_stp [kg m-3]
87 : real(kind_phys), intent(in) :: pi
88 :
89 : real(kind_phys), intent(in) :: icritc_in
90 : real(kind_phys), intent(in) :: icritw_in
91 : real(kind_phys), intent(in) :: conke_in
92 : real(kind_phys), intent(in) :: r3lcrit_in
93 : logical, intent(in) :: do_psrhmin_in
94 : real(kind_phys), intent(in) :: psrhmin_in
95 :
96 : ! Output arguments
97 : character(len=512), intent(out) :: errmsg ! error message
98 : integer, intent(out) :: errflg ! error flag
99 :
100 1024 : errmsg = ''
101 1024 : errflg = 0
102 :
103 : ! First populate tuning parameters in-module
104 1024 : icritc = icritc_in
105 1024 : icritw = icritw_in
106 1024 : conke = conke_in
107 1024 : r3lcrit = r3lcrit_in
108 1024 : do_psrhmin = do_psrhmin_in
109 1024 : psrhmin = psrhmin_in
110 :
111 1024 : if(amIRoot) then
112 2 : write(iulog,*) 'tuning parameters prognostic_cloud_water_init: icritw ',icritw,' icritc ',icritc,' conke ',conke,' r3lcrit ',r3lcrit
113 2 : write(iulog,*) 'prognostic_cloud_water_init: do_psrhmin = ', do_psrhmin
114 : endif
115 :
116 : !--------------------------------------------------
117 : ! Initialize constants used for prognostic condensate (inimc)
118 : ! Original author: P. Rasch, April 1997
119 : !--------------------------------------------------
120 :
121 1024 : rhonot = rhodair/1000.0_kind_phys ! convert from kg m-3 to g cm-3
122 :
123 : ! assumed densities of snow, water, ice [g cm-3]
124 1024 : rhos = 0.1_kind_phys
125 1024 : rhow = 1._kind_phys
126 1024 : rhoi = 1._kind_phys
127 :
128 1024 : esi = 1._kind_phys
129 1024 : esw = 0.1_kind_phys
130 1024 : t0 = tmelt
131 1024 : cldmin = 0.02_kind_phys
132 1024 : small = 1.e-22_kind_phys
133 1024 : c = 152.93_kind_phys
134 1024 : d = 0.25_kind_phys ! constant for graupel like snow
135 : ! values other than 0.25 are not supported.
136 1024 : if(d == 0.25_kind_phys) then
137 1024 : gam3pd = 2.549256966718531_kind_phys
138 1024 : gam4pd = 8.285085141835282_kind_phys
139 :
140 : ! on cray machines this can be calculated using gamma function:
141 : ! call gamma(3._kind_phys+d, signgam, gam3pd)
142 : ! gam3pd = sign(exp(gam3pd),signgam)
143 : ! call gamma(4._kind_phys+d, signgam, gam4pd)
144 : ! gam4pd = sign(exp(gam4pd),signgam)
145 : ! write(iulog,*) ' d, gamma(3+d), gamma(4+d) =', gam3pd, gam4pd
146 : else
147 0 : errflg = 1
148 0 : errmsg = 'prognostic_cloud_water_init: cannot use d /= 0.25'
149 : endif
150 :
151 1024 : nos = 3.e-2_kind_phys
152 :
153 : ! note: pi was originally calculated here as 4._kind_phys * atan(1.0_kind_phys)
154 : ! changed in ccpp-ization to use physconst value.
155 1024 : prhonos = pi*rhos*nos
156 1024 : thrpd = 3._kind_phys + d
157 :
158 1024 : mcon01 = pi*nos*c*gam3pd/4._kind_phys
159 1024 : mcon02 = 1._kind_phys/(c*gam4pd*sqrt(rhonot)/(6*prhonos**(d/4._kind_phys)))
160 1024 : mcon03 = -(0.5_kind_phys+d/4._kind_phys)
161 1024 : mcon04 = 4._kind_phys/(4._kind_phys+d)
162 1024 : mcon05 = (3+d)/(4+d)
163 1024 : mcon06 = (3+d)/4._kind_phys
164 1024 : mcon07 = mcon01*sqrt(rhonot)*mcon02**mcon05/prhonos**mcon06
165 1024 : mcon08 = -0.5_kind_phys/(4._kind_phys+d)
166 :
167 : ! Initialize parameters used by findmcnew
168 : ! Cloud particle densities [cm-3] for ...
169 1024 : capnw = 400._kind_phys ! ...warm continental
170 1024 : capnc = 150._kind_phys ! ...cold and oceanic
171 1024 : capnsi = 75._kind_phys ! ...sea ice
172 :
173 1024 : kconst = 1.18e6_kind_phys
174 :
175 : ! Autoconv collection efficiencies.
176 : ! Tripoli and Cotton (default) = 0.55
177 : ! Boucher 96 = 1.
178 : ! Baker 93 = 0.55*0.05
179 : ! Turn off = 0
180 1024 : effc = 0.55_kind_phys
181 :
182 1024 : alpha = 1.1_kind_phys ** 4
183 : ! constant for autoconversion
184 1024 : capc = pi**(-.333_kind_phys)*kconst*effc *(0.75_kind_phys)**(1.333_kind_phys)*alpha
185 :
186 : ! critical precip rate at which we assume the collector drops can change the
187 : ! drop size enough to enhance the auto-conversion process (mm/day)
188 1024 : critpr = 0.5_kind_phys
189 1024 : convfw = 1.94_kind_phys*2.13_kind_phys*sqrt(rhow*1000._kind_phys*9.81_kind_phys*2.7e-4_kind_phys)
190 :
191 : ! liquid microphysics - rain accreting water constant
192 : ! Tripoli and Cotton = default
193 : ! Beheng = 6.
194 1024 : cracw = .884_kind_phys*sqrt(9.81_kind_phys/(rhow*1000._kind_phys*2.7e-4_kind_phys)) ! tripoli and cotton
195 :
196 : ! ice microphysics autoconversion coefficient
197 1024 : ciautb = 5.e-4_kind_phys
198 :
199 1024 : if(amIRoot) then
200 2 : write(iulog,*) 'tuning parameters prognostic_cloud_water_init: capnw ',capnw,' capnc ',capnc,' capnsi ',capnsi,' kconst ',kconst
201 2 : write(iulog,*) 'tuning parameters prognostic_cloud_water_init: effc ',effc,' alpha ',alpha,' capc ',capc
202 2 : write(iulog,*) 'tuning parameters prognostic_cloud_water_init: critpr ',critpr,' convfw ',convfw,' cracw ',cracw,' ciautb ',ciautb
203 : endif
204 :
205 1024 : end subroutine prognostic_cloud_water_init
206 :
207 : ! Calculate prognostic condensate.
208 : ! Cloud water parameterization, returns tendencies to water vapor, temperature,
209 : ! and cloud water variables.
210 : !
211 : ! Rasch, P. J., and J. E. Kristjánsson, 1998: A Comparison of the CCM3 Model Climate Using Diagnosed and Predicted Condensate Parameterizations. J. Climate, 11, 1587–1614
212 : ! https://doi.org/10.1175/1520-0442(1998)011<1587:ACOTCM>2.0.CO;2
213 : !
214 : ! with modifications to improve the method of determining condensation/evaporation
215 : ! Zhang, M., W. Lin, C. Bretherton, J. Hack, and P. J. Rasch, 2003, A modified formulation of fractional stratiform condensation rate in the NCAR Community Atmospheric Model (CAM2), J. Geophys. Res., 108(D1), 4035
216 : ! https://doi.org/10.1029/2002JD002523
217 : !
218 : ! Original authors: M. Zhang, W. Lin, P. Rasch and J.E. Kristjansson
219 : ! B. A. Boville (latent heat of fusion)
220 : !> \section arg_table_prognostic_cloud_water_run Argument Table
221 : !! \htmlinclude arg_table_prognostic_cloud_water_run.html
222 70392 : subroutine prognostic_cloud_water_run( &
223 : ncol, pver, top_lev, deltat, &
224 : iulog, &
225 : pi, gravit, rh2o, epsilo, latvap, latice, cpair, &
226 70392 : dlat, &
227 140784 : pmid, pdel, &
228 70392 : zi, &
229 70392 : troplev, &
230 140784 : ttend, tn, &
231 140784 : qtend, qn, &
232 70392 : ltend, &
233 70392 : cldice, cldliq, &
234 : omega, &
235 70392 : cldn, &
236 140784 : fice, fsnow, &
237 140784 : rhdfda, rhu00, &
238 70392 : landm, seaicef, snowh, &
239 70392 : qme, &
240 140784 : prodprec, prodsnow, &
241 140784 : evapprec, evapsnow, &
242 70392 : evapheat, prfzheat, meltheat, &
243 70392 : precip, snowab, &
244 70392 : ice2pr, liq2pr, liq2snow, &
245 70392 : lsflxprc, &
246 70392 : lsflxsnw, &
247 70392 : pracwo, psacwo, psacio, &
248 70392 : fwaut, fsaut, fracw, fsacw, fsaci, &
249 0 : errmsg, errflg)
250 : ! Dependencies to-be-ccppized.
251 : use wv_saturation, only: qsat, estblf, svp_to_qsat, findsp_vc
252 :
253 : ! Input arguments
254 : integer, intent(in) :: ncol
255 : integer, intent(in) :: pver
256 : integer, intent(in) :: top_lev ! vertical_layer_index_of_troposphere_cloud_physics_top [index]
257 : real(kind_phys), intent(in) :: deltat ! timestep [s]
258 : integer, intent(in) :: iulog ! log output unit [1]
259 : real(kind_phys), intent(in) :: pi ! pi_constant [1]
260 : real(kind_phys), intent(in) :: gravit ! gravitational acceleration [m s-2]
261 : real(kind_phys), intent(in) :: rh2o ! gas_constant_of_water_vapor [J K-1 kg-1]
262 : real(kind_phys), intent(in) :: epsilo ! ratio_of_water_vapor_to_dry_air_molecular_weights [1]
263 : real(kind_phys), intent(in) :: latvap ! latent_heat_of_vaporization_of_water_at_0c [J kg-1]
264 : real(kind_phys), intent(in) :: latice ! latent_heat_of_fusion_of_water_at_0c [J kg-1]
265 : real(kind_phys), intent(in) :: cpair ! specific_heat_of_dry_air_at_constant_pressure [J K-1 kg-1]
266 : real(kind_phys), intent(in) :: dlat(:) ! latitude_degrees_north [degrees]
267 : real(kind_phys), intent(in) :: pmid(:,:) ! air_pressure [Pa]
268 : real(kind_phys), intent(in) :: pdel(:,:) ! air_pressure_thickness [Pa]
269 : real(kind_phys), intent(in) :: zi(:,:) ! geopotential_height_wrt_surface_at_interface [m]
270 : integer, intent(in) :: tropLev(:) ! tropopause_vertical_layer_index [index]
271 :
272 : real(kind_phys), intent(in) :: ttend(:,:) ! Temperature tendency [K s-1] -- from non-micro/macrophysics
273 : real(kind_phys), intent(in) :: tn(:,:) ! New Temperature [K]
274 : real(kind_phys), intent(in) :: qtend(:,:) ! Water vapor tendency [kg kg-1 s-1] -- from non-micro/macrophysics
275 : real(kind_phys), intent(in) :: qn(:,:) ! New Water vapor mixing ratio [kg kg-1]
276 : real(kind_phys), intent(in) :: ltend(:,:) ! Cloud liquid water tendency [kg kg-1 s-1] -- from non-micro/macrophysics
277 : real(kind_phys), intent(in) :: cldice(:,:) ! adv: cloud_ice_mixing_ratio_wrt_moist_air_and_condensed_water [kg kg-1]
278 : real(kind_phys), intent(in) :: cldliq(:,:) ! adv: cloud_liquid_water_mixing_ratio_wrt_moist_air_and_condensed_water [kg kg-1]
279 :
280 : real(kind_phys), intent(in) :: omega(:,:) ! lagrangian_tendency_of_air_pressure [Pa s-1]
281 : real(kind_phys), intent(in) :: cldn(:,:) ! New Cloud fraction [fraction]
282 : real(kind_phys), intent(in) :: fice(:,:) ! Ice fraction within cwat [fraction]
283 : real(kind_phys), intent(in) :: fsnow(:,:) ! Fraction of rain that freezes to snow [fraction]
284 :
285 : real(kind_phys), intent(in) :: rhdfda(:,:) ! dG(a)/da, rh=G(a), when rh>u00 [??]
286 : real(kind_phys), intent(in) :: rhu00(:,:) ! Rhlim for cloud [percent]
287 : real(kind_phys), intent(in) :: landm(:) ! smoothed_land_area_fraction [fraction]
288 : real(kind_phys), intent(in) :: seaicef(:) ! Sea ice fraction [1 (fraction)]
289 : real(kind_phys), intent(in) :: snowh(:) ! Snow depth over land, water equivalent [m]
290 :
291 : ! Output arguments
292 : real(kind_phys), intent(out) :: qme(:,:) ! Rate of condensation-evaporation of condensate (net_condensation_rate_due_to_microphysics) [kg kg-1 s-1]
293 : real(kind_phys), intent(out) :: prodprec(:,:) ! Conversion rate of condensate to precip (precipitation_production_due_to_microphysics) [kg kg-1 s-1]
294 : real(kind_phys), intent(out) :: prodsnow(:,:) ! Snow production rate (ignored in RK?) [kg kg-1 s-1]
295 : real(kind_phys), intent(out) :: evapprec(:,:) ! Falling precipitation evaporation rate (precipitation_evaporation_due_to_microphysics) [kg kg-1 s-1] -- & combined to apply q(wv) tendency
296 : real(kind_phys), intent(out) :: evapsnow(:,:) ! Falling snow evaporation rate [kg kg-1 s-1]
297 : real(kind_phys), intent(out) :: evapheat(:,:) ! heating rate due to evaporation of precipitation [J kg-1 s-1]
298 : real(kind_phys), intent(out) :: prfzheat(:,:) ! heating rate due to freezing of precipitation [J kg-1 s-1]
299 : real(kind_phys), intent(out) :: meltheat(:,:) ! heating rate due to snow melt [J kg-1 s-1]
300 : ! note -- these are precip units for atmosphere-surface exchange.
301 : real(kind_phys), intent(out) :: precip(:) ! Precipitation rate (lwe_stratiform_precipitation_rate_at_surface) [m s-1]
302 : real(kind_phys), intent(out) :: snowab(:) ! Snow rate (lwe_snow_precipitation_rate_at_surface_due_to_microphysics) [m s-1]
303 : ! / note
304 : real(kind_phys), intent(out) :: ice2pr(:,:) ! Conversion rate of ice to precip [kg kg-1 s-1] -- apply q(cldice) tendency
305 : real(kind_phys), intent(out) :: liq2pr(:,:) ! Conversion rate of liquid to precip [kg kg-1 s-1] -- apply q(cldliq) tendency
306 : real(kind_phys), intent(out) :: liq2snow(:,:) ! Conversion rate of liquid to snow [kg kg-1 s-1] -- discard??
307 : real(kind_phys), intent(out) :: lsflxprc(:,:) ! Grid-box mean, flux large scale cloud rain+snow at interfaces (stratiform_rain_and_snow_flux_at_interface) [kg m-2 s-1]
308 : real(kind_phys), intent(out) :: lsflxsnw(:,:) ! Grid-box mean, flux large scale cloud snow at interfaces (stratiform_snow_flux_at_interface) [kg m-2 s-1]
309 : real(kind_phys), intent(out) :: pracwo(:,:) ! accretion of cloud water by rain [s-1]
310 : real(kind_phys), intent(out) :: psacwo(:,:) ! accretion of cloud water by snow [s-1]
311 : real(kind_phys), intent(out) :: psacio(:,:) ! accretion of cloud ice by snow [s-1]
312 : real(kind_phys), intent(out) :: fwaut(ncol,pver) ! Relative importance of liquid autoconversion [fraction]
313 : real(kind_phys), intent(out) :: fsaut(ncol,pver) ! Relative importance of ice autoconversion [fraction]
314 : real(kind_phys), intent(out) :: fracw(ncol,pver) ! Relative importance of liquid collection by rain [fraction]
315 : real(kind_phys), intent(out) :: fsacw(ncol,pver) ! Relative importance of liquid collection by snow [fraction]
316 : real(kind_phys), intent(out) :: fsaci(ncol,pver) ! Relative importance of ice collection by snow [fraction]
317 :
318 : ! Output arguments
319 : character(len=512), intent(out) :: errmsg ! error message
320 : integer, intent(out) :: errflg ! error flag
321 :
322 : ! Local variables
323 : integer :: i, k, l ! Iteration index [1]
324 : integer :: iter ! # of iterations for precipitation calculation [1]
325 : logical :: error_found ! Flag for error detection [flag]
326 :
327 : ! Total cloud water (from cldice+cldliq)
328 140784 : real(kind_phys) :: cwat(ncol,pver) ! Cloud water mixing ratio [kg kg-1]
329 :
330 : ! Precipitation and conversion rates
331 : real(kind_phys) :: nice2pr ! Rate of conversion from ice to snow [kg kg-1 s-1]
332 : real(kind_phys) :: nliq2pr ! Rate of conversion from liquid to precipitation [kg kg-1 s-1]
333 : real(kind_phys) :: nliq2snow ! Rate of conversion from liquid to snow [kg kg-1 s-1]
334 140784 : real(kind_phys) :: precab(ncol) ! Rate of precipitation entering layer [kg m-2 s-1]
335 140784 : real(kind_phys) :: prprov(ncol) ! Provisional precipitation at bottom of layer [kg m-2 s-1]
336 :
337 : ! Cloud properties
338 140784 : real(kind_phys) :: cldm(ncol) ! Mean cloud fraction over timestep [1]
339 140784 : real(kind_phys) :: cldmax(ncol) ! Maximum cloud fraction above current level [1]
340 140784 : real(kind_phys) :: icwc(ncol) ! In-cloud water content [kg kg-1]
341 140784 : real(kind_phys) :: cwm(ncol) ! Cloud water mixing ratio at midpoint of timestep [kg kg-1]
342 140784 : real(kind_phys) :: cwn(ncol) ! Cloud water mixing ratio at end of timestep [kg kg-1]
343 140784 : real(kind_phys) :: coef(ncol) ! Conversion timescale for condensate to rain [s-1]
344 :
345 : ! Thermodynamic variables
346 140784 : real(kind_phys) :: t(ncol,pver) ! Temperature before timestep [K]
347 140784 : real(kind_phys) :: tsp(ncol,pver) ! Saturation point temperature [K]
348 140784 : real(kind_phys) :: es(ncol) ! Saturation vapor pressure [Pa]
349 140784 : real(kind_phys) :: q(ncol,pver) ! Water vapor mixing ratio [kg kg-1]
350 140784 : real(kind_phys) :: qs(ncol) ! Saturation specific humidity [kg kg-1]
351 140784 : real(kind_phys) :: qsp(ncol,pver) ! Saturation point mixing ratio [kg kg-1]
352 140784 : real(kind_phys) :: relhum(ncol) ! Relative humidity [1]
353 140784 : real(kind_phys) :: relhum1(ncol) ! Updated relative humidity [1]
354 140784 : real(kind_phys) :: rhu_adj(ncol,pver) ! Adjusted relative humidity limit for strat dehydration [1]
355 :
356 : ! Cloud-Evaporation scheme variables
357 70392 : real(kind_phys) :: calpha(ncol) ! Alpha term in C-E formulation [1]
358 140784 : real(kind_phys) :: cbeta(ncol) ! Beta term in C-E formulation [1]
359 140784 : real(kind_phys) :: cbetah(ncol) ! Beta-hat at saturation portion [1]
360 140784 : real(kind_phys) :: cgamma(ncol) ! Gamma term in C-E formulation [1]
361 140784 : real(kind_phys) :: cgamah(ncol) ! Gamma-hat at saturation portion [1]
362 140784 : real(kind_phys) :: rcgama(ncol) ! Ratio of gamma to gamma-hat [1]
363 140784 : real(kind_phys) :: csigma(ncol) ! Sigma term in C-E formulation [1]
364 140784 : real(kind_phys) :: qmeres(ncol) ! Residual condensation after C-E and evapprec [kg kg-1 s-1]
365 140784 : real(kind_phys) :: qmec1(ncol) ! Cloud condensation coefficient 1 C-E formulation [1]
366 140784 : real(kind_phys) :: qmec2(ncol) ! Cloud condensation coefficient 2 C-E formulation [1]
367 140784 : real(kind_phys) :: qmec3(ncol) ! Cloud condensation coefficient 3 C-E formulation [1]
368 140784 : real(kind_phys) :: qmec4(ncol) ! Cloud condensation coefficient 4 C-E formulation [1]
369 :
370 : ! Diagnostic arrays for cloud water budget
371 : ! Hardcoded 2 here = number of iterations (iter below)
372 140784 : real(kind_phys) :: rcwn(ncol,2,pver) ! Cloud water ratio evolution [kg kg-1]
373 140784 : real(kind_phys) :: rliq(ncol,2,pver) ! Liquid water ratio evolution [kg kg-1]
374 140784 : real(kind_phys) :: rice(ncol,2,pver) ! Ice ratio evolution [kg kg-1]
375 :
376 : ! Constants and parameters
377 : real(kind_phys), parameter :: omsm = 0.99999_kind_phys ! unity minus small number for rounding
378 : real(kind_phys) :: mincld ! Minimum cloud fraction [1]
379 : real(kind_phys) :: cpohl ! Ratio of specific heat to latent heat [K-1]
380 : real(kind_phys) :: hlocp ! Ratio of latent heat to specific heat [K]
381 : real(kind_phys) :: clrh2o ! Ratio of latent heat to water vapor gas constant [K]
382 : real(kind_phys) :: dto2 ! Half timestep [s]
383 :
384 : ! Work variables
385 : real(kind_phys) :: denom ! Denominator work variable [1]
386 : real(kind_phys) :: dqsdt ! Change in saturation specific humidity with temperature [kg kg-1 K-1]
387 140784 : real(kind_phys) :: gamma(ncol) ! Temperature derivative of saturation specific humidity [kg kg-1 K-1]
388 140784 : real(kind_phys) :: qtl(ncol) ! Saturation tendency [kg kg-1 s-1]
389 : real(kind_phys) :: qtmp ! Temporary mixing ratio [kg kg-1]
390 : real(kind_phys) :: ttmp ! Temporary temperature [K]
391 : real(kind_phys) :: qsn ! Updated saturation specific humidity [kg kg-1]
392 : real(kind_phys) :: esn ! Updated saturation vapor pressure [Pa]
393 : real(kind_phys) :: prtmp ! Temporary precipitation rate [kg m-2 s-1]
394 : real(kind_phys) :: ctmp ! Temporary condensation rate [kg kg-1 s-1]
395 : real(kind_phys) :: cdt ! Timestep factor [1]
396 : real(kind_phys) :: wtthick ! Layer thickness weight [1]
397 : real(kind_phys) :: pol ! Production over loss ratio [1]
398 :
399 70392 : errmsg = ''
400 70392 : errflg = 0
401 70392 : error_found = .false.
402 :
403 70392 : clrh2o = latvap/rh2o
404 70392 : cpohl = cpair/latvap
405 70392 : hlocp = latvap/cpair
406 70392 : dto2 = 0.5_kind_phys * deltat
407 :
408 : #ifdef PERGRO
409 : mincld = 1.e-4_kind_phys
410 : iter = 1 ! number of times to iterate the precipitation calculation
411 : #else
412 70392 : mincld = 1.e-4_kind_phys
413 70392 : iter = 2
414 : #endif
415 :
416 : ! initialize total cloud water [kg kg-1]
417 28436184 : cwat(:ncol,:pver) = cldice(:ncol,:pver) + cldliq(:ncol,:pver)
418 :
419 : ! initialize single level and multi level fields
420 1090992 : do i = 1, ncol
421 1020600 : precip(i) = 0.0_kind_phys
422 1020600 : precab(i) = 0.0_kind_phys
423 1020600 : snowab(i) = 0.0_kind_phys
424 1090992 : cldmax(i) = 0.0_kind_phys
425 : end do
426 :
427 1900584 : do k = 1, pver
428 28436184 : do i = 1, ncol
429 26535600 : q(i,k) = qn(i,k)
430 28365792 : t(i,k) = tn(i,k)
431 : ! q(i,k)=qn(i,k)-qtend(i,k)*deltat
432 : ! t(i,k)=tn(i,k)-ttend(i,k)*deltat
433 : end do
434 : end do
435 :
436 28436184 : qme (:ncol,:) = 0._kind_phys
437 28436184 : evapprec(:ncol,:) = 0._kind_phys
438 28436184 : prodprec(:ncol,:) = 0._kind_phys
439 28436184 : evapsnow(:ncol,:) = 0._kind_phys
440 28436184 : prodsnow(:ncol,:) = 0._kind_phys
441 28436184 : evapheat(:ncol,:) = 0._kind_phys
442 28436184 : meltheat(:ncol,:) = 0._kind_phys
443 28436184 : prfzheat(:ncol,:) = 0._kind_phys
444 28436184 : ice2pr(:ncol,:) = 0._kind_phys
445 28436184 : liq2pr(:ncol,:) = 0._kind_phys
446 28436184 : liq2snow(:ncol,:) = 0._kind_phys
447 28436184 : fwaut(:ncol,:) = 0._kind_phys
448 28436184 : fsaut(:ncol,:) = 0._kind_phys
449 28436184 : fracw(:ncol,:) = 0._kind_phys
450 28436184 : fsacw(:ncol,:) = 0._kind_phys
451 28436184 : fsaci(:ncol,:) = 0._kind_phys
452 29527176 : lsflxprc(:ncol,:) = 0._kind_phys
453 29527176 : lsflxsnw(:ncol,:) = 0._kind_phys
454 28436184 : pracwo(:ncol,:) = 0._kind_phys
455 28436184 : psacwo(:ncol,:) = 0._kind_phys
456 28436184 : psacio(:ncol,:) = 0._kind_phys
457 :
458 : ! reset diagnostic arrays
459 58632168 : rcwn(:,:,:) = 0._kind_phys
460 58632168 : rliq(:,:,:) = 0._kind_phys
461 58632168 : rice(:,:,:) = 0._kind_phys
462 :
463 : ! find the wet bulb temp and saturation value
464 : ! for the provisional t and q without condensation
465 1900584 : mp_level_loop: do k = top_lev, pver
466 : ! "True" means that ice will be taken into account.
467 0 : call findsp_vc(qn(:ncol,k), tn(:ncol,k), pmid(:ncol,k), .true., &
468 1830192 : tsp(:ncol,k), qsp(:ncol,k))
469 :
470 1830192 : call qsat(t(1:ncol,k), pmid(1:ncol,k), es(1:ncol), qs(1:ncol), ncol, gam=gamma(1:ncol))
471 :
472 28365792 : do i = 1, ncol
473 26535600 : relhum(i) = q(i,k)/qs(i)
474 26535600 : cldm(i) = max(cldn(i,k),mincld)
475 :
476 : ! the max cloud fraction above this level
477 26535600 : cldmax(i) = max(cldmax(i), cldm(i))
478 :
479 : ! define the coefficients for C - E calculation
480 26535600 : calpha(i) = 1.0_kind_phys/qs(i)
481 26535600 : cbeta (i) = q(i,k)/qs(i)**2*gamma(i)*cpohl
482 26535600 : cbetah(i) = 1.0_kind_phys/qs(i)*gamma(i)*cpohl
483 26535600 : cgamma(i) = calpha(i)+latvap*cbeta(i)/cpair
484 26535600 : cgamah(i) = calpha(i)+latvap*cbetah(i)/cpair
485 26535600 : rcgama(i) = cgamma(i)/cgamah(i)
486 :
487 26535600 : if(cldm(i) > mincld) then
488 6182808 : icwc(i) = max(0._kind_phys,cwat(i,k)/cldm(i))
489 : else
490 20352792 : icwc(i) = 0.0_kind_phys
491 : endif
492 : ! PJR the above logic give zero icwc with nonzero cwat, dont like it!
493 : ! PJR generates problems with csigma
494 : ! PJR set the icwc to a very small number, so we can start from zero cloud cover and make some clouds
495 : ! icwc(i) = max(1.e-8_kind_phys,cwat(i,k)/cldm(i))
496 :
497 : ! initial guess of evaporation, will be updated within iteration
498 106142400 : evapprec(i,k) = conke*(1._kind_phys - cldm(i))*sqrt(precab(i)) &
499 132678000 : *(1._kind_phys - min(relhum(i),1._kind_phys))
500 :
501 : ! zero qmeres before iteration for each level
502 28365792 : qmeres(i) = 0.0_kind_phys
503 : end do
504 :
505 28365792 : do i = 1, ncol
506 : ! calculate the cooling due to a phase change of the rainwater from above
507 28365792 : if (t(i,k) >= t0) then
508 4571995 : meltheat(i,k) = -latice * snowab(i) * gravit/pdel(i,k)
509 4571995 : snowab(i) = 0._kind_phys
510 : else
511 21963605 : meltheat(i,k) = 0._kind_phys
512 : endif
513 : end do
514 :
515 : ! calculate qme and formation of precip.
516 : !
517 : ! The cloud microphysics is highly nonlinear and coupled with qme
518 : ! Both rain processes and qme are calculated iteratively.
519 5490576 : qme_iter_loop: do l = 1, iter
520 56731584 : qme_update_loop: do i = 1, ncol
521 : ! calculation of qme has 4 scenarios
522 53071200 : call relhum_min_adj(ncol, pver, tropLev, dlat, rhu00, rhu_adj)
523 :
524 56731584 : if(relhum(i) > rhu_adj(i,k)) then
525 : ! 1. whole grid saturation
526 8991574 : if(relhum(i) >= 0.999_kind_phys .or. cldm(i) >= 0.999_kind_phys) then
527 1325880 : qme(i,k) = (calpha(i)*qtend(i,k)-cbetah(i)*ttend(i,k))/cgamah(i)
528 : ! 2. fractional saturation
529 : else
530 7665694 : if (rhdfda(i,k) .eq. 0._kind_phys .and. icwc(i) .eq. 0._kind_phys) then
531 0 : write(iulog,*) 'prognostic_cloud_water: empty rh cloud @ ', i, k
532 0 : write(iulog,*) 'relhum, iter ', relhum(i), l, rhu_adj(i,k), cldm(i), cldn(i,k)
533 :
534 0 : errflg = 1
535 0 : errmsg = 'prognostic_cloud_water: empty RH cloud'
536 0 : return
537 : endif
538 :
539 7665694 : csigma(i) = 1.0_kind_phys/(rhdfda(i,k)+cgamma(i)*icwc(i))
540 7665694 : qmec1(i) = (1.0_kind_phys-cldm(i))*csigma(i)*rhdfda(i,k)
541 45994164 : qmec2(i) = cldm(i)*calpha(i)/cgamah(i)+(1.0_kind_phys-rcgama(i)*cldm(i))* &
542 53659858 : csigma(i)*calpha(i)*icwc(i)
543 30662776 : qmec3(i) = cldm(i)*cbetah(i)/cgamah(i) + &
544 38328470 : (cbeta(i)-rcgama(i)*cldm(i)*cbetah(i))*csigma(i)*icwc(i)
545 7665694 : qmec4(i) = csigma(i)*cgamma(i)*icwc(i)
546 :
547 : ! Q = C-E = -C1*Al + C2*Aq - C3*At + C4*Er
548 61325552 : qme(i,k) = -qmec1(i)*ltend(i,k) + qmec2(i)*qtend(i,k) &
549 68991246 : -qmec3(i)*ttend(i,k) + qmec4(i)*evapprec(i,k)
550 : endif
551 : ! 3. when rh < rhu00, evaporate existing cloud water
552 44079626 : else if(cwat(i,k) > 0.0_kind_phys) then
553 : ! liquid water should be evaporated but not to exceed
554 : ! saturation point. if qn > qsp, not to evaporate cwat
555 28416630 : qme(i,k) = -min(max(0._kind_phys,qsp(i,k)-qn(i,k)),cwat(i,k))/deltat
556 : ! 4. no condensation nor evaporation
557 : else
558 15662996 : qme(i,k) = 0.0_kind_phys
559 : endif
560 : end do qme_update_loop
561 :
562 : ! Because of the finite time step, place a bound here not to exceed wet bulb point
563 : ! and not to evaporate more than available water
564 56731584 : qme_bound_loop: do i = 1, ncol
565 53071200 : qtmp = qn(i,k) - qme(i,k)*deltat
566 :
567 : ! possibilities to have qtmp > qsp
568 : !
569 : ! 1. if qn > qs(tn), it condenses;
570 : ! if after applying qme, qtmp > qsp, more condensation is applied.
571 : !
572 : ! 2. if qn < qs, evaporation should not exceed qsp,
573 53071200 : if(qtmp > qsp(i,k)) then
574 309736 : qme(i,k) = qme(i,k) + (qtmp-qsp(i,k))/deltat
575 : endif
576 :
577 : ! if net evaporation, it should not exceed available cwat
578 53071200 : if(qme(i,k) < -cwat(i,k)/deltat) then
579 744972 : qme(i,k) = -cwat(i,k)/deltat
580 : endif
581 :
582 : ! addition of residual condensation from previous step of iteration
583 53071200 : qme(i,k) = qme(i,k) + qmeres(i)
584 :
585 : ! limit qme for roundoff errors (multiply by slightly less than unity)
586 56731584 : qme(i,k) = qme(i,k) * omsm
587 : end do qme_bound_loop
588 :
589 56731584 : do i = 1, ncol
590 : ! as a safe limit, condensation should not reduce grid mean rh below rhu00
591 53071200 : if(qme(i,k) > 0.0_kind_phys .and. relhum(i) > rhu_adj(i,k)) then
592 5603377 : qme(i,k) = min(qme(i,k), (qn(i,k)-qs(i)*rhu_adj(i,k))/deltat)
593 : endif
594 :
595 : ! initial guess for cwm (mean cloud water over time step) if 1st iteration
596 56731584 : if(l < 2) then
597 26535600 : cwm(i) = max(cwat(i,k)+qme(i,k)*dto2, 0._kind_phys)
598 : endif
599 : enddo
600 :
601 : ! provisional precipitation falling through model layer
602 56731584 : prprov_update_loop: do i = 1, ncol
603 : ! prprov(i) = precab(i) + prodprec(i,k)*pdel(i,k)/gravit
604 : ! rain produced in this layer not too effective in collection process
605 53071200 : wtthick = max(0._kind_phys,min(0.5_kind_phys,((zi(i,k)-zi(i,k+1))/1000._kind_phys)**2))
606 56731584 : prprov(i) = precab(i) + wtthick*prodprec(i,k)*pdel(i,k)/gravit
607 : end do prprov_update_loop
608 :
609 : ! calculate conversion of condensate to precipitation by cloud microphysics
610 : call findmcnew( &
611 : ncol = ncol, &
612 : pver = pver, &
613 : pi = pi, &
614 : k = k, &
615 0 : precab = prprov(:ncol), &
616 0 : snowab = snowab(:ncol), &
617 0 : t = t(:ncol,:), &
618 0 : p = pmid(:ncol,:), &
619 0 : cwm = cwm(:ncol), &
620 0 : cldm = cldm(:ncol), &
621 0 : cldmax = cldmax(:ncol), &
622 0 : fice = fice(:ncol,k), &
623 0 : landm = landm(:ncol), &
624 0 : seaicef = seaicef(:ncol), &
625 0 : snowh = snowh(:ncol), &
626 : ! below output
627 0 : coef = coef(:ncol), &
628 0 : fwaut = fwaut(:ncol,k), &
629 0 : fsaut = fsaut(:ncol,k), &
630 0 : fracw = fracw(:ncol,k), &
631 0 : fsacw = fsacw(:ncol,k), &
632 0 : fsaci = fsaci(:ncol,k), &
633 0 : pracwo = pracwo(:ncol,k),&
634 0 : psacwo = psacwo(:ncol,k),&
635 3660384 : psacio = psacio(:ncol,k))
636 :
637 : ! calculate the precip rate
638 3660384 : error_found = .false.
639 56731584 : precip_update_loop: do i = 1, ncol
640 53071200 : if (cldm(i) > 0) then
641 : ! first predict the cloud water
642 53071200 : cdt = coef(i)*deltat
643 53071200 : if(cdt > 0.01_kind_phys) then
644 6526573 : pol = qme(i,k)/coef(i) ! production over loss
645 6526573 : cwn(i) = max(0._kind_phys,(cwat(i,k)-pol)*exp(-cdt)+ pol)
646 : else
647 46544627 : cwn(i) = max(0._kind_phys,(cwat(i,k) + qme(i,k)*deltat)/(1+cdt))
648 : endif
649 :
650 : ! now back out the tendency of net rain production
651 53071200 : prodprec(i,k) = max(0._kind_phys,qme(i,k)-(cwn(i)-cwat(i,k))/deltat)
652 : else
653 0 : prodprec(i,k) = 0.0_kind_phys
654 0 : cwn(i) = 0._kind_phys
655 : endif
656 :
657 : ! provisional calculation of conversion terms
658 53071200 : ice2pr(i,k) = prodprec(i,k)*(fsaut(i,k)+fsaci(i,k))
659 53071200 : liq2pr(i,k) = prodprec(i,k)*(fwaut(i,k)+fsacw(i,k)+fracw(i,k))
660 :
661 : ! revision suggested by Jim McCaa
662 : ! it controls the amount of snow hitting the sfc
663 : ! by forcing a lot of conversion of cloud liquid to snow phase
664 : ! it might be better done later by an explicit representation of
665 : ! rain accreting ice (and freezing), or by an explicit freezing of raindrops
666 : !
667 : ! old:
668 : ! liq2snow(i,k) = prodprec(i,k)*fsacw(i,k)
669 53071200 : liq2snow(i,k) = max(prodprec(i,k)*fsacw(i,k), fsnow(i,k)*liq2pr(i,k))
670 :
671 : ! bounds
672 53071200 : nice2pr = min(ice2pr(i,k),(cwat(i,k)+qme(i,k)*deltat)*fice(i,k)/deltat)
673 53071200 : nliq2pr = min(liq2pr(i,k),(cwat(i,k)+qme(i,k)*deltat)*(1._kind_phys-fice(i,k))/deltat)
674 :
675 53071200 : if (liq2pr(i,k) .ne. 0._kind_phys) then
676 13586986 : nliq2snow = liq2snow(i,k)*nliq2pr/liq2pr(i,k) ! correction
677 : else
678 39484214 : nliq2snow = liq2snow(i,k)
679 : endif
680 :
681 : ! avoid roundoff problems generating negatives
682 53071200 : nliq2snow = nliq2snow*omsm
683 53071200 : nliq2pr = nliq2pr*omsm
684 53071200 : nice2pr = nice2pr*omsm
685 :
686 : ! final estimates of conversion to precip and snow
687 53071200 : prodprec(i,k) = (nliq2pr + nice2pr)
688 53071200 : prodsnow(i,k) = (nice2pr + nliq2snow)
689 :
690 : ! compute diagnostics and sanity checks
691 53071200 : rcwn(i,l,k) = cwat(i,k) + (qme(i,k) - prodprec(i,k))*deltat
692 53071200 : rliq(i,l,k) = (cwat(i,k) + qme(i,k)*deltat)*(1._kind_phys-fice(i,k)) - nliq2pr*deltat
693 53071200 : rice(i,l,k) = (cwat(i,k) + qme(i,k)*deltat) * fice(i,k) - nice2pr*deltat
694 :
695 : ! Sanity checks
696 53071200 : if(abs(rcwn(i,l,k)) < 1.e-300_kind_phys) rcwn(i,l,k) = 0._kind_phys
697 53071200 : if(abs(rliq(i,l,k)) < 1.e-300_kind_phys) rliq(i,l,k) = 0._kind_phys
698 53071200 : if(abs(rice(i,l,k)) < 1.e-300_kind_phys) rice(i,l,k) = 0._kind_phys
699 53071200 : if(rcwn(i,l,k) < 0._kind_phys) error_found = .true.
700 53071200 : if(rliq(i,l,k) < 0._kind_phys) error_found = .true.
701 53071200 : if(rice(i,l,k) < 0._kind_phys) error_found = .true.
702 53071200 : if(error_found) then
703 0 : if(rcwn(i,l,k) < 0._kind_phys) then
704 0 : write(iulog,*) 'prognostic_cloud_water: prob with neg rcwn1 ', rcwn(i,l,k), cwn(i)
705 0 : write(iulog,*) ' cwat, qme*deltat, prodprec*deltat ', &
706 0 : cwat(i,k), qme(i,k)*deltat, &
707 0 : prodprec(i,k)*deltat, &
708 0 : (qme(i,k)-prodprec(i,k))*deltat
709 :
710 0 : errflg = 1
711 0 : errmsg = 'prognostic_cloud_water: negative rcwn1'
712 0 : return
713 : endif
714 :
715 0 : if (rliq(i,l,k) < 0._kind_phys) then
716 0 : write(iulog,*) 'prognostic_cloud_water: prob with neg rliq1 ', rliq(i,l,k)
717 0 : errflg = 1
718 0 : errmsg = 'prognostic_cloud_water: negative rliq1'
719 0 : return
720 : endif
721 :
722 0 : if (rice(i,l,k) < 0._kind_phys) then
723 0 : write(iulog,*) 'prognostic_cloud_water: prob with neg rice ', rice(i,l,k)
724 0 : errflg = 1
725 0 : errmsg = 'prognostic_cloud_water: negative rice'
726 0 : return
727 : endif
728 : endif
729 :
730 : ! Final version of condensate to precip terms
731 53071200 : liq2pr(i,k) = nliq2pr
732 53071200 : liq2snow(i,k) = nliq2snow
733 53071200 : ice2pr(i,k) = nice2pr
734 53071200 : cwn(i) = rcwn(i,l,k)
735 :
736 : ! update any remaining provisional values
737 53071200 : cwm(i) = (cwn(i) + cwat(i,k))*0.5_kind_phys
738 :
739 : ! update in cloud water
740 56731584 : if(cldm(i) > mincld) then
741 12365616 : icwc(i) = cwm(i)/cldm(i)
742 : else
743 40705584 : icwc(i) = 0.0_kind_phys
744 : endif
745 : ! PJR the above logic give zero icwc with nonzero cwat, dont like it!
746 : ! PJR generates problems with csigma
747 : ! PJR set the icwc to a very small number, so we can start from zero cloud cover and make some clouds
748 : ! icwc(i) = max(1.e-8_kind_phys,cwm(i)/cldm(i))
749 : end do precip_update_loop
750 :
751 : ! calculate provisional value of cloud water for
752 : ! evaporation of precipitate (evapprec) calculation
753 56731584 : qtl_update_loop: do i = 1, ncol
754 53071200 : qtmp = qn(i,k) - qme(i,k)*deltat
755 53071200 : ttmp = tn(i,k) + deltat/cpair * ( meltheat(i,k) + (latvap + latice*fice(i,k)) * qme(i,k) )
756 53071200 : esn = estblf(ttmp)
757 53071200 : qsn = svp_to_qsat(esn, pmid(i,k))
758 53071200 : qtl(i) = max((qsn - qtmp)/deltat,0._kind_phys)
759 56731584 : relhum1(i) = qtmp/qsn
760 : end do qtl_update_loop
761 :
762 56731584 : evap_update_loop: do i = 1, ncol
763 : #ifdef PERGRO
764 : evapprec(i,k) = conke*(1._kind_phys - max(cldm(i),mincld))* &
765 : sqrt(precab(i))*(1._kind_phys - min(relhum1(i),1._kind_phys))
766 : #else
767 212284800 : evapprec(i,k) = conke*(1._kind_phys - cldm(i))*sqrt(precab(i)) &
768 265356000 : *(1._kind_phys - min(relhum1(i),1._kind_phys))
769 : #endif
770 :
771 : ! limit the evaporation to the amount which is entering the box
772 : ! or saturates the box
773 53071200 : prtmp = precab(i)*gravit/pdel(i,k)
774 53071200 : evapprec(i,k) = min(evapprec(i,k), prtmp, qtl(i))*omsm
775 : #ifdef PERGRO
776 : ! zeroing needed for pert growth
777 : evapprec(i,k) = 0._kind_phys
778 : #endif
779 :
780 : ! Partition evaporation of precipitate between rain and snow using
781 : ! the fraction of snow falling into the box. Determine the heating
782 : ! due to evaporation. Note that evaporation is positive (loss of precip,
783 : ! gain of vapor) and that heating is negative.
784 53071200 : if (evapprec(i,k) > 0._kind_phys) then
785 14490359 : evapsnow(i,k) = evapprec(i,k) * snowab(i) / precab(i)
786 14490359 : evapheat(i,k) = -latvap * evapprec(i,k) - latice * evapsnow(i,k)
787 : else
788 38580841 : evapsnow(i,k) = 0._kind_phys
789 38580841 : evapheat(i,k) = 0._kind_phys
790 : end if
791 :
792 : ! Account for the latent heat of fusion for liquid drops collected by falling snow
793 56731584 : prfzheat(i,k) = latice * liq2snow(i,k)
794 : end do evap_update_loop
795 :
796 : ! now remove the residual of any over-saturation. Normally,
797 : ! the oversaturated water vapor should have been removed by
798 : ! qme formulation plus constraints by wet bulb tsp/qsp
799 : ! as computed above. However, because of non-linearity,
800 : ! addition of (qme-evapprec) to update t and q may still cause
801 : ! a very small amount of over saturation. It is called a
802 : ! residual of over-saturation because theoretically, qme
803 : ! should have taken care of all of large scale condensation.
804 58561776 : qmeres_update_loop: do i = 1, ncol
805 53071200 : qtmp = qn(i,k)-(qme(i,k)-evapprec(i,k))*deltat
806 424569600 : ttmp = tn(i,k) + deltat/cpair * ( meltheat(i,k) + evapheat(i,k) + prfzheat(i,k) &
807 477640800 : + (latvap + latice*fice(i,k)) * qme(i,k) )
808 :
809 53071200 : call qsat(ttmp, pmid(i,k), esn, qsn, dqsdt=dqsdt)
810 :
811 56731584 : if(qtmp > qsn) then
812 : ! now extra condensation to bring air to just saturation
813 90729 : ctmp = (qtmp-qsn)/(1._kind_phys+hlocp*dqsdt)/deltat
814 90729 : qme(i,k) = qme(i,k)+ctmp
815 : ! save residual on qmeres to addtion to qme on entering next iteration
816 : ! qme exit here contain the residual but overrided if back to iteration
817 90729 : qmeres(i) = ctmp
818 : else
819 52980471 : qmeres(i) = 0.0_kind_phys
820 : endif
821 : end do qmeres_update_loop
822 : end do qme_iter_loop ! loop over l (iteration)
823 :
824 : ! precipitation
825 28436184 : precip_snow_loop: do i = 1, ncol
826 26535600 : precip(i) = precip(i) + pdel(i,k)/gravit * (prodprec(i,k) - evapprec(i,k))
827 26535600 : precab(i) = precab(i) + pdel(i,k)/gravit * (prodprec(i,k) - evapprec(i,k))
828 26535600 : if(precab(i) < 0._kind_phys) then
829 0 : precab(i) = 0._kind_phys
830 : endif
831 :
832 : ! snowab(i) = snowab(i) + pdel(i,k)/gravit * (prodprec(i,k)*fice(i,k) - evapsnow(i,k))
833 26535600 : snowab(i) = snowab(i) + pdel(i,k)/gravit * (prodsnow(i,k) - evapsnow(i,k))
834 :
835 : ! If temperature above freezing, all precip is rain flux. if temperature below freezing, all precip is snow flux.
836 26535600 : lsflxprc(i,k+1) = precab(i) !! making this consistent with other precip fluxes. prc = rain + snow
837 : ! lsflxprc(i,k+1) = precab(i) - snowab(i)
838 28365792 : lsflxsnw(i,k+1) = snowab(i)
839 : end do precip_snow_loop
840 : end do mp_level_loop ! loop over k (level)
841 :
842 : ! Convert precip (lwe_stratiform_precipitation_rate_at_surface)
843 : ! and snowab (lwe_snow_precipitation_rate_at_surface_due_to_microphysics)
844 : ! from kg m-2 s-1 to m s-1 (precipitation units) for surface exchange.
845 : !
846 : ! If this conversion is removed in the future, the metadata needs to
847 : ! be updated.
848 1090992 : precip(:ncol) = precip(:ncol)/1000._kind_phys
849 1090992 : snowab(:ncol) = snowab(:ncol)/1000._kind_phys
850 : end subroutine prognostic_cloud_water_run
851 :
852 : ! Calculate the conversion of condensate to precipitate
853 : ! Rasch, P. J., and J. E. Kristjánsson, 1998: A Comparison of the CCM3 Model Climate Using Diagnosed and Predicted Condensate Parameterizations. J. Climate, 11, 1587–1614
854 : ! https://doi.org/10.1175/1520-0442(1998)011<1587:ACOTCM>2.0.CO;2
855 : !
856 : ! Original author: P. Rasch
857 3660384 : subroutine findmcnew( &
858 : ncol, pver, &
859 : pi, &
860 : k, &
861 3660384 : precab, snowab, t, p, cwm, cldm, cldmax, fice, &
862 3660384 : landm, seaicef, snowh, &
863 3660384 : coef, fwaut, fsaut, &
864 3660384 : fracw, fsacw, fsaci, &
865 3660384 : pracwo, psacwo, psacio)
866 :
867 : ! input arguments
868 : integer, intent(in) :: ncol ! number of atmospheric columns
869 : integer, intent(in) :: pver ! number of levels
870 : real(kind_phys), intent(in) :: pi ! pi_constant [1]
871 : integer, intent(in) :: k ! level index
872 :
873 : real(kind_phys), intent(in) :: precab(:) ! rate of precipitation from above [kg m-2 s-1]
874 : real(kind_phys), intent(in) :: t(:,:) ! temperature [K]
875 : real(kind_phys), intent(in) :: p(:,:) ! air_pressure [Pa]
876 : real(kind_phys), intent(in) :: cwm(:) ! condensate mixing ratio [kg kg-1]
877 : real(kind_phys), intent(in) :: cldm(:) ! cloud fraction
878 : real(kind_phys), intent(in) :: cldmax(:) ! max cloud fraction above this level
879 : real(kind_phys), intent(in) :: fice(:) ! fraction of cwat that is ice
880 : real(kind_phys), intent(in) :: landm(:) ! Land fraction ramped over water
881 : real(kind_phys), intent(in) :: seaicef(:) ! sea ice fraction
882 : real(kind_phys), intent(in) :: snowab(:) ! rate of snow from above [kg m-2 s-1]
883 : real(kind_phys), intent(in) :: snowh(:) ! Snow depth over land, water equivalent [m]
884 :
885 : ! output arguments
886 : real(kind_phys), intent(out) :: coef(:) ! conversion rate [s-1]
887 : real(kind_phys), intent(out) :: fwaut(:) ! relative importance of liquid autoconversion (a diagnostic)
888 : real(kind_phys), intent(out) :: fsaut(:) ! relative importance of ice autoconversion (a diagnostic)
889 : real(kind_phys), intent(out) :: fracw(:) ! relative importance of rain accreting liquid (a diagnostic)
890 : real(kind_phys), intent(out) :: fsacw(:) ! relative importance of snow accreting liquid (a diagnostic)
891 : real(kind_phys), intent(out) :: fsaci(:) ! relative importance of snow accreting ice (a diagnostic)
892 : real(kind_phys), intent(out) :: pracwo(:) ! accretion of cloud water by rain [s-1]
893 : real(kind_phys), intent(out) :: psacwo(:) ! accretion of cloud water by snow [s-1]
894 : real(kind_phys), intent(out) :: psacio(:) ! accretion of cloud ice by snow [s-1]
895 :
896 : ! local variables
897 : integer :: i, ii ! Loop index [index]
898 : integer :: ncols ! Number of active columns for microphysics (different from ncol!!) [count]
899 7320768 : integer :: ind(ncol) ! Active column indices [index]
900 : real(kind_phys) :: capn ! Local cloud particle number concentration [cm-3]
901 : real(kind_phys) :: capnoice ! Cloud particle concentration excluding sea ice [cm-3]
902 7320768 : real(kind_phys) :: cldloc(ncol) ! Non-zero cloud fraction [1]
903 7320768 : real(kind_phys) :: cldpr(ncol) ! Cloud fraction for precipitation [1]
904 7320768 : real(kind_phys) :: totmr(ncol) ! In-cloud total water mixing ratio [kg kg-1]
905 7320768 : real(kind_phys) :: icemr(ncol) ! In-cloud ice mixing ratio [kg kg-1]
906 7320768 : real(kind_phys) :: liqmr(ncol) ! In-cloud liquid water mixing ratio [kg kg-1]
907 7320768 : real(kind_phys) :: rainmr(ncol) ! In-cloud rain mixing ratio [kg kg-1]
908 : real(kind_phys) :: ciaut ! Ice autoconversion coefficient [s-1]
909 : real(kind_phys) :: pracw ! Rate of rain collecting cloud water [s-1]
910 : real(kind_phys) :: psaci ! Rate of snow collecting cloud ice [s-1]
911 : real(kind_phys) :: psacw ! Rate of snow collecting cloud water [s-1]
912 : real(kind_phys) :: psaut ! Rate of ice autoconversion [s-1]
913 : real(kind_phys) :: pwaut ! Rate of liquid water autoconversion [s-1]
914 : real(kind_phys) :: ptot ! Total conversion rate [s-1]
915 7320768 : real(kind_phys) :: prlloc(ncol) ! Local rain flux [mm day-1]
916 7320768 : real(kind_phys) :: prscgs(ncol) ! Local snow amount [g cm-2]
917 : real(kind_phys) :: snowfr ! Snow fraction of precipitation [1]
918 : real(kind_phys) :: vfallw ! Fall speed of liquid precipitation [m s-1]
919 7320768 : real(kind_phys) :: rho(ncol) ! Air density [kg m-3]
920 : real(kind_phys) :: rhocgs ! Air density in CGS units [g cm-3]
921 : real(kind_phys) :: r3l ! Cloud droplet volume radius [m]
922 : real(kind_phys) :: icrit ! Ice autoconversion threshold [kg kg-1]
923 : real(kind_phys) :: wsi ! Sea ice weight factor [1]
924 : real(kind_phys) :: wt ! Ice fraction weight [1]
925 : real(kind_phys) :: wland ! Land fraction weight [1]
926 : real(kind_phys) :: wp ! Pressure dependence weight [1]
927 : real(kind_phys) :: ftot ! Total fraction for conversion processes [1]
928 : real(kind_phys) :: con1 ! Work constant for radius calculation [m]
929 : real(kind_phys) :: con2 ! Work constant for density ratios [1]
930 : real(kind_phys) :: csacx ! Constant used for snow accreting liquid or ice [??]
931 : real(kind_phys) :: rat1 ! Density ratio work variable [1]
932 : real(kind_phys) :: rat2 ! Mixing ratio ratio work variable [1]
933 :
934 : ! find all the points where we need to do the microphysics
935 : ! and set the output variables to zero
936 3660384 : ncols = 0
937 56731584 : do i = 1, ncol
938 53071200 : coef(i) = 0._kind_phys
939 53071200 : fwaut(i) = 0._kind_phys
940 53071200 : fsaut(i) = 0._kind_phys
941 53071200 : fracw(i) = 0._kind_phys
942 53071200 : fsacw(i) = 0._kind_phys
943 53071200 : fsaci(i) = 0._kind_phys
944 53071200 : liqmr(i) = 0._kind_phys
945 53071200 : rainmr(i) = 0._kind_phys
946 56731584 : if (cwm(i) > 1.e-20_kind_phys) then
947 27471114 : ncols = ncols + 1
948 27471114 : ind(ncols) = i
949 : endif
950 : end do
951 :
952 31131498 : do ii = 1, ncols
953 : ! map from active microphysics columns to physics columns
954 27471114 : i = ind(ii)
955 :
956 : ! the local cloudiness at this level
957 27471114 : cldloc(i) = max(cldmin,cldm(i))
958 :
959 : ! a weighted mean between max cloudiness above, and this layer
960 27471114 : cldpr(i) = max(cldmin,(cldmax(i)+cldm(i))*0.5_kind_phys)
961 :
962 : ! decompose the suspended condensate into
963 : ! an incloud liquid and ice phase component
964 27471114 : totmr(i) = cwm(i)/cldloc(i)
965 27471114 : icemr(i) = totmr(i)*fice(i)
966 27471114 : liqmr(i) = totmr(i)*(1._kind_phys-fice(i))
967 :
968 : ! density
969 27471114 : rho(i) = p(i,k)/(287._kind_phys*t(i,k))
970 27471114 : rhocgs = rho(i)*1.e-3_kind_phys ! density in cgs units
971 :
972 : ! decompose the precipitate into a liquid and ice phase
973 27471114 : if (t(i,k) > t0) then
974 8255564 : vfallw = convfw/sqrt(rho(i))
975 8255564 : rainmr(i) = precab(i)/(rho(i)*vfallw*cldpr(i))
976 8255564 : snowfr = 0
977 : else
978 19215550 : snowfr = 1
979 19215550 : rainmr(i) = 0._kind_phys
980 : endif
981 :
982 : ! local snow amount in cgs units
983 27471114 : prscgs(i) = precab(i)/cldpr(i)*0.1_kind_phys*snowfr
984 :
985 : ! local rain amount in mm/day
986 31131498 : prlloc(i) = precab(i)*86400._kind_phys/cldpr(i)
987 : end do
988 :
989 3660384 : con1 = 1._kind_phys/(1.333_kind_phys*pi)**0.333_kind_phys * 0.01_kind_phys ! meters
990 :
991 : ! calculate the conversion terms
992 31131498 : do ii = 1, ncols
993 27471114 : i = ind(ii)
994 27471114 : rhocgs = rho(i)*1.e-3_kind_phys ! density in cgs units
995 :
996 : ! some temperature dependent constants
997 27471114 : wt = fice(i)
998 27471114 : icrit = icritc*wt + icritw*(1-wt)
999 :
1000 : ! jrm Reworked droplet number concentration algorithm
1001 : ! Start with pressure-dependent value appropriate for continental air
1002 : ! Note: reltab has a temperature dependence here
1003 27471114 : capn = capnw + (capnc-capnw) * min(1._kind_phys,max(0._kind_phys,1.0_kind_phys-(p(i,k)-0.8_kind_phys*p(i,pver))/(0.2_kind_phys*p(i,pver))))
1004 : ! Modify for snow depth over land
1005 27471114 : capn = capn + (capnc-capn) * min(1.0_kind_phys,max(0.0_kind_phys,snowh(i)*10._kind_phys))
1006 : ! Ramp between polluted value over land to clean value over ocean.
1007 27471114 : capn = capn + (capnc-capn) * min(1.0_kind_phys,max(0.0_kind_phys,1.0_kind_phys-landm(i)))
1008 : ! Ramp between the resultant value and a sea ice value in the presence of ice.
1009 27471114 : capn = capn + (capnsi-capn) * min(1.0_kind_phys,max(0.0_kind_phys,seaicef(i)))
1010 : ! end jrm
1011 :
1012 : ! useful terms in following calculations
1013 27471114 : rat1 = rhocgs/rhow
1014 27471114 : rat2 = liqmr(i)/capn
1015 27471114 : con2 = (rat1*rat2)**0.333_kind_phys
1016 :
1017 : ! volume radius
1018 27471114 : r3l = con1*con2
1019 :
1020 : ! critical threshold for autoconversion if modified for mixed phase
1021 : ! clouds to mimic a bergeron findeisen process
1022 : ! r3lc2 = r3lcrit*(1.-0.5*fice(i)*(1-fice(i)))
1023 : !
1024 : ! autoconversion of liquid
1025 : !
1026 : ! cwaut = 2.e-4
1027 : ! cwaut = 1.e-3
1028 : ! lcrit = 2.e-4
1029 : ! lcrit = 5.e-4
1030 : ! pwaut = max(0._kind_phys,liqmr(i)-lcrit)*cwaut
1031 : !
1032 : ! pwaut is following tripoli and cotton (and many others)
1033 : ! we reduce the autoconversion below critpr, because these are regions where
1034 : ! the drop size distribution is likely to imply much smaller collector drops than
1035 : ! those relevant for a cloud distribution corresponding to the value of effc = 0.55
1036 : ! suggested by cotton (see austin 1995 JAS, baker 1993)
1037 :
1038 : ! easy to follow form
1039 : ! pwaut = capc*liqmr(i)**2*rhocgs/rhow
1040 : ! *(liqmr(i)*rhocgs/(rhow*capn))**(.333)
1041 : ! *heavy(r3l,r3lcrit)
1042 : ! *max(0.10_kind_phys,min(1._kind_phys,prlloc(i)/critpr))
1043 : ! somewhat faster form
1044 :
1045 : ! using modified Heaviside function:
1046 27471114 : pwaut = capc*liqmr(i)**2*rat1*con2*heavymp(r3l,r3lcrit) * &
1047 54942228 : max(0.10_kind_phys,min(1._kind_phys,prlloc(i)/critpr))
1048 :
1049 : ! not using modified Heaviside function:
1050 : ! pwaut = capc*liqmr(i)**2*rat1*con2*heavym(r3l,r3lcrit)* &
1051 : ! max(0.10_kind_phys,min(1._kind_phys,prlloc(i)/critpr))
1052 :
1053 : ! autoconversion of ice
1054 27471114 : ciaut = ciautb
1055 :
1056 : ! autoconversion of ice condensate
1057 : #ifdef PERGRO
1058 : psaut = heavyp(icemr(i),icrit)*icemr(i)*ciaut
1059 : #else
1060 27471114 : psaut = max(0._kind_phys,icemr(i)-icrit)*ciaut
1061 : #endif
1062 :
1063 : ! collection of liquid by rain
1064 : ! pracw = cracw*rho(i)*liqmr(i)*rainmr(i) !(beheng 1994)
1065 27471114 : pracw = cracw*rho(i)*sqrt(rho(i))*liqmr(i)*rainmr(i) !(tripoli and cotton)
1066 27471114 : pracwo(i) = pracw
1067 :
1068 : ! the following lines calculate the slope parameter and snow mixing ratio
1069 : ! from the precip rate using the equations found in lin et al 83
1070 : ! in the most natural form, but it is expensive, so after some tedious
1071 : ! algebraic manipulation you can use the cheaper form found below
1072 : ! vfalls = c*gam4pd/(6*lamdas**d)*sqrt(rhonot/rhocgs)
1073 : ! $ *0.01 ! convert from cm/s to m/s
1074 : ! snowmr(i) = snowfr*precab(i)/(rho(i)*vfalls*cldpr(i))
1075 : ! snowmr(i) = ( prscgs(i)*mcon02 * (rhocgs**mcon03) )**mcon04
1076 : ! lamdas = (prhonos/max(rhocgs*snowmr(i),small))**0.25
1077 : ! csacw = mcon01*sqrt(rhonot/rhocgs)/(lamdas**thrpd)
1078 : !
1079 : ! coefficient for collection by snow independent of phase
1080 27471114 : csacx = mcon07*rhocgs**mcon08*prscgs(i)**mcon05
1081 :
1082 : ! collection of liquid by snow (lin et al 1983)
1083 27471114 : psacw = csacx*liqmr(i)*esw
1084 :
1085 : #ifdef PERGRO
1086 : ! this is necessary for pergro
1087 : psacw = 0._kind_phys
1088 : #endif
1089 :
1090 27471114 : psacwo(i) = psacw
1091 :
1092 : ! collection of ice by snow (lin et al 1983)
1093 27471114 : psaci = csacx*icemr(i)*esi
1094 27471114 : psacio(i) = psaci
1095 :
1096 : ! total conversion of condensate to precipitate
1097 27471114 : ptot = pwaut + psaut + pracw + psacw + psaci
1098 :
1099 : ! the reciprocal of cloud water amount (or zero if no cloud water)
1100 : ! rcwm = totmr(i)/(max(totmr(i),small)**2)
1101 :
1102 : ! turn the tendency back into a loss rate [s-1]
1103 27471114 : if (totmr(i) > 0._kind_phys) then
1104 27471114 : coef(i) = ptot/totmr(i)
1105 : else
1106 0 : coef(i) = 0._kind_phys
1107 : endif
1108 :
1109 27471114 : if (ptot.gt.0._kind_phys) then
1110 18508370 : fwaut(i) = pwaut/ptot
1111 18508370 : fsaut(i) = psaut/ptot
1112 18508370 : fracw(i) = pracw/ptot
1113 18508370 : fsacw(i) = psacw/ptot
1114 18508370 : fsaci(i) = psaci/ptot
1115 : else
1116 8962744 : fwaut(i) = 0._kind_phys
1117 8962744 : fsaut(i) = 0._kind_phys
1118 8962744 : fracw(i) = 0._kind_phys
1119 8962744 : fsacw(i) = 0._kind_phys
1120 8962744 : fsaci(i) = 0._kind_phys
1121 : endif
1122 :
1123 31131498 : ftot = fwaut(i)+fsaut(i)+fracw(i)+fsacw(i)+fsaci(i)
1124 : end do
1125 :
1126 3660384 : end subroutine findmcnew
1127 :
1128 : ! For special WACCM/CAM-chem cases with CAM4 physics:
1129 : ! Sets rhu to a different value poleward of +/- 50 deg latitude and
1130 : ! levels above the tropopause if polstrat_rhmin is specified
1131 53071200 : subroutine relhum_min_adj(ncol, pver, tropLev, dlat, rhu, rhu_adj)
1132 : integer, intent(in) :: ncol
1133 : integer, intent(in) :: pver
1134 : integer, intent(in) :: tropLev(:)
1135 : real(kind_phys), intent(in) :: dlat(:) ! latitudes in degrees
1136 : real(kind_phys), intent(in) :: rhu(:,:)
1137 : real(kind_phys), intent(out) :: rhu_adj(:,:)
1138 :
1139 : integer :: i, k
1140 :
1141 22805284320 : rhu_adj(:,:) = rhu(:,:)
1142 53071200 : if (.not. do_psrhmin) return
1143 :
1144 0 : do k = 1, pver
1145 0 : do i = 1, ncol
1146 : ! assumption that up is lower index
1147 0 : if ((k .lt. tropLev(i)) .and. (abs(dlat(i)) .gt. 50._kind_phys)) then
1148 0 : rhu_adj(i,k) = psrhmin
1149 : endif
1150 : enddo
1151 : enddo
1152 : end subroutine relhum_min_adj
1153 :
1154 : ! Heavyside step functions used in findmcnew
1155 : ! Heavyside step function
1156 0 : pure function heavy(a1, a2) result(h)
1157 : ! h: 1 if a1 > a2, 0 otherwise
1158 : real(kind_phys), intent(in) :: a1, a2
1159 : real(kind_phys) :: h
1160 0 : h = max(0.0_kind_phys, sign(1.0_kind_phys, a1-a2))
1161 0 : end function heavy
1162 :
1163 : ! Modified Heavyside function with minimum value of 0.01
1164 0 : pure function heavym(a1, a2) result(h)
1165 : ! h: 1 if a1 > a2, 0.01 otherwise
1166 : real(kind_phys), intent(in) :: a1, a2
1167 : real(kind_phys) :: h
1168 0 : h = max(0.01_kind_phys, sign(1.0_kind_phys, a1-a2))
1169 0 : end function heavym
1170 :
1171 : ! Smoothed Heavyside function using ratio formulation
1172 0 : pure function heavyp(a1, a2) result(h)
1173 : ! h: Ratio of a1/(a1+a2+eps)
1174 : real(kind_phys), intent(in) :: a1, a2
1175 : real(kind_phys) :: h
1176 : real(kind_phys), parameter :: eps = 1.0e-36_kind_phys
1177 0 : h = a1/(a2 + a1 + eps)
1178 0 : end function heavyp
1179 :
1180 : ! Modified smooth Heavyside with offset term
1181 27471114 : pure function heavymp(a1, a2) result(h)
1182 : ! h: Ratio (a1 + 0.01*a2)/(a1+a2+eps)
1183 : real(kind_phys), intent(in) :: a1, a2
1184 : real(kind_phys) :: h
1185 : real(kind_phys), parameter :: eps = 1.0e-36_kind_phys
1186 27471114 : h = (a1 + 0.01_kind_phys*a2)/(a2 + a1 + eps)
1187 27471114 : end function heavymp
1188 :
1189 : end module prognostic_cloud_water
|