Line data Source code
1 : module macrop_driver
2 :
3 : !-------------------------------------------------------------------------------------------------------
4 : ! Purpose:
5 : !
6 : ! Provides the CAM interface to the prognostic cloud macrophysics
7 : !
8 : ! Author: Andrew Gettelman, Cheryl Craig October 2010
9 : ! Origin: modified from stratiform.F90 elements
10 : ! (Boville 2002, Coleman 2004, Park 2009, Kay 2010)
11 : !-------------------------------------------------------------------------------------------------------
12 :
13 : use shr_kind_mod, only: r8=>shr_kind_r8
14 : use spmd_utils, only: masterproc
15 : use ppgrid, only: pcols, pver, pverp
16 : use physconst, only: latice, latvap
17 : use phys_control, only: phys_getopts
18 : use constituents, only: cnst_get_ind, pcnst
19 : use physics_buffer, only: physics_buffer_desc, pbuf_set_field, pbuf_get_field, pbuf_old_tim_idx
20 : use time_manager, only: is_first_step
21 : use cldwat2m_macro, only: ini_macro
22 : use perf_mod, only: t_startf, t_stopf
23 : use cam_logfile, only: iulog
24 : use cam_abortutils, only: endrun
25 :
26 : implicit none
27 : private
28 : save
29 :
30 : public :: macrop_driver_readnl
31 : public :: macrop_driver_register
32 : public :: macrop_driver_init
33 : public :: macrop_driver_tend
34 : public :: liquid_macro_tend
35 :
36 : logical, public :: do_cldice ! .true., park macrophysics is prognosing cldice
37 : logical, public :: do_cldliq ! .true., park macrophysics is prognosing cldliq
38 : logical, public :: do_detrain ! .true., park macrophysics is detraining ice into stratiform
39 :
40 : ! ------------------------- !
41 : ! Private Module Parameters !
42 : ! ------------------------- !
43 :
44 : ! 'cu_det_st' : If .true. (.false.), detrain cumulus liquid condensate into the pre-existing liquid stratus
45 : ! (environment) without (with) macrophysical evaporation. If there is no pre-esisting stratus,
46 : ! evaporate cumulus liquid condensate. This option only influences the treatment of cumulus
47 : ! liquid condensate, not cumulus ice condensate.
48 :
49 : logical, parameter :: cu_det_st = .false.
50 :
51 : ! Parameters used for selecting generalized critical RH for liquid and ice stratus
52 : integer :: rhminl_opt = 0
53 : integer :: rhmini_opt = 0
54 :
55 :
56 : character(len=16) :: shallow_scheme
57 : logical :: use_shfrc ! Local copy of flag from convect_shallow_use_shfrc
58 :
59 : integer :: &
60 : ixcldliq, &! cloud liquid amount index
61 : ixcldice, &! cloud ice amount index
62 : ixnumliq, &! cloud liquid number index
63 : ixnumice, &! cloud ice water index
64 : qcwat_idx, &! qcwat index in physics buffer
65 : lcwat_idx, &! lcwat index in physics buffer
66 : iccwat_idx, &! iccwat index in physics buffer
67 : nlwat_idx, &! nlwat index in physics buffer
68 : niwat_idx, &! niwat index in physics buffer
69 : tcwat_idx, &! tcwat index in physics buffer
70 : CC_T_idx, &!
71 : CC_qv_idx, &!
72 : CC_ql_idx, &!
73 : CC_qi_idx, &!
74 : CC_nl_idx, &!
75 : CC_ni_idx, &!
76 : CC_qlst_idx, &!
77 : cld_idx, &! cld index in physics buffer
78 : ast_idx, &! stratiform cloud fraction index in physics buffer
79 : aist_idx, &! ice stratiform cloud fraction index in physics buffer
80 : alst_idx, &! liquid stratiform cloud fraction index in physics buffer
81 : qist_idx, &! ice stratiform in-cloud IWC
82 : qlst_idx, &! liquid stratiform in-cloud LWC
83 : concld_idx, &! concld index in physics buffer
84 : fice_idx, &
85 : cmeliq_idx, &
86 : shfrc_idx
87 :
88 : ! Physics buffer indices for convective_cloud_cover
89 : integer :: sh_frac_idx = 0
90 : integer :: dp_frac_idx = 0
91 :
92 : integer :: &
93 : dlfzm_idx = -1, & ! ZM detrained convective cloud water mixing ratio.
94 : dnlfzm_idx = -1, & ! ZM detrained convective cloud water num concen.
95 : dnifzm_idx = -1 ! ZM detrained convective cloud ice num concen.
96 :
97 :
98 : integer :: &
99 : tke_idx = -1, &! tke defined at the model interfaces
100 : qtl_flx_idx = -1, &! overbar(w'qtl' where qtl = qv + ql) from the PBL scheme
101 : qti_flx_idx = -1, &! overbar(w'qti' where qti = qv + qi) from the PBL scheme
102 : cmfr_det_idx = -1, &! detrained convective mass flux from UNICON
103 : qlr_det_idx = -1, &! detrained convective ql from UNICON
104 : qir_det_idx = -1, &! detrained convective qi from UNICON
105 : cmfmc_sh_idx = -1
106 :
107 : contains
108 :
109 : ! ===============================================================================
110 1024 : subroutine macrop_driver_readnl(nlfile)
111 :
112 : use namelist_utils, only: find_group_name
113 : use units, only: getunit, freeunit
114 : use mpishorthand
115 :
116 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
117 :
118 : ! Namelist variables
119 : logical :: macro_park_do_cldice = .true. ! do_cldice = .true., park macrophysics is prognosing cldice
120 : logical :: macro_park_do_cldliq = .true. ! do_cldliq = .true., park macrophysics is prognosing cldliq
121 : logical :: macro_park_do_detrain = .true. ! do_detrain = .true., park macrophysics is detraining ice into stratiform
122 :
123 : ! Local variables
124 : integer :: unitn, ierr
125 : character(len=*), parameter :: subname = 'macrop_driver_readnl'
126 :
127 : namelist /macro_park_nl/ macro_park_do_cldice, macro_park_do_cldliq, macro_park_do_detrain
128 : !-----------------------------------------------------------------------------
129 :
130 1026 : if (masterproc) then
131 2 : unitn = getunit()
132 2 : open( unitn, file=trim(nlfile), status='old' )
133 2 : call find_group_name(unitn, 'macro_park_nl', status=ierr)
134 2 : if (ierr == 0) then
135 0 : read(unitn, macro_park_nl, iostat=ierr)
136 0 : if (ierr /= 0) then
137 0 : call endrun(subname // ':: ERROR reading namelist')
138 : end if
139 : end if
140 2 : close(unitn)
141 2 : call freeunit(unitn)
142 :
143 : ! set local variables
144 :
145 2 : do_cldice = macro_park_do_cldice
146 2 : do_cldliq = macro_park_do_cldliq
147 2 : do_detrain = macro_park_do_detrain
148 :
149 : end if
150 :
151 : #ifdef SPMD
152 : ! Broadcast namelist variables
153 1024 : call mpibcast(do_cldice, 1, mpilog, 0, mpicom)
154 1024 : call mpibcast(do_cldliq, 1, mpilog, 0, mpicom)
155 1024 : call mpibcast(do_detrain, 1, mpilog, 0, mpicom)
156 : #endif
157 :
158 1024 : end subroutine macrop_driver_readnl
159 :
160 : !================================================================================================
161 :
162 0 : subroutine macrop_driver_register
163 :
164 : !---------------------------------------------------------------------- !
165 : ! !
166 : ! Register the constituents (cloud liquid and cloud ice) and the fields !
167 : ! in the physics buffer. !
168 : ! !
169 : !---------------------------------------------------------------------- !
170 :
171 :
172 : use physics_buffer, only : pbuf_add_field, dtype_r8, dyn_time_lvls
173 :
174 : !-----------------------------------------------------------------------
175 :
176 0 : call phys_getopts(shallow_scheme_out=shallow_scheme)
177 :
178 0 : call pbuf_add_field('AST', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), ast_idx)
179 0 : call pbuf_add_field('AIST', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), aist_idx)
180 0 : call pbuf_add_field('ALST', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), alst_idx)
181 0 : call pbuf_add_field('QIST', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), qist_idx)
182 0 : call pbuf_add_field('QLST', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), qlst_idx)
183 0 : call pbuf_add_field('CLD', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), cld_idx)
184 0 : call pbuf_add_field('CONCLD', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), concld_idx)
185 :
186 0 : call pbuf_add_field('QCWAT', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), qcwat_idx)
187 0 : call pbuf_add_field('LCWAT', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), lcwat_idx)
188 0 : call pbuf_add_field('ICCWAT', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), iccwat_idx)
189 0 : call pbuf_add_field('NLWAT', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), nlwat_idx)
190 0 : call pbuf_add_field('NIWAT', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), niwat_idx)
191 0 : call pbuf_add_field('TCWAT', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), tcwat_idx)
192 :
193 0 : call pbuf_add_field('FICE', 'physpkg', dtype_r8, (/pcols,pver/), fice_idx)
194 :
195 0 : call pbuf_add_field('CMELIQ', 'physpkg', dtype_r8, (/pcols,pver/), cmeliq_idx)
196 :
197 0 : end subroutine macrop_driver_register
198 :
199 : !============================================================================ !
200 : ! !
201 : !============================================================================ !
202 :
203 0 : subroutine macrop_driver_init(pbuf2d)
204 :
205 : !-------------------------------------------- !
206 : ! !
207 : ! Initialize the cloud water parameterization !
208 : ! !
209 : !-------------------------------------------- !
210 0 : use physics_buffer, only : pbuf_get_index
211 : use cam_history, only: addfld, add_default
212 : use convect_shallow, only: convect_shallow_use_shfrc
213 :
214 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
215 :
216 : logical :: history_aerosol ! Output the MAM aerosol tendencies
217 : logical :: history_budget ! Output tendencies and state variables for CAM4
218 : ! temperature, water vapor, cloud ice and cloud
219 : ! liquid budgets.
220 : integer :: history_budget_histfile_num ! output history file number for budget fields
221 : integer :: istat
222 :
223 : character(len=*), parameter :: subname = 'macrop_driver_init'
224 : !-----------------------------------------------------------------------
225 :
226 : ! Initialization routine for cloud macrophysics
227 0 : if (shallow_scheme .eq. 'UNICON') rhminl_opt = 1
228 0 : call ini_macro(rhminl_opt, rhmini_opt)
229 :
230 : call phys_getopts(history_aerosol_out = history_aerosol , &
231 : history_budget_out = history_budget , &
232 0 : history_budget_histfile_num_out = history_budget_histfile_num )
233 :
234 : ! Find out whether shfrc from convect_shallow will be used in cldfrc
235 :
236 0 : if( convect_shallow_use_shfrc() ) then
237 0 : use_shfrc = .true.
238 0 : shfrc_idx = pbuf_get_index('shfrc')
239 : else
240 0 : use_shfrc = .false.
241 : endif
242 :
243 0 : call addfld ('DPDLFLIQ', (/ 'lev' /), 'A', 'kg/kg/s', 'Detrained liquid water from deep convection', sampled_on_subcycle=.true.)
244 0 : call addfld ('DPDLFICE', (/ 'lev' /), 'A', 'kg/kg/s', 'Detrained ice from deep convection', sampled_on_subcycle=.true.)
245 0 : call addfld ('SHDLFLIQ', (/ 'lev' /), 'A', 'kg/kg/s', 'Detrained liquid water from shallow convection', sampled_on_subcycle=.true.)
246 0 : call addfld ('SHDLFICE', (/ 'lev' /), 'A', 'kg/kg/s', 'Detrained ice from shallow convection', sampled_on_subcycle=.true.)
247 0 : call addfld ('DPDLFT', (/ 'lev' /), 'A', 'K/s', 'T-tendency due to deep convective detrainment', sampled_on_subcycle=.true.)
248 0 : call addfld ('SHDLFT', (/ 'lev' /), 'A', 'K/s', 'T-tendency due to shallow convective detrainment', sampled_on_subcycle=.true.)
249 :
250 0 : call addfld ('ZMDLF', (/ 'lev' /), 'A', 'kg/kg/s', 'Detrained liquid water from ZM convection', sampled_on_subcycle=.true.)
251 :
252 0 : call addfld ('MACPDT', (/ 'lev' /), 'A', 'W/kg', 'Heating tendency - Revised macrophysics', sampled_on_subcycle=.true.)
253 0 : call addfld ('MACPDQ', (/ 'lev' /), 'A', 'kg/kg/s', 'Q tendency - Revised macrophysics', sampled_on_subcycle=.true.)
254 0 : call addfld ('MACPDLIQ', (/ 'lev' /), 'A', 'kg/kg/s', 'CLDLIQ tendency - Revised macrophysics', sampled_on_subcycle=.true.)
255 0 : call addfld ('MACPDICE', (/ 'lev' /), 'A', 'kg/kg/s', 'CLDICE tendency - Revised macrophysics', sampled_on_subcycle=.true.)
256 :
257 : call addfld ('CLDVAPADJ', (/ 'lev' /), 'A', 'kg/kg/s', &
258 0 : 'Q tendency associated with liq/ice adjustment - Revised macrophysics', sampled_on_subcycle=.true.)
259 0 : call addfld ('CLDLIQADJ', (/ 'lev' /), 'A', 'kg/kg/s', 'CLDLIQ adjustment tendency - Revised macrophysics', sampled_on_subcycle=.true.)
260 0 : call addfld ('CLDICEADJ', (/ 'lev' /), 'A', 'kg/kg/s', 'CLDICE adjustment tendency - Revised macrophysics', sampled_on_subcycle=.true.)
261 : call addfld ('CLDLIQDET', (/ 'lev' /), 'A', 'kg/kg/s', &
262 0 : 'Detrainment of conv cld liq into envrionment - Revised macrophysics', sampled_on_subcycle=.true.)
263 : call addfld ('CLDICEDET', (/ 'lev' /), 'A', 'kg/kg/s', &
264 0 : 'Detrainment of conv cld ice into envrionment - Revised macrophysics', sampled_on_subcycle=.true.)
265 0 : call addfld ('CLDLIQLIM', (/ 'lev' /), 'A', 'kg/kg/s', 'CLDLIQ limiting tendency - Revised macrophysics', sampled_on_subcycle=.true.)
266 0 : call addfld ('CLDICELIM', (/ 'lev' /), 'A', 'kg/kg/s', 'CLDICE limiting tendency - Revised macrophysics', sampled_on_subcycle=.true.)
267 :
268 0 : call addfld ('AST', (/ 'lev' /), 'A', '1', 'Stratus cloud fraction', sampled_on_subcycle=.true.)
269 0 : call addfld ('LIQCLDF', (/ 'lev' /), 'A', '1', 'Stratus Liquid cloud fraction', sampled_on_subcycle=.true.)
270 0 : call addfld ('ICECLDF', (/ 'lev' /), 'A', '1', 'Stratus ICE cloud fraction', sampled_on_subcycle=.true.)
271 :
272 0 : call addfld ('CLDST', (/ 'lev' /), 'A', 'fraction', 'Stratus cloud fraction', sampled_on_subcycle=.true.)
273 0 : call addfld ('CONCLD', (/ 'lev' /), 'A', 'fraction', 'Convective cloud cover', sampled_on_subcycle=.true.)
274 :
275 0 : call addfld ('CLR_LIQ', (/ 'lev' /), 'A', 'fraction', 'Clear sky fraction for liquid stratus', sampled_on_subcycle=.true.)
276 0 : call addfld ('CLR_ICE', (/ 'lev' /), 'A', 'fraction', 'Clear sky fraction for ice stratus', sampled_on_subcycle=.true.)
277 :
278 0 : call addfld ('CLDLIQSTR', (/ 'lev' /), 'A', 'kg/kg', 'Stratiform CLDLIQ', sampled_on_subcycle=.true.)
279 0 : call addfld ('CLDICESTR', (/ 'lev' /), 'A', 'kg/kg', 'Stratiform CLDICE', sampled_on_subcycle=.true.)
280 0 : call addfld ('CLDLIQCON', (/ 'lev' /), 'A', 'kg/kg', 'Convective CLDLIQ', sampled_on_subcycle=.true.)
281 0 : call addfld ('CLDICECON', (/ 'lev' /), 'A', 'kg/kg', 'Convective CLDICE', sampled_on_subcycle=.true.)
282 :
283 0 : call addfld ('CLDSICE', (/ 'lev' /), 'A', 'kg/kg', 'CloudSat equivalent ice mass mixing ratio', sampled_on_subcycle=.true.)
284 0 : call addfld ('CMELIQ', (/ 'lev' /), 'A', 'kg/kg/s', 'Rate of cond-evap of liq within the cloud', sampled_on_subcycle=.true.)
285 :
286 0 : call addfld ('TTENDICE', (/ 'lev' /), 'A', 'K/s', 'T tendency from Ice Saturation Adjustment', sampled_on_subcycle=.true.)
287 0 : call addfld ('QVTENDICE', (/ 'lev' /), 'A', 'kg/kg/s', 'Q tendency from Ice Saturation Adjustment', sampled_on_subcycle=.true.)
288 0 : call addfld ('QITENDICE', (/ 'lev' /), 'A', 'kg/kg/s', 'CLDICE tendency from Ice Saturation Adjustment', sampled_on_subcycle=.true.)
289 0 : call addfld ('NITENDICE', (/ 'lev' /), 'A', 'kg/kg/s', 'NUMICE tendency from Ice Saturation Adjustment', sampled_on_subcycle=.true.)
290 0 : if ( history_budget ) then
291 :
292 0 : call add_default ('DPDLFLIQ ', history_budget_histfile_num, ' ')
293 0 : call add_default ('DPDLFICE ', history_budget_histfile_num, ' ')
294 0 : call add_default ('SHDLFLIQ ', history_budget_histfile_num, ' ')
295 0 : call add_default ('SHDLFICE ', history_budget_histfile_num, ' ')
296 0 : call add_default ('DPDLFT ', history_budget_histfile_num, ' ')
297 0 : call add_default ('SHDLFT ', history_budget_histfile_num, ' ')
298 0 : call add_default ('ZMDLF ', history_budget_histfile_num, ' ')
299 :
300 0 : call add_default ('MACPDT ', history_budget_histfile_num, ' ')
301 0 : call add_default ('MACPDQ ', history_budget_histfile_num, ' ')
302 0 : call add_default ('MACPDLIQ ', history_budget_histfile_num, ' ')
303 0 : call add_default ('MACPDICE ', history_budget_histfile_num, ' ')
304 :
305 0 : call add_default ('CLDVAPADJ', history_budget_histfile_num, ' ')
306 0 : call add_default ('CLDLIQLIM', history_budget_histfile_num, ' ')
307 0 : call add_default ('CLDLIQDET', history_budget_histfile_num, ' ')
308 0 : call add_default ('CLDLIQADJ', history_budget_histfile_num, ' ')
309 0 : call add_default ('CLDICELIM', history_budget_histfile_num, ' ')
310 0 : call add_default ('CLDICEDET', history_budget_histfile_num, ' ')
311 0 : call add_default ('CLDICEADJ', history_budget_histfile_num, ' ')
312 :
313 0 : call add_default ('CMELIQ ', history_budget_histfile_num, ' ')
314 :
315 : end if
316 :
317 : ! Get constituent indices
318 0 : call cnst_get_ind('CLDLIQ', ixcldliq)
319 0 : call cnst_get_ind('CLDICE', ixcldice)
320 0 : call cnst_get_ind('NUMLIQ', ixnumliq)
321 0 : call cnst_get_ind('NUMICE', ixnumice)
322 :
323 : ! Get physics buffer indices
324 0 : CC_T_idx = pbuf_get_index('CC_T')
325 0 : CC_qv_idx = pbuf_get_index('CC_qv')
326 0 : CC_ql_idx = pbuf_get_index('CC_ql')
327 0 : CC_qi_idx = pbuf_get_index('CC_qi')
328 0 : CC_nl_idx = pbuf_get_index('CC_nl')
329 0 : CC_ni_idx = pbuf_get_index('CC_ni')
330 0 : CC_qlst_idx = pbuf_get_index('CC_qlst')
331 0 : cmfmc_sh_idx = pbuf_get_index('CMFMC_SH')
332 :
333 0 : sh_frac_idx = pbuf_get_index('SH_FRAC')
334 0 : dp_frac_idx = pbuf_get_index('DP_FRAC')
335 :
336 0 : if (rhminl_opt > 0 .or. rhmini_opt > 0) then
337 0 : cmfr_det_idx = pbuf_get_index('cmfr_det', istat)
338 0 : if (istat < 0) call endrun(subname//': macrop option requires cmfr_det in pbuf')
339 0 : if (rhminl_opt > 0) then
340 0 : qlr_det_idx = pbuf_get_index('qlr_det', istat)
341 0 : if (istat < 0) call endrun(subname//': macrop option requires qlr_det in pbuf')
342 : end if
343 0 : if (rhmini_opt > 0) then
344 0 : qir_det_idx = pbuf_get_index('qir_det', istat)
345 0 : if (istat < 0) call endrun(subname//': macrop option requires qir_det in pbuf')
346 : end if
347 : end if
348 :
349 0 : if (rhminl_opt == 2 .or. rhmini_opt == 2) then
350 0 : tke_idx = pbuf_get_index('tke')
351 0 : if (rhminl_opt == 2) then
352 0 : qtl_flx_idx = pbuf_get_index('qtl_flx', istat)
353 0 : if (istat < 0) call endrun(subname//': macrop option requires qtl_flx in pbuf')
354 : end if
355 0 : if (rhmini_opt == 2) then
356 0 : qti_flx_idx = pbuf_get_index('qti_flx', istat)
357 0 : if (istat < 0) call endrun(subname//': macrop option requires qti_flx in pbuf')
358 : end if
359 : end if
360 :
361 : ! Init pbuf fields. Note that the fields CLD, CONCLD, QCWAT, LCWAT,
362 : ! ICCWAT, and TCWAT are initialized in phys_inidat.
363 0 : if (is_first_step()) then
364 0 : call pbuf_set_field(pbuf2d, ast_idx, 0._r8)
365 0 : call pbuf_set_field(pbuf2d, aist_idx, 0._r8)
366 0 : call pbuf_set_field(pbuf2d, alst_idx, 0._r8)
367 0 : call pbuf_set_field(pbuf2d, qist_idx, 0._r8)
368 0 : call pbuf_set_field(pbuf2d, qlst_idx, 0._r8)
369 0 : call pbuf_set_field(pbuf2d, nlwat_idx, 0._r8)
370 0 : call pbuf_set_field(pbuf2d, niwat_idx, 0._r8)
371 : end if
372 :
373 : ! the following are physpkg, so they need to be init every time
374 0 : call pbuf_set_field(pbuf2d, fice_idx, 0._r8)
375 0 : call pbuf_set_field(pbuf2d, cmeliq_idx, 0._r8)
376 :
377 0 : end subroutine macrop_driver_init
378 :
379 : !============================================================================ !
380 : ! !
381 : !============================================================================ !
382 :
383 :
384 0 : subroutine macrop_driver_tend( &
385 : state, ptend, dtime, landfrac, &
386 : ocnfrac, snowh, &
387 : dlf, dlf2, cmfmc, ts, &
388 : sst, zdu, &
389 : pbuf, &
390 : det_s, det_ice)
391 :
392 : !-------------------------------------------------------- !
393 : ! !
394 : ! Purpose: !
395 : ! !
396 : ! Interface to detrain, cloud fraction and !
397 : ! cloud macrophysics subroutines !
398 : ! !
399 : ! Author: A. Gettelman, C. Craig, Oct 2010 !
400 : ! based on stratiform_tend by D.B. Coleman 4/2010 !
401 : ! !
402 : !-------------------------------------------------------- !
403 :
404 0 : use cloud_fraction_fice, only: cloud_fraction_fice_run
405 : use physics_types, only: physics_state, physics_ptend
406 : use physics_types, only: physics_ptend_init, physics_update
407 : use physics_types, only: physics_ptend_sum, physics_state_copy
408 : use physics_types, only: physics_state_dealloc
409 : use cam_history, only: outfld
410 : use constituents, only: cnst_get_ind, pcnst
411 : use cldwat2m_macro, only: mmacro_pcond
412 : use physconst, only: cpair, tmelt, gravit, cappa, rair, pref, lapse_rate
413 : use time_manager, only: get_nstep
414 :
415 : use ref_pres, only: top_lev => trop_cloud_top_lev
416 :
417 : ! CCPPized schemes
418 : use convective_cloud_cover, only: convective_cloud_cover_run
419 : use compute_cloud_fraction, only: compute_cloud_fraction_run
420 : use ref_pres, only: trop_cloud_top_lev ! bring in explicitly for CCPPized scheme
421 :
422 : !
423 : ! Input arguments
424 : !
425 :
426 : type(physics_state), intent(in) :: state ! State variables
427 : type(physics_ptend), intent(out) :: ptend ! macrophysics parameterization tendencies
428 : type(physics_buffer_desc), pointer :: pbuf(:) ! Physics buffer
429 :
430 : real(r8), intent(in) :: dtime ! Timestep
431 : real(r8), intent(in) :: landfrac(pcols) ! Land fraction (fraction)
432 : real(r8), intent(in) :: ocnfrac (pcols) ! Ocean fraction (fraction)
433 : real(r8), intent(in) :: snowh(pcols) ! Snow depth over land, water equivalent (m)
434 : real(r8), intent(in) :: dlf(pcols,pver) ! Detrained water from convection schemes
435 : real(r8), intent(in) :: dlf2(pcols,pver) ! Detrained water from shallow convection scheme
436 : real(r8), intent(in) :: cmfmc(pcols,pverp) ! Deep + Shallow Convective mass flux [ kg /s/m^2 ]
437 :
438 : real(r8), intent(in) :: ts(pcols) ! Surface temperature
439 : real(r8), intent(in) :: sst(pcols) ! Sea surface temperature
440 : real(r8), intent(in) :: zdu(pcols,pver) ! Detrainment rate from deep convection
441 :
442 :
443 : ! These two variables are needed for energy check
444 : real(r8), intent(out) :: det_s(pcols) ! Integral of detrained static energy from ice
445 : real(r8), intent(out) :: det_ice(pcols) ! Integral of detrained ice for energy check
446 :
447 : !
448 : ! Local variables
449 : !
450 :
451 0 : type(physics_state) :: state_loc ! Local copy of the state variable
452 0 : type(physics_ptend) :: ptend_loc ! Local parameterization tendencies
453 :
454 : integer i,k
455 : integer :: lchnk ! Chunk identifier
456 : integer :: ncol ! Number of atmospheric columns
457 :
458 : ! Physics buffer fields
459 :
460 : integer itim_old
461 0 : real(r8), pointer, dimension(:,:) :: qcwat ! Cloud water old q
462 0 : real(r8), pointer, dimension(:,:) :: tcwat ! Cloud water old temperature
463 0 : real(r8), pointer, dimension(:,:) :: lcwat ! Cloud liquid water old q
464 0 : real(r8), pointer, dimension(:,:) :: iccwat ! Cloud ice water old q
465 0 : real(r8), pointer, dimension(:,:) :: nlwat ! Cloud liquid droplet number condentration. old.
466 0 : real(r8), pointer, dimension(:,:) :: niwat ! Cloud ice droplet number condentration. old.
467 0 : real(r8), pointer, dimension(:,:) :: CC_T ! Grid-mean microphysical tendency
468 0 : real(r8), pointer, dimension(:,:) :: CC_qv ! Grid-mean microphysical tendency
469 0 : real(r8), pointer, dimension(:,:) :: CC_ql ! Grid-mean microphysical tendency
470 0 : real(r8), pointer, dimension(:,:) :: CC_qi ! Grid-mean microphysical tendency
471 0 : real(r8), pointer, dimension(:,:) :: CC_nl ! Grid-mean microphysical tendency
472 0 : real(r8), pointer, dimension(:,:) :: CC_ni ! Grid-mean microphysical tendency
473 0 : real(r8), pointer, dimension(:,:) :: CC_qlst ! In-liquid stratus microphysical tendency
474 0 : real(r8), pointer, dimension(:,:) :: cld ! Total cloud fraction
475 0 : real(r8), pointer, dimension(:,:) :: ast ! Relative humidity cloud fraction
476 0 : real(r8), pointer, dimension(:,:) :: aist ! Physical ice stratus fraction
477 0 : real(r8), pointer, dimension(:,:) :: alst ! Physical liquid stratus fraction
478 0 : real(r8), pointer, dimension(:,:) :: qist ! Physical in-cloud IWC
479 0 : real(r8), pointer, dimension(:,:) :: qlst ! Physical in-cloud LWC
480 0 : real(r8), pointer, dimension(:,:) :: concld ! Convective cloud fraction
481 :
482 0 : real(r8), pointer, dimension(:,:) :: shfrc ! Cloud fraction from shallow convection scheme
483 0 : real(r8), pointer, dimension(:,:) :: cmfmc_sh ! Shallow convective mass flux (pcols,pverp) [ kg/s/m^2 ]
484 :
485 0 : real(r8), pointer, dimension(:,:) :: cmeliq
486 :
487 0 : real(r8), pointer, dimension(:,:) :: tke
488 0 : real(r8), pointer, dimension(:,:) :: qtl_flx
489 0 : real(r8), pointer, dimension(:,:) :: qti_flx
490 0 : real(r8), pointer, dimension(:,:) :: cmfr_det
491 0 : real(r8), pointer, dimension(:,:) :: qlr_det
492 0 : real(r8), pointer, dimension(:,:) :: qir_det
493 :
494 : ! Convective cloud to the physics buffer for purposes of ql contrib. to radn.
495 :
496 0 : real(r8), pointer, dimension(:,:) :: fice_ql ! Cloud ice/water partitioning ratio.
497 :
498 : ! ZM microphysics
499 : real(r8), pointer :: dlfzm(:,:) ! ZM detrained convective cloud water mixing ratio.
500 : real(r8), pointer :: dnlfzm(:,:) ! ZM detrained convective cloud water num concen.
501 : real(r8), pointer :: dnifzm(:,:) ! ZM detrained convective cloud ice num concen.
502 :
503 : real(r8) :: latsub
504 :
505 : ! tendencies for ice saturation adjustment
506 : real(r8) :: stend(pcols,pver)
507 : real(r8) :: qvtend(pcols,pver)
508 : real(r8) :: qitend(pcols,pver)
509 : real(r8) :: initend(pcols,pver)
510 :
511 : ! Local variables for cldfrc
512 :
513 : real(r8) cldst(pcols,pver) ! Stratus cloud fraction
514 : real(r8) rhcloud(pcols,pver) ! Relative humidity cloud (last timestep)
515 : real(r8) rhu00(pcols,pver) ! RH threshold for cloud
516 : real(r8) icecldf(pcols,pver) ! Ice cloud fraction
517 : real(r8) liqcldf(pcols,pver) ! Liquid cloud fraction (combined into cloud)
518 : real(r8) relhum(pcols,pver) ! RH, output to determine drh/da
519 :
520 : ! Local variables for macrophysics
521 :
522 : real(r8) rdtime ! 1./dtime
523 : real(r8) qtend(pcols,pver) ! Moisture tendencies
524 : real(r8) ttend(pcols,pver) ! Temperature tendencies
525 : real(r8) ltend(pcols,pver) ! Cloud liquid water tendencies
526 : real(r8) fice(pcols,pver) ! Fractional ice content within cloud
527 : real(r8) fsnow(pcols,pver) ! Fractional snow production
528 : real(r8) homoo(pcols,pver)
529 : real(r8) qcreso(pcols,pver)
530 : real(r8) prcio(pcols,pver)
531 : real(r8) praio(pcols,pver)
532 : real(r8) qireso(pcols,pver)
533 : real(r8) ftem(pcols,pver)
534 : real(r8) pracso (pcols,pver)
535 : real(r8) dpdlfliq(pcols,pver)
536 : real(r8) dpdlfice(pcols,pver)
537 : real(r8) shdlfliq(pcols,pver)
538 : real(r8) shdlfice(pcols,pver)
539 : real(r8) dpdlft (pcols,pver)
540 : real(r8) shdlft (pcols,pver)
541 :
542 : real(r8) dum1
543 : real(r8) qc(pcols,pver)
544 : real(r8) qi(pcols,pver)
545 : real(r8) nc(pcols,pver)
546 : real(r8) ni(pcols,pver)
547 :
548 : logical lq(pcnst)
549 :
550 : ! Output from mmacro_pcond
551 :
552 : real(r8) tlat(pcols,pver)
553 : real(r8) qvlat(pcols,pver)
554 : real(r8) qcten(pcols,pver)
555 : real(r8) qiten(pcols,pver)
556 : real(r8) ncten(pcols,pver)
557 : real(r8) niten(pcols,pver)
558 :
559 : ! Output from mmacro_pcond
560 :
561 : real(r8) qvadj(pcols,pver) ! Macro-physics adjustment tendency from "positive_moisture" call (vapor)
562 : real(r8) qladj(pcols,pver) ! Macro-physics adjustment tendency from "positive_moisture" call (liquid)
563 : real(r8) qiadj(pcols,pver) ! Macro-physics adjustment tendency from "positive_moisture" call (ice)
564 : real(r8) qllim(pcols,pver) ! Macro-physics tendency from "instratus_condensate" call (liquid)
565 : real(r8) qilim(pcols,pver) ! Macro-physics tendency from "instratus_condensate" call (ice)
566 :
567 : ! For revised macophysics, mmacro_pcond
568 :
569 : real(r8) itend(pcols,pver)
570 : real(r8) lmitend(pcols,pver)
571 : real(r8) zeros(pcols,pver)
572 : real(r8) t_inout(pcols,pver)
573 : real(r8) qv_inout(pcols,pver)
574 : real(r8) ql_inout(pcols,pver)
575 : real(r8) qi_inout(pcols,pver)
576 : real(r8) concld_old(pcols,pver)
577 :
578 : ! Note that below 'clr_old' is defined using 'alst_old' not 'ast_old' for full consistency with the
579 : ! liquid condensation process which is using 'alst' not 'ast'.
580 : ! For microconsistency use 'concld_old', since 'alst_old' was computed using 'concld_old'.
581 : ! Since convective updraft fractional area is small, it does not matter whether 'concld' or 'concld_old' is used.
582 : ! Note also that 'clri_old' is defined using 'ast_old' since current microphysics is operating on 'ast_old'
583 : real(r8) clrw_old(pcols,pver) ! (1 - concld_old - alst_old)
584 : real(r8) clri_old(pcols,pver) ! (1 - concld_old - ast_old)
585 :
586 : real(r8) nl_inout(pcols,pver)
587 : real(r8) ni_inout(pcols,pver)
588 :
589 : real(r8) nltend(pcols,pver)
590 : real(r8) nitend(pcols,pver)
591 :
592 :
593 : ! For detraining cumulus condensate into the 'stratus' without evaporation
594 : ! This is for use in mmacro_pcond
595 :
596 : real(r8) dlf_T(pcols,pver)
597 : real(r8) dlf_qv(pcols,pver)
598 : real(r8) dlf_ql(pcols,pver)
599 : real(r8) dlf_qi(pcols,pver)
600 : real(r8) dlf_nl(pcols,pver)
601 : real(r8) dlf_ni(pcols,pver)
602 :
603 : ! Local variables for CFMIP calculations
604 : real(r8) :: mr_lsliq(pcols,pver) ! mixing_ratio_large_scale_cloud_liquid (kg/kg)
605 : real(r8) :: mr_lsice(pcols,pver) ! mixing_ratio_large_scale_cloud_ice (kg/kg)
606 : real(r8) :: mr_ccliq(pcols,pver) ! mixing_ratio_convective_cloud_liquid (kg/kg)
607 : real(r8) :: mr_ccice(pcols,pver) ! mixing_ratio_convective_cloud_ice (kg/kg)
608 :
609 : ! CloudSat equivalent ice mass mixing ratio (kg/kg)
610 : real(r8) :: cldsice(pcols,pver)
611 :
612 : ! shallowcu, deepcu pbuf fields
613 0 : real(r8), pointer, dimension(:,:) :: deepcu ! deep convection cloud fraction
614 0 : real(r8), pointer, dimension(:,:) :: shallowcu ! shallow convection cloud fraction
615 :
616 : ! For CCPPized schemes
617 : character(len=512) :: errmsg
618 : integer :: errflg
619 :
620 : ! ======================================================================
621 :
622 0 : lchnk = state%lchnk
623 0 : ncol = state%ncol
624 :
625 0 : call physics_state_copy(state, state_loc) ! Copy state to local state_loc.
626 :
627 : ! Associate pointers with physics buffer fields
628 :
629 0 : itim_old = pbuf_old_tim_idx()
630 :
631 0 : call pbuf_get_field(pbuf, qcwat_idx, qcwat, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
632 0 : call pbuf_get_field(pbuf, tcwat_idx, tcwat, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
633 0 : call pbuf_get_field(pbuf, lcwat_idx, lcwat, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
634 0 : call pbuf_get_field(pbuf, iccwat_idx, iccwat, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
635 0 : call pbuf_get_field(pbuf, nlwat_idx, nlwat, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
636 0 : call pbuf_get_field(pbuf, niwat_idx, niwat, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
637 :
638 0 : call pbuf_get_field(pbuf, cc_t_idx, cc_t, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
639 0 : call pbuf_get_field(pbuf, cc_qv_idx, cc_qv, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
640 0 : call pbuf_get_field(pbuf, cc_ql_idx, cc_ql, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
641 0 : call pbuf_get_field(pbuf, cc_qi_idx, cc_qi, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
642 0 : call pbuf_get_field(pbuf, cc_nl_idx, cc_nl, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
643 0 : call pbuf_get_field(pbuf, cc_ni_idx, cc_ni, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
644 0 : call pbuf_get_field(pbuf, cc_qlst_idx, cc_qlst, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
645 :
646 0 : call pbuf_get_field(pbuf, cld_idx, cld, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
647 0 : call pbuf_get_field(pbuf, concld_idx, concld, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
648 0 : call pbuf_get_field(pbuf, ast_idx, ast, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
649 0 : call pbuf_get_field(pbuf, aist_idx, aist, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
650 0 : call pbuf_get_field(pbuf, alst_idx, alst, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
651 0 : call pbuf_get_field(pbuf, qist_idx, qist, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
652 0 : call pbuf_get_field(pbuf, qlst_idx, qlst, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
653 :
654 0 : call pbuf_get_field(pbuf, cmeliq_idx, cmeliq)
655 :
656 : ! For purposes of convective ql.
657 :
658 0 : call pbuf_get_field(pbuf, fice_idx, fice_ql )
659 :
660 0 : call pbuf_get_field(pbuf, cmfmc_sh_idx, cmfmc_sh)
661 :
662 : ! check that qcwat and tcwat were initialized; if not then do it now.
663 0 : if (qcwat(1,1) == huge(1._r8)) then
664 0 : qcwat(:ncol,:) = state%q(:ncol,:,1)
665 : end if
666 0 : if (tcwat(1,1) == huge(1._r8)) then
667 0 : tcwat(:ncol,:) = state%t(:ncol,:)
668 : end if
669 :
670 : ! Initialize convective detrainment tendency
671 :
672 0 : dlf_T(:,:) = 0._r8
673 0 : dlf_qv(:,:) = 0._r8
674 0 : dlf_ql(:,:) = 0._r8
675 0 : dlf_qi(:,:) = 0._r8
676 0 : dlf_nl(:,:) = 0._r8
677 0 : dlf_ni(:,:) = 0._r8
678 :
679 : ! ------------------------------------- !
680 : ! From here, process computation begins !
681 : ! ------------------------------------- !
682 :
683 : ! ----------------------------------------------------------------------------- !
684 : ! Detrainment of convective condensate into the environment or stratiform cloud !
685 : ! ----------------------------------------------------------------------------- !
686 :
687 0 : lq(:) = .FALSE.
688 0 : lq(ixcldliq) = .TRUE.
689 0 : lq(ixcldice) = .TRUE.
690 0 : lq(ixnumliq) = .TRUE.
691 0 : lq(ixnumice) = .TRUE.
692 0 : call physics_ptend_init(ptend_loc, state%psetcols, 'pcwdetrain', ls=.true., lq=lq) ! Initialize local physics_ptend object
693 :
694 : ! Procedures :
695 : ! (1) Partition detrained convective cloud water into liquid and ice based on T.
696 : ! This also involves heating.
697 : ! If convection scheme can handle this internally, this step is not necssary.
698 : ! (2) Assuming a certain effective droplet radius, computes number concentration
699 : ! of detrained convective cloud liquid and ice.
700 : ! (3) If 'cu_det_st = .true' ('false'), detrain convective cloud 'liquid' into
701 : ! the pre-existing 'liquid' stratus ( mean environment ). The former does
702 : ! not involve any macrophysical evaporation while the latter does. This is
703 : ! a kind of 'targetted' deposition. Then, force in-stratus LWC to be bounded
704 : ! by qcst_min and qcst_max in mmacro_pcond.
705 : ! (4) In contrast to liquid, convective ice is detrained into the environment
706 : ! and involved in the sublimation. Similar bounds as liquid stratus are imposed.
707 : ! This is the key procesure generating upper-level cirrus clouds.
708 : ! The unit of dlf : [ kg/kg/s ]
709 :
710 0 : det_s(:) = 0._r8
711 0 : det_ice(:) = 0._r8
712 :
713 0 : dpdlfliq = 0._r8
714 0 : dpdlfice = 0._r8
715 0 : shdlfliq = 0._r8
716 0 : shdlfice = 0._r8
717 0 : dpdlft = 0._r8
718 0 : shdlft = 0._r8
719 :
720 0 : do k = top_lev, pver
721 0 : do i = 1, state_loc%ncol
722 0 : if( state_loc%t(i,k) > 268.15_r8 ) then
723 0 : dum1 = 0.0_r8
724 0 : elseif( state_loc%t(i,k) < 238.15_r8 ) then
725 0 : dum1 = 1.0_r8
726 : else
727 0 : dum1 = ( 268.15_r8 - state_loc%t(i,k) ) / 30._r8
728 : endif
729 :
730 : ! If detrainment was done elsewhere, still update the variables used for output
731 : ! assuming that the temperature split between liquid and ice is the same as assumed
732 : ! here.
733 0 : if (do_detrain) then
734 0 : ptend_loc%q(i,k,ixcldliq) = dlf(i,k) * ( 1._r8 - dum1 )
735 0 : ptend_loc%q(i,k,ixcldice) = dlf(i,k) * dum1
736 : ! dum2 = dlf(i,k) * ( 1._r8 - dum1 )
737 0 : ptend_loc%q(i,k,ixnumliq) = 3._r8 * ( max(0._r8, ( dlf(i,k) - dlf2(i,k) )) * ( 1._r8 - dum1 ) ) / &
738 : (4._r8*3.14_r8* 8.e-6_r8**3*997._r8) + & ! Deep Convection
739 0 : 3._r8 * ( dlf2(i,k) * ( 1._r8 - dum1 ) ) / &
740 0 : (4._r8*3.14_r8*10.e-6_r8**3*997._r8) ! Shallow Convection
741 : ! dum2 = dlf(i,k) * dum1
742 0 : ptend_loc%q(i,k,ixnumice) = 3._r8 * ( max(0._r8, ( dlf(i,k) - dlf2(i,k) )) * dum1 ) / &
743 : (4._r8*3.14_r8*25.e-6_r8**3*500._r8) + & ! Deep Convection
744 0 : 3._r8 * ( dlf2(i,k) * dum1 ) / &
745 0 : (4._r8*3.14_r8*50.e-6_r8**3*500._r8) ! Shallow Convection
746 0 : ptend_loc%s(i,k) = dlf(i,k) * dum1 * latice
747 : else
748 0 : ptend_loc%q(i,k,ixcldliq) = 0._r8
749 0 : ptend_loc%q(i,k,ixcldice) = 0._r8
750 0 : ptend_loc%q(i,k,ixnumliq) = 0._r8
751 0 : ptend_loc%q(i,k,ixnumice) = 0._r8
752 0 : ptend_loc%s(i,k) = 0._r8
753 : end if
754 :
755 : ! Only rliq is saved from deep convection, which is the reserved liquid. We need to keep
756 : ! track of the integrals of ice and static energy that is effected from conversion to ice
757 : ! so that the energy checker doesn't complain.
758 0 : det_s(i) = det_s(i) + ptend_loc%s(i,k)*state_loc%pdel(i,k)/gravit
759 0 : det_ice(i) = det_ice(i) - ptend_loc%q(i,k,ixcldice)*state_loc%pdel(i,k)/gravit
760 :
761 : ! Targetted detrainment of convective liquid water either directly into the
762 : ! existing liquid stratus or into the environment.
763 0 : if( cu_det_st ) then
764 : dlf_T(i,k) = ptend_loc%s(i,k)/cpair
765 : dlf_qv(i,k) = 0._r8
766 : dlf_ql(i,k) = ptend_loc%q(i,k,ixcldliq)
767 : dlf_qi(i,k) = ptend_loc%q(i,k,ixcldice)
768 : dlf_nl(i,k) = ptend_loc%q(i,k,ixnumliq)
769 : dlf_ni(i,k) = ptend_loc%q(i,k,ixnumice)
770 : ptend_loc%q(i,k,ixcldliq) = 0._r8
771 : ptend_loc%q(i,k,ixcldice) = 0._r8
772 : ptend_loc%q(i,k,ixnumliq) = 0._r8
773 : ptend_loc%q(i,k,ixnumice) = 0._r8
774 : ptend_loc%s(i,k) = 0._r8
775 : dpdlfliq(i,k) = 0._r8
776 : dpdlfice(i,k) = 0._r8
777 : shdlfliq(i,k) = 0._r8
778 : shdlfice(i,k) = 0._r8
779 : dpdlft (i,k) = 0._r8
780 : shdlft (i,k) = 0._r8
781 : else
782 0 : dpdlfliq(i,k) = ( dlf(i,k) - dlf2(i,k) ) * ( 1._r8 - dum1 )
783 0 : dpdlfice(i,k) = ( dlf(i,k) - dlf2(i,k) ) * ( dum1 )
784 0 : dpdlft (i,k) = ( dlf(i,k) - dlf2(i,k) ) * dum1 * latice/cpair
785 :
786 0 : shdlfliq(i,k) = dlf2(i,k) * ( 1._r8 - dum1 )
787 0 : shdlfice(i,k) = dlf2(i,k) * ( dum1 )
788 0 : shdlft (i,k) = dlf2(i,k) * dum1 * latice/cpair
789 : endif
790 : end do
791 : end do
792 :
793 0 : call outfld( 'DPDLFLIQ ', dpdlfliq, pcols, lchnk )
794 0 : call outfld( 'DPDLFICE ', dpdlfice, pcols, lchnk )
795 0 : call outfld( 'SHDLFLIQ ', shdlfliq, pcols, lchnk )
796 0 : call outfld( 'SHDLFICE ', shdlfice, pcols, lchnk )
797 0 : call outfld( 'DPDLFT ', dpdlft , pcols, lchnk )
798 0 : call outfld( 'SHDLFT ', shdlft , pcols, lchnk )
799 :
800 0 : call outfld( 'ZMDLF', dlf , pcols, state_loc%lchnk )
801 :
802 0 : det_ice(:ncol) = det_ice(:ncol)/1000._r8 ! divide by density of water
803 :
804 : ! Add the detrainment tendency to the output tendency
805 0 : call physics_ptend_init(ptend, state%psetcols, 'macrop')
806 0 : call physics_ptend_sum(ptend_loc, ptend, ncol)
807 :
808 : ! update local copy of state with the detrainment tendency
809 : ! ptend_loc is reset to zero by this call
810 0 : call physics_update(state_loc, ptend_loc, dtime)
811 :
812 : ! -------------------------------------- !
813 : ! Computation of Various Cloud Fractions !
814 : ! -------------------------------------- !
815 :
816 : ! ----------------------------------------------------------------------------- !
817 : ! Treatment of cloud fraction in CAM4 and CAM5 differs !
818 : ! (1) CAM4 !
819 : ! . Cumulus AMT = Deep Cumulus AMT ( empirical fcn of mass flux ) + !
820 : ! Shallow Cumulus AMT ( empirical fcn of mass flux ) !
821 : ! . Stratus AMT = max( RH stratus AMT, Stability Stratus AMT ) !
822 : ! . Cumulus and Stratus are 'minimally' overlapped without hierarchy. !
823 : ! . Cumulus LWC,IWC is assumed to be the same as Stratus LWC,IWC !
824 : ! (2) CAM5 !
825 : ! . Cumulus AMT = Deep Cumulus AMT ( empirical fcn of mass flux ) + !
826 : ! Shallow Cumulus AMT ( internally fcn of mass flux and w ) !
827 : ! . Stratus AMT = fcn of environmental-mean RH ( no Stability Stratus ) !
828 : ! . Cumulus and Stratus are non-overlapped with higher priority on Cumulus !
829 : ! . Cumulus ( both Deep and Shallow ) has its own LWC and IWC. !
830 : ! ----------------------------------------------------------------------------- !
831 :
832 0 : concld_old(:ncol,top_lev:pver) = concld(:ncol,top_lev:pver)
833 :
834 0 : nullify(tke, qtl_flx, qti_flx, cmfr_det, qlr_det, qir_det)
835 0 : if (tke_idx > 0) call pbuf_get_field(pbuf, tke_idx, tke)
836 0 : if (qtl_flx_idx > 0) call pbuf_get_field(pbuf, qtl_flx_idx, qtl_flx)
837 0 : if (qti_flx_idx > 0) call pbuf_get_field(pbuf, qti_flx_idx, qti_flx)
838 0 : if (cmfr_det_idx > 0) call pbuf_get_field(pbuf, cmfr_det_idx, cmfr_det)
839 0 : if (qlr_det_idx > 0) call pbuf_get_field(pbuf, qlr_det_idx, qlr_det)
840 0 : if (qir_det_idx > 0) call pbuf_get_field(pbuf, qir_det_idx, qir_det)
841 :
842 0 : clrw_old(:ncol,:top_lev-1) = 0._r8
843 0 : clri_old(:ncol,:top_lev-1) = 0._r8
844 0 : do k = top_lev, pver
845 0 : do i = 1, ncol
846 0 : clrw_old(i,k) = max( 0._r8, min( 1._r8, 1._r8 - concld(i,k) - alst(i,k) ) )
847 0 : clri_old(i,k) = max( 0._r8, min( 1._r8, 1._r8 - concld(i,k) - ast(i,k) ) )
848 : end do
849 : end do
850 :
851 0 : if( use_shfrc ) then
852 0 : call pbuf_get_field(pbuf, shfrc_idx, shfrc )
853 : else
854 0 : allocate(shfrc(pcols,pver))
855 0 : shfrc(:,:) = 0._r8
856 : endif
857 :
858 : ! CAM5 only uses 'concld' output from the below subroutine.
859 : ! Stratus ('ast' = max(alst,aist)) and total cloud fraction ('cld = ast + concld')
860 : ! will be computed using this updated 'concld' in the stratiform macrophysics
861 : ! scheme (mmacro_pcond) later below.
862 :
863 0 : call t_startf("cldfrc")
864 :
865 0 : call pbuf_get_field(pbuf, sh_frac_idx, shallowcu )
866 0 : call pbuf_get_field(pbuf, dp_frac_idx, deepcu )
867 :
868 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
869 0 : shallowcu(:,:) = 0._r8
870 0 : deepcu(:,:) = 0._r8
871 0 : concld(:,:) = 0._r8
872 : !REMOVECAM_END
873 :
874 : ! compute convective cloud fraction using CCPP-ized subroutine
875 : call convective_cloud_cover_run( &
876 : ncol = ncol, &
877 : pver = pver, &
878 : top_lev_cloudphys = trop_cloud_top_lev, & ! CAM5 macrophysics.
879 : use_shfrc = use_shfrc, &
880 0 : shfrc = shfrc(:ncol,:), &
881 0 : cmfmc_total = cmfmc(:ncol,:), &
882 0 : cmfmc_sh = cmfmc_sh(:ncol,:), &
883 0 : shallowcu = shallowcu(:ncol,:), &
884 0 : deepcu = deepcu(:ncol,:), &
885 0 : concld = concld(:ncol,:), &
886 : errmsg = errmsg, &
887 0 : errflg = errflg)
888 :
889 : ! write out convective cloud fraction diagnostic.
890 0 : call outfld( 'SH_CLD ', shallowcu , pcols, lchnk )
891 0 : call outfld( 'DP_CLD ', deepcu , pcols, lchnk )
892 :
893 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
894 0 : cld(:,:) = 0._r8
895 0 : rhcloud(:,:) = 0._r8
896 0 : cldst(:,:) = 0._r8
897 0 : rhu00(:,:) = 0._r8
898 0 : icecldf(:,:) = 0._r8
899 0 : liqcldf(:,:) = 0._r8
900 0 : relhum(:,:) = 0._r8
901 : !REMOVECAM_END
902 :
903 : ! call CCPPized cloud fraction scheme
904 : call compute_cloud_fraction_run( &
905 : ncol = ncol, &
906 : pver = pver, &
907 : cappa = cappa, &
908 : gravit = gravit, &
909 : rair = rair, &
910 : tmelt = tmelt, &
911 : pref = pref, &
912 : lapse_rate = lapse_rate, &
913 : top_lev_cloudphys = trop_cloud_top_lev, & ! CAM5 macrophysics.
914 0 : pmid = state_loc%pmid(:ncol,:), &
915 0 : ps = state_loc%pint(:ncol,pverp), &
916 0 : temp = state_loc%t(:ncol,:), &
917 0 : sst = sst(:ncol), &
918 0 : q = state_loc%q(:ncol,:,1), & ! note: assumes wv is at index 1.
919 0 : cldice = state_loc%q(:ncol,:,ixcldice), &
920 0 : phis = state_loc%phis(:ncol), &
921 0 : shallowcu = shallowcu(:ncol,:), &
922 0 : deepcu = deepcu(:ncol,:), &
923 0 : concld = concld(:ncol,:), &
924 0 : landfrac = landfrac(:ncol), &
925 0 : ocnfrac = ocnfrac(:ncol), &
926 0 : snowh = snowh(:ncol), &
927 : rhpert_flag = .false., & ! below output:
928 0 : cloud = cld(:ncol, :), &
929 0 : rhcloud = rhcloud(:ncol, :), &
930 0 : cldst = cldst(:ncol,:), &
931 0 : rhu00 = rhu00(:ncol,:), &
932 0 : icecldf = icecldf(:ncol,:), &
933 0 : liqcldf = liqcldf(:ncol,:), &
934 0 : relhum = relhum(:ncol,:), &
935 : errmsg = errmsg, &
936 0 : errflg = errflg)
937 :
938 0 : call t_stopf("cldfrc")
939 :
940 : ! ---------------------------------------------- !
941 : ! Stratiform Cloud Macrophysics and Microphysics !
942 : ! ---------------------------------------------- !
943 :
944 0 : lchnk = state_loc%lchnk
945 0 : ncol = state_loc%ncol
946 0 : rdtime = 1._r8/dtime
947 :
948 : ! Define fractional amount of stratus condensate and precipitation in ice phase.
949 : ! This uses a ramp ( -30 ~ -10 for fice, -5 ~ 0 for fsnow ).
950 : ! The ramp within convective cloud may be different
951 :
952 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
953 0 : fice(:,:) = 0._r8
954 0 : fsnow(:,:) = 0._r8
955 : !REMOVECAM_END
956 :
957 0 : call cloud_fraction_fice_run(ncol, state_loc%t(:ncol,:), tmelt, top_lev, pver, fice(:ncol,:), fsnow(:ncol,:), errmsg, errflg)
958 :
959 0 : lq(:) = .FALSE.
960 :
961 0 : lq(1) = .true.
962 0 : lq(ixcldice) = .true.
963 0 : lq(ixcldliq) = .true.
964 :
965 0 : lq(ixnumliq) = .true.
966 0 : lq(ixnumice) = .true.
967 :
968 : ! Initialize local physics_ptend object again
969 : call physics_ptend_init(ptend_loc, state%psetcols, 'macro_park', &
970 0 : ls=.true., lq=lq )
971 :
972 : ! --------------------------------- !
973 : ! Liquid Macrop_Driver Macrophysics !
974 : ! --------------------------------- !
975 :
976 0 : call t_startf('mmacro_pcond')
977 :
978 0 : zeros(:ncol,top_lev:pver) = 0._r8
979 0 : qc(:ncol,top_lev:pver) = state_loc%q(:ncol,top_lev:pver,ixcldliq)
980 0 : qi(:ncol,top_lev:pver) = state_loc%q(:ncol,top_lev:pver,ixcldice)
981 0 : nc(:ncol,top_lev:pver) = state_loc%q(:ncol,top_lev:pver,ixnumliq)
982 0 : ni(:ncol,top_lev:pver) = state_loc%q(:ncol,top_lev:pver,ixnumice)
983 :
984 : ! In CAM5, 'microphysical forcing' ( CC_... ) and 'the other advective forcings' ( ttend, ... )
985 : ! are separately provided into the prognostic microp_driver macrophysics scheme. This is an
986 : ! attempt to resolve in-cloud and out-cloud forcings.
987 :
988 0 : if( get_nstep() .le. 1 ) then
989 0 : tcwat(:ncol,top_lev:pver) = state_loc%t(:ncol,top_lev:pver)
990 0 : qcwat(:ncol,top_lev:pver) = state_loc%q(:ncol,top_lev:pver,1)
991 0 : lcwat(:ncol,top_lev:pver) = qc(:ncol,top_lev:pver) + qi(:ncol,top_lev:pver)
992 0 : iccwat(:ncol,top_lev:pver) = qi(:ncol,top_lev:pver)
993 0 : nlwat(:ncol,top_lev:pver) = nc(:ncol,top_lev:pver)
994 0 : niwat(:ncol,top_lev:pver) = ni(:ncol,top_lev:pver)
995 0 : ttend(:ncol,:) = 0._r8
996 0 : qtend(:ncol,:) = 0._r8
997 0 : ltend(:ncol,:) = 0._r8
998 0 : itend(:ncol,:) = 0._r8
999 0 : nltend(:ncol,:) = 0._r8
1000 0 : nitend(:ncol,:) = 0._r8
1001 0 : CC_T(:ncol,:) = 0._r8
1002 0 : CC_qv(:ncol,:) = 0._r8
1003 0 : CC_ql(:ncol,:) = 0._r8
1004 0 : CC_qi(:ncol,:) = 0._r8
1005 0 : CC_nl(:ncol,:) = 0._r8
1006 0 : CC_ni(:ncol,:) = 0._r8
1007 0 : CC_qlst(:ncol,:) = 0._r8
1008 : else
1009 0 : ttend(:ncol,top_lev:pver) = ( state_loc%t(:ncol,top_lev:pver) - tcwat(:ncol,top_lev:pver)) * rdtime &
1010 0 : - CC_T(:ncol,top_lev:pver)
1011 0 : qtend(:ncol,top_lev:pver) = ( state_loc%q(:ncol,top_lev:pver,1) - qcwat(:ncol,top_lev:pver)) * rdtime &
1012 0 : - CC_qv(:ncol,top_lev:pver)
1013 0 : ltend(:ncol,top_lev:pver) = ( qc(:ncol,top_lev:pver) + qi(:ncol,top_lev:pver) - lcwat(:ncol,top_lev:pver) ) * rdtime &
1014 0 : - (CC_ql(:ncol,top_lev:pver) + CC_qi(:ncol,top_lev:pver))
1015 0 : itend(:ncol,top_lev:pver) = ( qi(:ncol,top_lev:pver) - iccwat(:ncol,top_lev:pver)) * rdtime &
1016 0 : - CC_qi(:ncol,top_lev:pver)
1017 0 : nltend(:ncol,top_lev:pver) = ( nc(:ncol,top_lev:pver) - nlwat(:ncol,top_lev:pver)) * rdtime &
1018 0 : - CC_nl(:ncol,top_lev:pver)
1019 0 : nitend(:ncol,top_lev:pver) = ( ni(:ncol,top_lev:pver) - niwat(:ncol,top_lev:pver)) * rdtime &
1020 0 : - CC_ni(:ncol,top_lev:pver)
1021 : endif
1022 0 : lmitend(:ncol,top_lev:pver) = ltend(:ncol,top_lev:pver) - itend(:ncol,top_lev:pver)
1023 :
1024 0 : t_inout(:ncol,top_lev:pver) = tcwat(:ncol,top_lev:pver)
1025 0 : qv_inout(:ncol,top_lev:pver) = qcwat(:ncol,top_lev:pver)
1026 0 : ql_inout(:ncol,top_lev:pver) = lcwat(:ncol,top_lev:pver) - iccwat(:ncol,top_lev:pver)
1027 0 : qi_inout(:ncol,top_lev:pver) = iccwat(:ncol,top_lev:pver)
1028 0 : nl_inout(:ncol,top_lev:pver) = nlwat(:ncol,top_lev:pver)
1029 0 : ni_inout(:ncol,top_lev:pver) = niwat(:ncol,top_lev:pver)
1030 :
1031 : ! Liquid Microp_Driver Macrophysics.
1032 : ! The main roles of this subroutines are
1033 : ! (1) compute net condensation rate of stratiform liquid ( cmeliq )
1034 : ! (2) compute liquid stratus and ice stratus fractions.
1035 : ! Note 'ttend...' are advective tendencies except microphysical process while
1036 : ! 'CC...' are microphysical tendencies.
1037 :
1038 : call mmacro_pcond( lchnk, ncol, dtime, state_loc%pmid, state_loc%pdel, &
1039 : t_inout, qv_inout, ql_inout, qi_inout, nl_inout, ni_inout, &
1040 : ttend, qtend, lmitend, itend, nltend, nitend, &
1041 : CC_T, CC_qv, CC_ql, CC_qi, CC_nl, CC_ni, CC_qlst, &
1042 : dlf_T, dlf_qv, dlf_ql, dlf_qi, dlf_nl, dlf_ni, &
1043 : concld_old, concld, clrw_old, clri_old, landfrac, snowh, &
1044 : tke, qtl_flx, qti_flx, cmfr_det, qlr_det, qir_det, &
1045 : tlat, qvlat, qcten, qiten, ncten, niten, &
1046 : cmeliq, qvadj, qladj, qiadj, qllim, qilim, &
1047 0 : cld, alst, aist, qlst, qist, do_cldice )
1048 :
1049 : ! Copy of concld/fice to put in physics buffer
1050 : ! Below are used only for convective cloud.
1051 :
1052 0 : fice_ql(:ncol,:top_lev-1) = 0._r8
1053 0 : fice_ql(:ncol,top_lev:pver) = fice(:ncol,top_lev:pver)
1054 :
1055 :
1056 : ! Compute net stratus fraction using maximum over-lapping assumption
1057 0 : ast(:ncol,:top_lev-1) = 0._r8
1058 0 : ast(:ncol,top_lev:pver) = max( alst(:ncol,top_lev:pver), aist(:ncol,top_lev:pver) )
1059 :
1060 0 : call t_stopf('mmacro_pcond')
1061 :
1062 0 : do k = top_lev, pver
1063 0 : do i = 1, ncol
1064 0 : ptend_loc%s(i,k) = tlat(i,k)
1065 0 : ptend_loc%q(i,k,1) = qvlat(i,k)
1066 0 : ptend_loc%q(i,k,ixcldliq) = qcten(i,k)
1067 0 : ptend_loc%q(i,k,ixcldice) = qiten(i,k)
1068 0 : ptend_loc%q(i,k,ixnumliq) = ncten(i,k)
1069 0 : ptend_loc%q(i,k,ixnumice) = niten(i,k)
1070 :
1071 : ! Check to make sure that the macrophysics code is respecting the flags that control
1072 : ! whether cldwat should be prognosing cloud ice and cloud liquid or not.
1073 0 : if ((.not. do_cldice) .and. (qiten(i,k) /= 0.0_r8)) then
1074 : call endrun("macrop_driver:ERROR - "// &
1075 0 : "Cldwat is configured not to prognose cloud ice, but mmacro_pcond has ice mass tendencies.")
1076 : end if
1077 0 : if ((.not. do_cldice) .and. (niten(i,k) /= 0.0_r8)) then
1078 : call endrun("macrop_driver:ERROR -"// &
1079 0 : " Cldwat is configured not to prognose cloud ice, but mmacro_pcond has ice number tendencies.")
1080 : end if
1081 :
1082 0 : if ((.not. do_cldliq) .and. (qcten(i,k) /= 0.0_r8)) then
1083 : call endrun("macrop_driver:ERROR - "// &
1084 0 : "Cldwat is configured not to prognose cloud liquid, but mmacro_pcond has liquid mass tendencies.")
1085 : end if
1086 0 : if ((.not. do_cldliq) .and. (ncten(i,k) /= 0.0_r8)) then
1087 : call endrun("macrop_driver:ERROR - "// &
1088 0 : "Cldwat is configured not to prognose cloud liquid, but mmacro_pcond has liquid number tendencies.")
1089 : end if
1090 : end do
1091 : end do
1092 :
1093 : ! update the output tendencies with the mmacro_pcond tendencies
1094 0 : call physics_ptend_sum(ptend_loc, ptend, ncol)
1095 :
1096 : ! state_loc is the equlibrium state after macrophysics
1097 0 : call physics_update(state_loc, ptend_loc, dtime)
1098 :
1099 0 : call outfld('CLR_LIQ', clrw_old, pcols, lchnk)
1100 0 : call outfld('CLR_ICE', clri_old, pcols, lchnk)
1101 :
1102 0 : call outfld( 'MACPDT ', tlat , pcols, lchnk )
1103 0 : call outfld( 'MACPDQ ', qvlat, pcols, lchnk )
1104 0 : call outfld( 'MACPDLIQ ', qcten, pcols, lchnk )
1105 0 : call outfld( 'MACPDICE ', qiten, pcols, lchnk )
1106 0 : call outfld( 'CLDVAPADJ', qvadj, pcols, lchnk )
1107 0 : call outfld( 'CLDLIQADJ', qladj, pcols, lchnk )
1108 0 : call outfld( 'CLDICEADJ', qiadj, pcols, lchnk )
1109 0 : call outfld( 'CLDLIQDET', dlf_ql, pcols, lchnk )
1110 0 : call outfld( 'CLDICEDET', dlf_qi, pcols, lchnk )
1111 0 : call outfld( 'CLDLIQLIM', qllim, pcols, lchnk )
1112 0 : call outfld( 'CLDICELIM', qilim, pcols, lchnk )
1113 :
1114 0 : call outfld( 'ICECLDF ', aist, pcols, lchnk )
1115 0 : call outfld( 'LIQCLDF ', alst, pcols, lchnk )
1116 0 : call outfld( 'AST', ast, pcols, lchnk )
1117 :
1118 0 : call outfld( 'CONCLD ', concld, pcols, lchnk )
1119 0 : call outfld( 'CLDST ', cldst, pcols, lchnk )
1120 :
1121 0 : call outfld( 'CMELIQ' , cmeliq, pcols, lchnk )
1122 :
1123 :
1124 : ! calculations and outfld calls for CLDLIQSTR, CLDICESTR, CLDLIQCON, CLDICECON for CFMIP
1125 :
1126 : ! initialize local variables
1127 0 : mr_ccliq = 0._r8 !! not seen by radiation, so setting to 0
1128 0 : mr_ccice = 0._r8 !! not seen by radiation, so setting to 0
1129 0 : mr_lsliq = 0._r8
1130 0 : mr_lsice = 0._r8
1131 :
1132 0 : do k=top_lev,pver
1133 0 : do i=1,ncol
1134 0 : if (cld(i,k) .gt. 0._r8) then
1135 0 : mr_lsliq(i,k) = state_loc%q(i,k,ixcldliq)
1136 0 : mr_lsice(i,k) = state_loc%q(i,k,ixcldice)
1137 : else
1138 0 : mr_lsliq(i,k) = 0._r8
1139 0 : mr_lsice(i,k) = 0._r8
1140 : end if
1141 : end do
1142 : end do
1143 :
1144 0 : call outfld( 'CLDLIQSTR ', mr_lsliq, pcols, lchnk )
1145 0 : call outfld( 'CLDICESTR ', mr_lsice, pcols, lchnk )
1146 0 : call outfld( 'CLDLIQCON ', mr_ccliq, pcols, lchnk )
1147 0 : call outfld( 'CLDICECON ', mr_ccice, pcols, lchnk )
1148 :
1149 : ! ------------------------------------------------- !
1150 : ! Save equilibrium state variables for macrophysics !
1151 : ! at the next time step !
1152 : ! ------------------------------------------------- !
1153 0 : cldsice = 0._r8
1154 0 : do k = top_lev, pver
1155 0 : tcwat(:ncol,k) = state_loc%t(:ncol,k)
1156 0 : qcwat(:ncol,k) = state_loc%q(:ncol,k,1)
1157 0 : lcwat(:ncol,k) = state_loc%q(:ncol,k,ixcldliq) + state_loc%q(:ncol,k,ixcldice)
1158 0 : iccwat(:ncol,k) = state_loc%q(:ncol,k,ixcldice)
1159 0 : nlwat(:ncol,k) = state_loc%q(:ncol,k,ixnumliq)
1160 0 : niwat(:ncol,k) = state_loc%q(:ncol,k,ixnumice)
1161 0 : cldsice(:ncol,k) = lcwat(:ncol,k) * min(1.0_r8, max(0.0_r8, (tmelt - tcwat(:ncol,k)) / 20._r8))
1162 : end do
1163 :
1164 0 : call outfld( 'CLDSICE' , cldsice, pcols, lchnk )
1165 :
1166 : ! ptend_loc is deallocated in physics_update above
1167 0 : call physics_state_dealloc(state_loc)
1168 :
1169 0 : end subroutine macrop_driver_tend
1170 :
1171 :
1172 : ! Saturation adjustment for liquid
1173 : !
1174 : ! With CLUBB, we are seeing relative humidity with respect to water
1175 : ! greater than 1. This should not be happening and is not what the
1176 : ! microphsyics expects from the macrophysics. As a work around while
1177 : ! this issue is investigated in CLUBB, this routine will enfornce a
1178 : ! maximum RHliq of 1 everywhere in the atmosphere. Any excess water will
1179 : ! be converted into cloud drops.
1180 0 : subroutine liquid_macro_tend(npccn,t,p,qv,qc,nc,xxlv,deltat,stend,qvtend,qctend,nctend,vlen)
1181 :
1182 0 : use wv_sat_methods, only: wv_sat_qsat_ice_vect, wv_sat_qsat_water_vect
1183 : use micro_pumas_utils, only: rhow
1184 : use physconst, only: rair
1185 : use cldfrc2m, only: rhmini_const, rhmaxi_const
1186 :
1187 : integer, intent(in) :: vlen
1188 : real(r8), dimension(vlen), intent(in) :: npccn !Activated number of cloud condensation nuclei
1189 : real(r8), dimension(vlen), intent(in) :: t !temperature (k)
1190 : real(r8), dimension(vlen), intent(in) :: p !pressure (pa)
1191 : real(r8), dimension(vlen), intent(in) :: qv !water vapor mixing ratio
1192 : real(r8), dimension(vlen), intent(in) :: qc !liquid mixing ratio
1193 : real(r8), dimension(vlen), intent(in) :: nc !liquid number concentration
1194 : real(r8), intent(in) :: xxlv !latent heat of vaporization
1195 : real(r8), intent(in) :: deltat !timestep
1196 : real(r8), dimension(vlen), intent(out) :: stend ! 'temperature' tendency
1197 : real(r8), dimension(vlen), intent(out) :: qvtend !vapor tendency
1198 : real(r8), dimension(vlen), intent(out) :: qctend !liquid mass tendency
1199 : real(r8), dimension(vlen), intent(out) :: nctend !liquid number tendency
1200 :
1201 0 : real(r8) :: ESL(vlen)
1202 0 : real(r8) :: QSL(vlen)
1203 : real(r8) :: drop_size_param
1204 : integer :: i
1205 :
1206 0 : drop_size_param = 3._r8/(4._r8*3.14_r8*6.e-6_r8**3*rhow)
1207 :
1208 0 : do i = 1, vlen
1209 0 : stend(i) = 0._r8
1210 0 : qvtend(i) = 0._r8
1211 0 : qctend(i) = 0._r8
1212 0 : nctend(i) = 0._r8
1213 : end do
1214 :
1215 : ! calculate qsatl from t,p,q
1216 : !$acc data copyin(t,p) copyout(ESL,QSL)
1217 0 : call wv_sat_qsat_water_vect(t, p, ESL, QSL, vlen)
1218 : !$acc end data
1219 :
1220 0 : do i = 1, vlen
1221 : ! Don't allow supersaturation with respect to liquid.
1222 0 : if (qv(i) > QSL(i)) then
1223 :
1224 0 : qctend(i) = (qv(i) - QSL(i)) / deltat
1225 0 : qvtend(i) = 0._r8 - qctend(i)
1226 0 : stend(i) = qctend(i) * xxlv ! moist static energy tend...[J/kg/s] !
1227 :
1228 : ! If drops exists (more than 1 L-1) and there is condensation,
1229 : ! do not add to number (= growth), otherwise add 6um drops.
1230 : !
1231 : ! This is somewhat arbitrary, but ensures that some reasonable droplet
1232 : ! size is created to remove the excess water. This could be enhanced to
1233 : ! look at npccn, but ideally this entire routine should go away.
1234 0 : if ((nc(i)*p(i)/rair/t(i) < 1e3_r8) .and. (qc(i)+qctend(i)*deltat > 1e-18_r8)) then
1235 0 : nctend(i) = nctend(i) + qctend(i)*drop_size_param
1236 : end if
1237 : end if
1238 : end do
1239 :
1240 0 : end subroutine liquid_macro_tend
1241 :
1242 : end module macrop_driver
|