Line data Source code
1 : !===========================================================================
2 : ! CAM interface to the UNIFIED CONVECTION SCHEME (UNICON)
3 : !
4 : ! The USE_UNICON macro converts this module to a stub interface which allows
5 : ! CAM to be built without the unicon and unicon_utils modules.
6 : !
7 : !===========================================================================
8 :
9 : module unicon_cam
10 :
11 : use shr_kind_mod, only: r8 => shr_kind_r8, i4 => shr_kind_i4
12 : use spmd_utils, only: masterproc
13 : use ppgrid, only: pcols, pver, pverp, begchunk, endchunk
14 : use physconst, only: rair, cpair, gravit, latvap, latice, zvir, mwdry
15 :
16 : use constituents, only: pcnst, cnst_add, qmin, cnst_get_type_byind, cnst_get_ind, cnst_name
17 : use rad_constituents, only: rad_cnst_get_info, rad_cnst_get_mode_num_idx, rad_cnst_get_mam_mmr_idx
18 : use physics_types, only: physics_state, physics_ptend, physics_ptend_init
19 : use camsrfexch, only: cam_in_t
20 : use physics_buffer, only: pbuf_add_field, dtype_r8, dyn_time_lvls, pbuf_old_tim_idx, &
21 : physics_buffer_desc, pbuf_get_index, pbuf_get_field, pbuf_set_field
22 :
23 : use phys_control, only: phys_getopts
24 : use cam_history, only: outfld, addfld, horiz_only
25 :
26 : use time_manager, only: is_first_step
27 : use cam_abortutils, only: endrun
28 :
29 : #ifdef USE_UNICON
30 : use unicon, only: unicon_init, compute_unicon
31 : use unicon_utils, only: unicon_utils_init, positive_moisture, positive_tracer
32 : #endif
33 :
34 : implicit none
35 : private
36 : save
37 :
38 : public :: &
39 : unicon_cam_readnl, &
40 : unicon_cam_register, &
41 : unicon_cam_init, &
42 : unicon_implements_cnst, &
43 : unicon_init_cnst, &
44 : unicon_out_t, &
45 : unicon_cam_tend, &
46 : unicon_cam_org_diags
47 :
48 : ! namelist variables
49 : logical :: unicon_offline_dat_out = .false.
50 : integer :: unicon_offline_dat_hfile = 2
51 :
52 : ! properties
53 : real(r8) :: xlv ! Latent heat of vaporization
54 : real(r8) :: xlf ! Latent heat of fusion
55 : real(r8) :: xls ! Latent heat of sublimation
56 : real(r8) :: cp ! Specific heat of dry air
57 :
58 : integer, parameter :: &
59 : nseg = 1, &! Number of updraft segments [ # ]
60 : mix = pcols, &! Maximum number of columns
61 : mkx = pver, &! Number of vertical layers
62 : ncnst = pcnst ! Number of advected constituents
63 :
64 : ! For advecting organization-related variables
65 : integer, parameter :: n_org = 5 ! Number of constituents
66 : character(len=8), dimension(n_org), parameter :: & ! Constituent names
67 : cnst_names = (/'ORGawk ','ORGthl ','ORGqto ','ORGuoo ','ORGvoo '/)
68 :
69 : integer :: awk_cnst_ind, thl_cnst_ind, qt_cnst_ind, u_cnst_ind, v_cnst_ind
70 :
71 : ! fields added to physics buffer by this module
72 : integer :: &
73 : cushavg_idx, &
74 : cuorg_idx, &
75 : awk_PBL_idx, &
76 : delta_thl_PBL_idx, &
77 : delta_qt_PBL_idx, &
78 : delta_u_PBL_idx, &
79 : delta_v_PBL_idx, &
80 : delta_tr_PBL_idx, &
81 : cu_cmfr_idx, &
82 : cu_thlr_idx, &
83 : cu_qtr_idx, &
84 : cu_ur_idx, &
85 : cu_vr_idx, &
86 : cu_qlr_idx, &
87 : cu_qir_idx, &
88 : cu_trr_idx, &
89 : cmfr_det_idx, &
90 : qlr_det_idx, &
91 : qir_det_idx
92 :
93 : ! fields expected to be in the physics buffer
94 : integer :: &
95 : sgh30_idx = -1, &
96 : ast_idx = -1, &
97 : tke_idx = -1, &
98 : bprod_idx = -1, &
99 : kpblh_idx = -1, &
100 : pblh_idx = -1, &
101 : went_idx = -1, &
102 : cush_idx = -1, &
103 : shfrc_idx = -1, &
104 : icwmrsh_idx = -1, &
105 : rprdsh_idx = -1, &
106 : prec_sh_idx = -1, &
107 : snow_sh_idx = -1, &
108 : nevapr_shcu_idx = -1, &
109 : am_evp_st_idx = -1, & ! Evaporation area of stratiform precipitation [fraction]
110 : evprain_st_idx = -1, & ! Grid-mean evaporation rate of stratiform rain [kg/kg/s] >= 0.
111 : evpsnow_st_idx = -1 ! Grid-mean evaporation rate of stratiform snow [kg/kg/s] >= 0.
112 :
113 : ! constituent indices
114 : integer :: ixcldliq, ixcldice, ixnumliq, ixnumice
115 :
116 : ! unicon output fields
117 : type unicon_out_t
118 : real(r8) :: cmfmc(mix,mkx+1) ! Upward convective mass flux at the interface [ kg / s / m2 ]
119 : real(r8) :: slflx(mix,mkx+1) ! Net upward convective flux of liquid static energy [ J / s / m2 ]
120 : real(r8) :: qtflx(mix,mkx+1) ! Net upward convective flux of total specific humidity [ kg / s / m2 ]
121 : real(r8) :: rqc(mix,mkx) ! Prod of suspended LWC+IWC by expel of excessive in-cumulus condensate [ kg / kg / s ] > 0
122 : real(r8) :: rliq(mix) ! Vertical integral of 'rqc' in flux unit [ m / s ]
123 : real(r8) :: cnt(mix) ! Cloud top interface index ( ki = kpen )
124 : real(r8) :: cnb(mix) ! Cloud base interface index ( ki = krel-1 )
125 : end type unicon_out_t
126 :
127 : ! logical array to identify constituents that are mode number concentrations
128 : logical :: cnst_is_mam_num(ncnst)
129 : ! logical array to identify constituents that are mode specie mass mixing ratios
130 : logical :: cnst_is_mam_mmr(ncnst)
131 :
132 : !==================================================================================================
133 : contains
134 : !==================================================================================================
135 :
136 : !> \brief Read namelist group unicon_nl
137 : !!
138 : !! \param[in] nlfile ! filepath for file containing namelist input
139 :
140 1536 : subroutine unicon_cam_readnl(nlfile)
141 :
142 : use namelist_utils, only: find_group_name
143 : use units, only: getunit, freeunit
144 : use mpishorthand
145 :
146 : character(len=*), intent(in) :: nlfile
147 :
148 : ! Local variables
149 : integer :: unitn, ierr
150 : character(len=*), parameter :: subname = 'unicon_cam_readnl'
151 :
152 : namelist /unicon_nl/ unicon_offline_dat_out, unicon_offline_dat_hfile
153 :
154 : !-----------------------------------------------------------------------------
155 :
156 1536 : if (masterproc) then
157 2 : unitn = getunit()
158 2 : open( unitn, file=trim(nlfile), status='old' )
159 2 : call find_group_name(unitn, 'unicon_nl', status=ierr)
160 2 : if (ierr == 0) then
161 0 : read(unitn, unicon_nl, iostat=ierr)
162 0 : if (ierr /= 0) then
163 0 : call endrun(subname // ':: ERROR reading namelist')
164 : end if
165 : end if
166 2 : close(unitn)
167 2 : call freeunit(unitn)
168 : end if
169 :
170 : #ifdef SPMD
171 : ! Broadcast namelist variables
172 1536 : call mpibcast (unicon_offline_dat_out, 1, mpilog, 0, mpicom)
173 1536 : call mpibcast (unicon_offline_dat_hfile, 1, mpiint, 0, mpicom)
174 : #endif
175 :
176 1536 : end subroutine unicon_cam_readnl
177 :
178 : !================================================================================================
179 :
180 0 : subroutine unicon_cam_register
181 :
182 : ! Register fields in the constituent array and the physics buffer.
183 :
184 :
185 : ! Jun.02.2012. Sungsu for advecting organization-related horizontal heterogeneity
186 : ! within PBL.
187 : ! For the time being, in order to save computation time, advection of aerosol perturbations
188 : ! are simply neglected.
189 :
190 : #ifdef USE_UNICON
191 :
192 : call cnst_add(cnst_names(1), mwdry, cpair, 0._r8, awk_cnst_ind, &
193 : 'Wake area within PBL associated with organization', readiv=.false., mixtype = 'dry')
194 : call cnst_add(cnst_names(2), mwdry, cpair, 0._r8, thl_cnst_ind, &
195 : 'Perturbation of thl associated with organization', readiv=.false., mixtype = 'dry')
196 : call cnst_add(cnst_names(3), mwdry, cpair, 0._r8, qt_cnst_ind, &
197 : 'Perturbation of qt associated with organization', readiv=.false., mixtype = 'dry')
198 : call cnst_add(cnst_names(4), mwdry, cpair, 0._r8, u_cnst_ind, &
199 : 'Perturbation of u associated with organization', readiv=.false., mixtype = 'dry')
200 : call cnst_add(cnst_names(5), mwdry, cpair, 0._r8, v_cnst_ind, &
201 : 'Perturbation of v associated with organization', readiv=.false., mixtype = 'dry')
202 :
203 :
204 : call pbuf_add_field('cushavg', 'global', dtype_r8, (/pcols,dyn_time_lvls/), cushavg_idx)
205 : call pbuf_add_field('cuorg', 'global', dtype_r8, (/pcols,dyn_time_lvls/), cuorg_idx)
206 : call pbuf_add_field('awk_PBL', 'global', dtype_r8, (/pcols,dyn_time_lvls/), awk_PBL_idx)
207 : call pbuf_add_field('delta_thl_PBL', 'global', dtype_r8, (/pcols,dyn_time_lvls/), delta_thl_PBL_idx)
208 : call pbuf_add_field('delta_qt_PBL', 'global', dtype_r8, (/pcols,dyn_time_lvls/), delta_qt_PBL_idx)
209 : call pbuf_add_field('delta_u_PBL', 'global', dtype_r8, (/pcols,dyn_time_lvls/), delta_u_PBL_idx)
210 : call pbuf_add_field('delta_v_PBL', 'global', dtype_r8, (/pcols,dyn_time_lvls/), delta_v_PBL_idx)
211 : call pbuf_add_field('delta_tr_PBL', 'global', dtype_r8, (/pcols,pcnst,dyn_time_lvls/), delta_tr_PBL_idx)
212 : call pbuf_add_field('cu_cmfr', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), cu_cmfr_idx)
213 : call pbuf_add_field('cu_thlr', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), cu_thlr_idx)
214 : call pbuf_add_field('cu_qtr', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), cu_qtr_idx)
215 : call pbuf_add_field('cu_ur', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), cu_ur_idx)
216 : call pbuf_add_field('cu_vr', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), cu_vr_idx)
217 : call pbuf_add_field('cu_qlr', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), cu_qlr_idx)
218 : call pbuf_add_field('cu_qir', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), cu_qir_idx)
219 : call pbuf_add_field('cu_trr', 'global', dtype_r8, (/pcols,pver,pcnst,dyn_time_lvls/), cu_trr_idx)
220 : call pbuf_add_field('cmfr_det', 'global', dtype_r8, (/pcols,pver/), cmfr_det_idx)
221 : call pbuf_add_field('qlr_det', 'global', dtype_r8, (/pcols,pver/), qlr_det_idx)
222 : call pbuf_add_field('qir_det', 'global', dtype_r8, (/pcols,pver/), qir_det_idx)
223 :
224 : #endif
225 :
226 0 : end subroutine unicon_cam_register
227 :
228 : !==================================================================================================
229 :
230 0 : subroutine unicon_cam_init(pbuf2d)
231 :
232 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
233 :
234 : ! local variables
235 : character(len=*), parameter :: sub='unicon_cam_init: '
236 : integer :: i, icnst, j, l, m, nmodes, nspec
237 :
238 : character(len=8) :: units
239 : character(len=30) :: varname
240 : character(len=60) :: surname
241 : character(len=2) :: numcha
242 : integer :: msfc
243 : ! ------------------------------------------------------------------------------------------- !
244 :
245 : #ifdef USE_UNICON
246 :
247 : ! constants
248 : xlv = latvap
249 : xlf = latice
250 : xls = xlv + xlf
251 : cp = cpair
252 :
253 : call unicon_init(latvap, cpair, latice, zvir, rair, gravit)
254 :
255 : ! save some constituent indices
256 : call cnst_get_ind('CLDLIQ', ixcldliq)
257 : call cnst_get_ind('CLDICE', ixcldice)
258 : call cnst_get_ind('NUMLIQ', ixnumliq)
259 : call cnst_get_ind('NUMICE', ixnumice)
260 :
261 : ! save some physics buffer indices
262 : sgh30_idx = pbuf_get_index('SGH30')
263 : ast_idx = pbuf_get_index('AST')
264 : tke_idx = pbuf_get_index('tke')
265 : bprod_idx = pbuf_get_index('bprod')
266 : kpblh_idx = pbuf_get_index('kpblh')
267 : pblh_idx = pbuf_get_index('pblh')
268 : went_idx = pbuf_get_index('went')
269 : cush_idx = pbuf_get_index('cush')
270 : shfrc_idx = pbuf_get_index('shfrc')
271 : icwmrsh_idx = pbuf_get_index('ICWMRSH')
272 : rprdsh_idx = pbuf_get_index('RPRDSH')
273 : prec_sh_idx = pbuf_get_index('PREC_SH')
274 : snow_sh_idx = pbuf_get_index('SNOW_SH')
275 : nevapr_shcu_idx = pbuf_get_index('NEVAPR_SHCU')
276 : am_evp_st_idx = pbuf_get_index('am_evp_st')
277 : evprain_st_idx = pbuf_get_index('evprain_st')
278 : evpsnow_st_idx = pbuf_get_index('evpsnow_st')
279 :
280 : ! physics buffer fields that need initializers -- these are only
281 : ! fields that have been added to pbuf by this module, or by the
282 : ! convection driver layer.
283 : if (is_first_step()) then
284 :
285 : if (cush_idx > 0) then
286 : call pbuf_set_field(pbuf2d, cush_idx, 1.e3_r8)
287 : else
288 : call endrun(sub//'cush not in pbuf')
289 : end if
290 : if (cushavg_idx > 0) then
291 : call pbuf_set_field(pbuf2d, cushavg_idx, 1.e3_r8)
292 : else
293 : call endrun(sub//'cushavg not in pbuf')
294 : end if
295 : if (cuorg_idx > 0) then
296 : call pbuf_set_field(pbuf2d, cuorg_idx, 0.0_r8)
297 : else
298 : call endrun(sub//'cuorg not in pbuf')
299 : end if
300 : if (awk_pbl_idx > 0) then
301 : call pbuf_set_field(pbuf2d, awk_pbl_idx, 0.0_r8)
302 : else
303 : call endrun(sub//'awk_PBL not in pbuf')
304 : end if
305 : if (delta_thl_PBL_idx > 0) then
306 : call pbuf_set_field(pbuf2d, delta_thl_PBL_idx, 0.0_r8)
307 : else
308 : call endrun(sub//'delta_thl_PBL not in pbuf')
309 : end if
310 : if (delta_qt_PBL_idx > 0) then
311 : call pbuf_set_field(pbuf2d, delta_qt_PBL_idx, 0.0_r8)
312 : else
313 : call endrun(sub//'delta_qt_PBL not in pbuf')
314 : end if
315 : if (delta_u_PBL_idx > 0) then
316 : call pbuf_set_field(pbuf2d, delta_u_PBL_idx, 0.0_r8)
317 : else
318 : call endrun(sub//'delta_u_PBL not in pbuf')
319 : end if
320 : if (delta_v_PBL_idx > 0) then
321 : call pbuf_set_field(pbuf2d, delta_v_PBL_idx, 0.0_r8)
322 : else
323 : call endrun(sub//'delta_v_PBL not in pbuf')
324 : end if
325 : if (delta_tr_PBL_idx > 0) then
326 : call pbuf_set_field(pbuf2d, delta_tr_PBL_idx, 0.0_r8)
327 : else
328 : call endrun(sub//'delta_tr_PBL not in pbuf')
329 : end if
330 : if (cu_cmfr_idx > 0) then
331 : call pbuf_set_field(pbuf2d, cu_cmfr_idx, 0.0_r8)
332 : else
333 : call endrun(sub//'cu_cmfr not in pbuf')
334 : end if
335 : if (cu_thlr_idx > 0) then
336 : call pbuf_set_field(pbuf2d, cu_thlr_idx, 0.0_r8)
337 : else
338 : call endrun(sub//'cu_thlr not in pbuf')
339 : end if
340 : if (cu_qtr_idx > 0) then
341 : call pbuf_set_field(pbuf2d, cu_qtr_idx, 0.0_r8)
342 : else
343 : call endrun(sub//'cu_qtr not in pbuf')
344 : end if
345 : if (cu_ur_idx > 0) then
346 : call pbuf_set_field(pbuf2d, cu_ur_idx, 0.0_r8)
347 : else
348 : call endrun(sub//'cu_ur not in pbuf')
349 : end if
350 : if (cu_vr_idx > 0) then
351 : call pbuf_set_field(pbuf2d, cu_vr_idx, 0.0_r8)
352 : else
353 : call endrun(sub//'cu_vr not in pbuf')
354 : end if
355 : if (cu_qlr_idx > 0) then
356 : call pbuf_set_field(pbuf2d, cu_qlr_idx, 0.0_r8)
357 : else
358 : call endrun(sub//'cu_qlr not in pbuf')
359 : end if
360 : if (cu_qir_idx > 0) then
361 : call pbuf_set_field(pbuf2d, cu_qir_idx, 0.0_r8)
362 : else
363 : call endrun(sub//'cu_qir not in pbuf')
364 : end if
365 : if (cu_trr_idx > 0) then
366 : call pbuf_set_field(pbuf2d, cu_trr_idx, 0.0_r8)
367 : else
368 : call endrun(sub//'cu_trr not in pbuf')
369 : end if
370 : if (cmfr_det_idx > 0) then
371 : call pbuf_set_field(pbuf2d, cmfr_det_idx, 0.0_r8)
372 : else
373 : call endrun(sub//'cmfr_det not in pbuf')
374 : end if
375 : if (qlr_det_idx > 0) then
376 : call pbuf_set_field(pbuf2d, qlr_det_idx, 0.0_r8)
377 : else
378 : call endrun(sub//'qlr_det not in pbuf')
379 : end if
380 : if (qir_det_idx > 0) then
381 : call pbuf_set_field(pbuf2d, qir_det_idx, 0.0_r8)
382 : else
383 : call endrun(sub//'qir_det not in pbuf')
384 : end if
385 :
386 : end if
387 :
388 : ! Set arrays to identify the modal aerosol constituents
389 : cnst_is_mam_num = .false.
390 : cnst_is_mam_mmr = .false.
391 : call rad_cnst_get_info(0, nmodes=nmodes)
392 : do i = 1, nmodes
393 : call rad_cnst_get_mode_num_idx(i, icnst)
394 : cnst_is_mam_num(icnst) = .true.
395 : call rad_cnst_get_info(0, i, nspec=nspec)
396 : do j = 1, nspec
397 : call rad_cnst_get_mam_mmr_idx(i, j, icnst)
398 : cnst_is_mam_mmr(icnst) = .true.
399 : end do
400 : end do
401 :
402 : ! ------------------------- !
403 : ! Internal Output Variables !
404 : ! ------------------------- !
405 :
406 : ! Sungsu for advection of convective organization
407 :
408 : call addfld('ORGawk', (/ 'lev' /), 'A', 'no', 'Convective Organization - Wake Area' )
409 : call addfld('ORGthl', (/ 'lev' /), 'A', 'K', 'Convective Organization - Perturbation of thl in the non-wake area' )
410 : call addfld('ORGqto', (/ 'lev' /), 'A', 'kg/kg', 'Convective Organization - Perturbation of qt in the non-wake area' )
411 : call addfld('ORGuoo', (/ 'lev' /), 'A', 'm/s', 'Convective Organization - Perturbation of u in the non-wake area' )
412 : call addfld('ORGvoo', (/ 'lev' /), 'A', 'm/s', 'Convective Organization - Perturbation of v in the non-wake area' )
413 :
414 : ! From the main unified convection scheme
415 :
416 : call addfld('flxrain_SP', (/ 'ilev' /), 'A', 'kg/m/s2', 'Convective net rain flux')
417 : call addfld('flxsnow_SP', (/ 'ilev' /), 'A', 'kg/m/s2', 'Convective net snow flux')
418 :
419 : call addfld('cmf_SP', (/ 'ilev' /), 'A', 'kg/m2/s', 'Convective net mass flux')
420 : call addfld('qtflx_SP', (/ 'ilev' /), 'A', 'kg/m2/s', 'Convective net qt flux')
421 : call addfld('slflx_SP', (/ 'ilev' /), 'A', 'J/m2/s' , 'Convective net sl flux')
422 : call addfld('uflx_SP', (/ 'ilev' /), 'A', 'kg/m/s2', 'Convective net u flux')
423 : call addfld('vflx_SP', (/ 'ilev' /), 'A', 'kg/m/s2', 'Convective net v flux')
424 :
425 : call addfld('cmf_u_SP', (/ 'ilev' /), 'A', 'kg/m2/s', 'Convective updraft mass flux')
426 : call addfld('qtflx_u_SP', (/ 'ilev' /), 'A', 'kg/m2/s', 'Convective updraft qt flux')
427 : call addfld('slflx_u_SP', (/ 'ilev' /), 'A', 'J/m2/s' , 'Convective updraft sl flux')
428 : call addfld('uflx_u_SP', (/ 'ilev' /), 'A', 'kg/m/s2', 'Convective updraft u flux')
429 : call addfld('vflx_u_SP', (/ 'ilev' /), 'A', 'kg/m/s2', 'Convective updraft v flux')
430 :
431 : call addfld('cmf_d_SP', (/ 'ilev' /), 'A', 'kg/m2/s', 'Convective downdraft mass flux')
432 : call addfld('qtflx_d_SP', (/ 'ilev' /), 'A', 'kg/m2/s', 'Convective downdraft qt flux')
433 : call addfld('slflx_d_SP', (/ 'ilev' /), 'A', 'J/m2/s' , 'Convective downdraft sl flux')
434 : call addfld('uflx_d_SP', (/ 'ilev' /), 'A', 'kg/m/s2', 'Convective downdraft u flux')
435 : call addfld('vflx_d_SP', (/ 'ilev' /), 'A', 'kg/m/s2', 'Convective downdraft v flux')
436 :
437 : call addfld('qtten_u_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qt by updraft')
438 : call addfld('slten_u_SP', (/ 'lev' /), 'A', 'J/kg/s', 'Convective tendency sl by updraft')
439 : call addfld('uten_u_SP', (/ 'lev' /), 'A', 'm/s/s', 'Convective tendency u by updraft')
440 : call addfld('vten_u_SP', (/ 'lev' /), 'A', 'm/s/s' , 'Convective tendency v by updraft')
441 : call addfld('sten_u_SP', (/ 'lev' /), 'A', 'J/kg/s', 'Convective tendency s by updraft')
442 : call addfld('qvten_u_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qv by updraft')
443 : call addfld('qlten_u_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency ql by updraft')
444 : call addfld('qiten_u_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qi by updraft')
445 : call addfld('nlten_u_SP', (/ 'lev' /), 'A', '1/kg/s', 'Convective tendency nl by updraft')
446 : call addfld('niten_u_SP', (/ 'lev' /), 'A', '1/kg/s', 'Convective tendency ni by updraft')
447 :
448 : call addfld('qtten_d_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qt by downdraft')
449 : call addfld('slten_d_SP', (/ 'lev' /), 'A', 'J/kg/s', 'Convective tendency sl by downdraft')
450 : call addfld('uten_d_SP', (/ 'lev' /), 'A', 'm/s/s', 'Convective tendency u by downdraft')
451 : call addfld('vten_d_SP', (/ 'lev' /), 'A', 'm/s/s', 'Convective tendency v by downdraft')
452 : call addfld('sten_d_SP', (/ 'lev' /), 'A', 'J/kg/s', 'Convective tendency s by downdraft')
453 : call addfld('qvten_d_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qv by downdraft')
454 : call addfld('qlten_d_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency ql by downdraft')
455 : call addfld('qiten_d_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qi by downdraft')
456 : call addfld('nlten_d_SP', (/ 'lev' /), 'A', '1/kg/s', 'Convective tendency nl by downdraft')
457 : call addfld('niten_d_SP', (/ 'lev' /), 'A', '1/kg/s', 'Convective tendency ni by downdraft')
458 :
459 : call addfld('qtten_evp_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qt by evap of precip within environment')
460 : call addfld('slten_evp_SP', (/ 'lev' /), 'A', 'J/kg/s', 'Convective tendency sl by evap of precip within environment')
461 : call addfld('uten_evp_SP', (/ 'lev' /), 'A', 'm/s/s', 'Convective tendency u by evap of precip within environment')
462 : call addfld('vten_evp_SP', (/ 'lev' /), 'A', 'm/s/s', 'Convective tendency v by evap of precip within environment')
463 : call addfld('sten_evp_SP', (/ 'lev' /), 'A', 'J/kg/s', 'Convective tendency s by evap of precip within environment')
464 : call addfld('qvten_evp_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qv by evap of precip within environment')
465 : call addfld('qlten_evp_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency ql by evap of precip within environment')
466 : call addfld('qiten_evp_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qi by evap of precip within environment')
467 : call addfld('nlten_evp_SP', (/ 'lev' /), 'A', '#/kg/s', 'Convective tendency nl by evap of precip within environment')
468 : call addfld('niten_evp_SP', (/ 'lev' /), 'A', '#/kg/s', 'Convective tendency ni by evap of precip within environment')
469 :
470 : call addfld('qtten_dis_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qt by dissipative heating')
471 : call addfld('slten_dis_SP', (/ 'lev' /), 'A', 'J/kg/s', 'Convective tendency sl by dissipative heating')
472 : call addfld('uten_dis_SP', (/ 'lev' /), 'A', 'm/s/s', 'Convective tendency u by dissipative heating')
473 : call addfld('vten_dis_SP', (/ 'lev' /), 'A', 'm/s/s', 'Convective tendency v by dissipative heating')
474 : call addfld('sten_dis_SP', (/ 'lev' /), 'A', 'J/kg/s', 'Convective tendency s by dissipative heating')
475 : call addfld('qvten_dis_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qv by dissipative heating')
476 : call addfld('qlten_dis_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency ql by dissipative heating')
477 : call addfld('qiten_dis_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qi by dissipative heating')
478 : call addfld('nlten_dis_SP', (/ 'lev' /), 'A', '1/kg/s', 'Convective tendency nl by dissipative heating')
479 : call addfld('niten_dis_SP', (/ 'lev' /), 'A', '1/kg/s', 'Convective tendency ni by dissipative heating')
480 :
481 : call addfld('qtten_pos_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qt by positive tracer constraints')
482 : call addfld('slten_pos_SP', (/ 'lev' /), 'A', 'J/kg/s', 'Convective tendency sl by positive tracer constraints')
483 : call addfld('uten_pos_SP', (/ 'lev' /), 'A', 'm/s/s', 'Convective tendency u by positive tracer constraints')
484 : call addfld('vten_pos_SP', (/ 'lev' /), 'A', 'm/s/s', 'Convective tendency v by positive tracer constraints')
485 : call addfld('sten_pos_SP', (/ 'lev' /), 'A', 'J/kg/s', 'Convective tendency s by positive tracer constraints')
486 : call addfld('qvten_pos_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qv by positive tracer constraints')
487 : call addfld('qlten_pos_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency ql by positive tracer constraints')
488 : call addfld('qiten_pos_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qi by positive tracer constraints')
489 : call addfld('nlten_pos_SP', (/ 'lev' /), 'A', '1/kg/s', 'Convective tendency nl by positive tracer constraints')
490 : call addfld('niten_pos_SP', (/ 'lev' /), 'A', '1/kg/s', 'Convective tendency ni by positive tracer constraints')
491 :
492 : call addfld('qlten_sub_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency ql by compensating subsidence')
493 : call addfld('qiten_sub_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qi by compensating subsidence')
494 :
495 : call addfld('qlten_det_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency ql by detrainment')
496 : call addfld('qiten_det_SP', (/ 'lev' /), 'A', 'kg/kg/s', 'Convective tendency qi by detrainment')
497 :
498 : call addfld('CMFR_DET', (/ 'lev' /), 'A', 'kg/m2/s', 'Detrained convective mass flux' )
499 : call addfld('QLR_DET', (/ 'lev' /), 'A', 'kg/kg', 'Detrained convective LWC' )
500 : call addfld('QIR_DET', (/ 'lev' /), 'A', 'kg/kg', 'Detrained convective IWC' )
501 :
502 : !f call addfld('exit_Cu_SP', horiz_only, 'A', '1', 'Exit identifier of UNICON')
503 : call addfld('cush_SP', horiz_only, 'A', 'm', 'Cumulus top height from UNICON')
504 : call addfld('cushavg_SP', horiz_only, 'A', 'm', 'Mean cumulus top height from UNICON')
505 : call addfld('cuorg_SP', horiz_only, 'A', '1', 'Cumulus organization parameter from UNICON')
506 : call addfld('Radius_SP', horiz_only, 'A', 'm', 'Cumulus plume radius at surface from UNICON')
507 : !d call addfld('orgforce_SP', horiz_only, 'A', 'kg/m/s^2', 'Various organization forcing of UNICON')
508 : !d call addfld('tau_org_SP', horiz_only, 'A', 's', 'Damping time scale of convective organization of UNICON')
509 : !d call addfld('tau_TKE_SP', horiz_only, 'A', 's', 'Damping time scale of meso-scale TKE of UNICON')
510 : call addfld('sgh30_SP', horiz_only, 'A', 'm', 'Standard deviation of subgrid topography at 30 sec')
511 :
512 :
513 : call addfld('kw_SP', horiz_only, 'A', 'no', 'Internally computed kw from SPARKCONV' )
514 :
515 : call addfld('sigma_w_SP', horiz_only, 'A', 'm/s', 'Standard deviation of updraft w at surface from UNICON')
516 : call addfld('sigma_thl_SP', horiz_only, 'A', 'K', 'Standard deviation of updraft thl at surface from UNICON')
517 : call addfld('sigma_qt_SP', horiz_only, 'A', 'g/kg', 'Standard deviation of updraft qt at surface from UNICON')
518 : call addfld('sigma_u_SP', horiz_only, 'A', 'm/s', 'Standard deviation of updraft u at surface from UNICON')
519 : call addfld('sigma_v_SP', horiz_only, 'A', 'm/s', 'Standard deviation of updraft v at surface from UNICON')
520 :
521 : call addfld('w_org_SP', horiz_only, 'A', 'm2/s2', 'Organization-generated additional w at surface from UNICON')
522 : call addfld('thl_org_SP', horiz_only, 'A', 'K', 'Organization-generated additional thl at surface from UNICON')
523 : call addfld('qt_org_SP', horiz_only, 'A', 'g/kg', 'Organization-generated additional qt at surface from UNICON')
524 : call addfld('u_org_SP', horiz_only, 'A', 'm/s', 'Organization-generated additional u at surface from UNICON')
525 : call addfld('v_org_SP', horiz_only, 'A', 'm/s', 'Organization-generated additional v at surface from UNICON')
526 :
527 : call addfld('tkes_SP', horiz_only, 'A', 'm2/s2', 'TKE at surface within UNICON')
528 : call addfld('went_SP', horiz_only, 'A', 'm/s', 'Entrainment rate at the PBL top interface from UW PBL')
529 : call addfld('went_eff_SP', horiz_only, 'A', 'm/s', 'Effective entrainment rate at the PBL top interface in UNICON')
530 :
531 : call addfld('am_u_SP', (/ 'lev' /), 'A', '1', 'Convective updraft fractional area at mid-layer')
532 : call addfld('am_d_SP', (/ 'lev' /), 'A', '1', 'Convective downdraft fractional area at mid-layer')
533 :
534 : call addfld('qlm_u_SP', (/ 'lev' /), 'A', 'kg/kg', 'Area-weighted in-cumulus LWC condensate of updraft at mid-layer')
535 : call addfld('qlm_d_SP', (/ 'lev' /), 'A', 'kg/kg', 'Area-weighted in-cumulus IWC condensate of downdraft at mid-layer')
536 :
537 : call addfld('qim_u_SP', (/ 'lev' /), 'A', 'kg/kg', 'Area-weighted in-cumulus LWC condensate of updraft at mid-layer')
538 : call addfld('qim_d_SP', (/ 'lev' /), 'A', 'kg/kg', 'Area-weighted in-cumulus IWC condensate of downdraft at mid-layer')
539 :
540 : call addfld('thl_u_SP', (/ 'ilev' /), 'A', 'K', 'Mass-flux weighted updraft mean thl')
541 : call addfld('qt_u_SP', (/ 'ilev' /), 'A', 'kg/kg', 'Mass-flux weighted updraft mean qt')
542 : call addfld('u_u_SP', (/ 'ilev' /), 'A', 'm/s', 'Mass-flux weighted updraft mean u')
543 : call addfld('v_u_SP', (/ 'ilev' /), 'A', 'm/s', 'Mass-flux weighted updraft mean v')
544 : call addfld('w_u_SP', (/ 'ilev' /), 'A', 'm/s', 'Mass-flux weighted updraft mean w')
545 : call addfld('ql_u_SP', (/ 'ilev' /), 'A', 'kg/kg', 'Mass-flux weighted updraft mean ql')
546 : call addfld('qi_u_SP', (/ 'ilev' /), 'A', 'kg/kg', 'Mass-flux weighted updraft mean qi')
547 : call addfld('wa_u_SP', (/ 'ilev' /), 'A', 'm/s', 'Area weighted updraft mean w')
548 : call addfld('qla_u_SP', (/ 'ilev' /), 'A', 'kg/kg', 'Area weighted updraft mean ql')
549 : call addfld('qia_u_SP', (/ 'ilev' /), 'A', 'kg/kg', 'Area weighted updraft mean qi')
550 : call addfld('a_u_SP', (/ 'ilev' /), 'A', '1', 'Convective updraft fractional area')
551 : call addfld('rad_u_SP', (/ 'ilev' /), 'A', 'm', 'Number-weighted effective radius of updraft plumes')
552 : call addfld('num_u_SP', (/ 'ilev' /), 'A', '1/m^2', 'Number concentration of updraft plumes')
553 : call addfld('gamw_u_SP', (/ 'ilev' /), 'A', 'ratio', 'The ratio of w_u to wa_u')
554 : call addfld('nl_u_SP', (/ 'ilev' /), 'A', '1/kg', 'Mass-flux weighted updraft mean nl')
555 : call addfld('ni_u_SP', (/ 'ilev' /), 'A', '1/kg', 'Mass-flux weighted updraft mean ni')
556 : call addfld('thva_u_SP', (/ 'ilev' /), 'A', 'K', 'Area weighted updraft mean thv')
557 :
558 : call addfld('thl_d_SP', (/ 'ilev' /), 'A', 'K', 'Mass-flux weighted downdraft mean thl')
559 : call addfld('qt_d_SP', (/ 'ilev' /), 'A', 'kg/kg', 'Mass-flux weighted downdraft mean qt')
560 : call addfld('u_d_SP', (/ 'ilev' /), 'A', 'm/s', 'Mass-flux weighted downdraft mean u')
561 : call addfld('v_d_SP', (/ 'ilev' /), 'A', 'm/s', 'Mass-flux weighted downdraft mean v')
562 : call addfld('w_d_SP', (/ 'ilev' /), 'A', 'm/s', 'Mass-flux weighted downdraft mean w')
563 : call addfld('ql_d_SP', (/ 'ilev' /), 'A', 'kg/kg', 'Mass-flux weighted downdraft mean ql')
564 : call addfld('qi_d_SP', (/ 'ilev' /), 'A', 'kg/kg', 'Mass-flux weighted downdraft mean qi')
565 : call addfld('wa_d_SP' , (/ 'ilev' /), 'A', 'm/s', 'Area weighted downdraft mean w')
566 : call addfld('qla_d_SP', (/ 'ilev' /), 'A', 'kg/kg', 'Area weighted downdraft mean ql')
567 : call addfld('qia_d_SP', (/ 'ilev' /), 'A', 'kg/kg', 'Area weighted downdraft mean qi')
568 : call addfld('a_d_SP', (/ 'ilev' /), 'A', '1', 'Convective downdraft fractional area')
569 : call addfld('nl_d_SP', (/ 'ilev' /), 'A', '1/kg', 'Mass-flux weighted downdraft mean nl')
570 : call addfld('ni_d_SP', (/ 'ilev' /), 'A', '1/kg', 'Mass-flux weighted downdraft mean ni')
571 :
572 : call addfld('thv_b_SP', (/ 'ilev' /), 'A', 'K', 'thv_b : Environmental buoyancy at the base interface')
573 : call addfld('thv_t_SP', (/ 'ilev' /), 'A', 'K', 'thv_t : Environmental buoyancy at the top interface')
574 : call addfld('thv_mt_SP', (/ 'ilev' /), 'A', 'K', 'thv_mt : Environmental buoyancy at the top interface of lower layer')
575 : call addfld('thv_min_SP', (/ 'ilev' /), 'A', 'K', 'thv_min : Minimum environmental buoyancy for downdraft sorting')
576 :
577 : !a call addfld('CFL_SP', (/ 'lev' /), 'A', '1', 'Numerical stability parameter of UNICON')
578 :
579 : call addfld('cu_cmfr_SP', (/ 'lev' /), 'A', 'kg/m2/s', 'Mass flux of mixing environmental airs')
580 : call addfld('cu_thlr_SP', (/ 'lev' /), 'A', 'K', 'thl of mixing environmental airs')
581 : call addfld('cu_qtr_SP', (/ 'lev' /), 'A', 'kg/kg', 'qt of mixing environmental airs')
582 : call addfld('cu_qlr_SP', (/ 'lev' /), 'A', 'kg/kg', 'ql of mixing environmental airs')
583 : call addfld('cu_qir_SP', (/ 'lev' /), 'A', 'kg/kg', 'qi of mixing environmental airs')
584 : call addfld('cu_ur_SP', (/ 'lev' /), 'A', 'm/s', 'u of mixing environmental airs')
585 : call addfld('cu_vr_SP', (/ 'lev' /), 'A', 'm/s', 'v of mixing environmental airs')
586 : call addfld('cu_thvr_SP', (/ 'lev' /), 'A', 'K', 'thv of mixing environmental airs')
587 : call addfld('cu_rhr_SP', (/ 'lev' /), 'A', '1', 'rh of mixing environmental airs')
588 : call addfld('cu_nlr_SP', (/ 'lev' /), 'A', '1/kg', 'nl of mixing environmental airs')
589 : call addfld('cu_nir_SP', (/ 'lev' /), 'A', '1/kg', 'ni of mixing environmental airs')
590 :
591 : call addfld('a_p_SP', (/ 'ilev' /), 'A', '1', 'Convective precipitation area')
592 : call addfld('am_evp_SP', (/ 'lev' /), 'A', '1', 'Convective evaporation area')
593 : call addfld('am_pu_SP', (/ 'lev' /), 'A', '1', 'Overlapping area between conv precipitation and sat updraft area')
594 : call addfld('x_um_SP', (/ 'lev' /), 'A', 'm', 'Zonal displacement of the updraft area from the surface')
595 : call addfld('y_um_SP', (/ 'lev' /), 'A', 'm', 'Meridional displacement of the updraft area from the surface')
596 : call addfld('x_p_SP', (/ 'ilev' /), 'A', 'm', 'Zonal displacement of the precipitation area from the surface')
597 : call addfld('y_p_SP', (/ 'ilev' /), 'A', 'm', 'Meridional displacement of the precipitation area from the surface')
598 :
599 : do l = 1, pcnst
600 :
601 : if (cnst_is_mam_num(l) .or. cnst_is_mam_mmr(l)) then
602 :
603 : units = '1/kg/s'
604 : if (cnst_is_mam_mmr(l)) units = 'kg/kg/s'
605 :
606 : varname = trim(cnst_name(l))//'_u_SP'
607 : surname = trim(cnst_name(l))//' tendency by convective updraft from UNICON'
608 : call addfld(trim(varname), (/ 'lev' /), 'A', units, trim(surname))
609 :
610 : varname = trim(cnst_name(l))//'_d_SP'
611 : surname = trim(cnst_name(l))//' tendency by convective downdraft from UNICON'
612 : call addfld(trim(varname), (/ 'lev' /), 'A', units, trim(surname))
613 :
614 : varname = trim(cnst_name(l))//'_evp_SP'
615 : surname = trim(cnst_name(l))//' tendency by evap. of precip in env. from UNICON'
616 : call addfld(trim(varname), (/ 'lev' /), 'A', units, trim(surname))
617 :
618 : varname = trim(cnst_name(l))//'_dis_SP'
619 : surname = trim(cnst_name(l))//' tendency by dissipative heating from UNICON'
620 : call addfld(trim(varname), (/ 'lev' /), 'A', units, trim(surname))
621 :
622 : varname = trim(cnst_name(l))//'_pos_SP'
623 : surname = trim(cnst_name(l))//' tendency by positive moisture from UNICON'
624 : call addfld(trim(varname), (/ 'lev' /), 'A', units, trim(surname))
625 :
626 : end if
627 : end do
628 :
629 :
630 : ! Nov.15.2012. Below output corresponding to individual updraft segment is designed to write out individual
631 : ! segment values for writing UNICON_II paper.
632 :
633 : do msfc = 1, nseg
634 : write(numcha,'(i2.2)') msfc
635 :
636 : ! The properties of individual updraft segment
637 :
638 : call addfld('thl_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'K', numcha//' updraft segment : updraft thl' )
639 : call addfld('qt_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'kg/kg', numcha//' updraft segment : updraft qt' )
640 : call addfld('u_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'm/s', numcha//' updraft segment : updraft u' )
641 : call addfld('v_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'm/s', numcha//' updraft segment : updraft v' )
642 : call addfld('w_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'm/s', numcha//' updraft segment : updraft w' )
643 : call addfld('ql_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'kg/kg', numcha//' updraft segment : updraft ql' )
644 : call addfld('qi_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'kg/kg', numcha//' updraft segment : updraft qi' )
645 : call addfld('cmf_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'kg/s/m^2', numcha//' updraft segment : updraft cmf' )
646 : call addfld('a_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1', numcha//' updraft segment : updraft fractional area')
647 : call addfld('num_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1/m^2', numcha//' updraft segment : updraft number density' )
648 : call addfld('rad_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'm', numcha//' updraft segment : updraft plume radius' )
649 : call addfld('nl_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1/kg', numcha//' updraft segment : updraft nl' )
650 : call addfld('ni_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1/kg', numcha//' updraft segment : updraft ni' )
651 :
652 : call addfld('eps0_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1/Pa', numcha//' updraft segment : updraft eps0' )
653 : call addfld('eps_u'//numcha//'_SP' , (/ 'ilev' /), 'A', '1/Pa', numcha//' updraft segment : updraft eps' )
654 : call addfld('del_u'//numcha//'_SP' , (/ 'ilev' /), 'A', '1/Pa', numcha//' updraft segment : updraft del' )
655 : call addfld('eeps_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1', numcha//' updraft segment : updraft eeps' )
656 : call addfld('ddel_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1', numcha//' updraft segment : updraft ddel' )
657 : call addfld('xc_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1', numcha//' updraft segment : updraft xc' )
658 : call addfld('xs_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1', numcha//' updraft segment : updraft xs' )
659 : call addfld('xemin_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1', numcha//' updraft segment : updraft xemin' )
660 : call addfld('xemax_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1', numcha//' updraft segment : updraft xemax' )
661 : call addfld('cridis_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'm', numcha//' updraft segment : updraft cridis' )
662 : call addfld('thvcuenv_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'K', numcha//' updraft segment : updraft thvcuenv')
663 : call addfld('thvegenv_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'K', numcha//' updraft segment : updraft thvegenv')
664 : call addfld('thvxsenv_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'K', numcha//' updraft segment : updraft thvxsenv')
665 : call addfld('fmix_u'//numcha//'_SP', (/ 'ilev' /), 'A', '1', numcha//' updraft segment : updraft fmix' )
666 : call addfld('cmfumix_u'//numcha//'_SP', (/ 'ilev' /), 'A', 'kg/s/m^2', numcha//' updraft segment : updraft cmfumix' )
667 :
668 : ! call addfld('ktop'//numcha//'_SP', horiz_only, 'A', '1', numcha//' updraft segment : top layer index')
669 : call addfld('ptop'//numcha//'_SP', horiz_only, 'A', 'Pa', numcha//' updraft segment : updraft top pressure')
670 : call addfld('ztop'//numcha//'_SP', horiz_only, 'A', 'm', numcha//' updraft segment : updraft top height')
671 :
672 : ! The properties of mass flux weighted ( or area-weighted or net=sum ) downdraft properties for individual updraft segment
673 :
674 : call addfld('thl_d'//numcha//'_SP', (/ 'ilev' /), 'A', 'K',&
675 : 'Mass-flux weighted mean downdraft thl for '// numcha//' updraft segment')
676 : call addfld('qt_d'//numcha//'_SP' , (/ 'ilev' /), 'A', 'kg/kg',&
677 : 'Mass-flux weighted mean downdraft qt for '// numcha//' updraft segment')
678 : call addfld('u_d'//numcha//'_SP' , (/ 'ilev' /), 'A', 'm/s',&
679 : 'Mass-flux weighted mean downdraft u for '// numcha//' updraft segment')
680 : call addfld('v_d'//numcha//'_SP' , (/ 'ilev' /), 'A', 'm/s',&
681 : 'Mass-flux weighted mean downdraft v for '// numcha//' updraft segment')
682 : call addfld('w_d'//numcha//'_SP' , (/ 'ilev' /), 'A', 'm/s',&
683 : 'Mass-flux weighted mean downdraft w for '// numcha//' updraft segment')
684 : call addfld('ql_d'//numcha//'_SP' , (/ 'ilev' /), 'A', 'kg/kg',&
685 : 'Mass-flux weighted mean downdraft ql for '// numcha//' updraft segment')
686 : call addfld('qi_d'//numcha//'_SP' , (/ 'ilev' /), 'A', 'kg/kg',&
687 : 'Mass-flux weighted mean downdraft qi for '// numcha//' updraft segment')
688 : call addfld('wa_d'//numcha//'_SP' , (/ 'ilev' /), 'A', 'm/s',&
689 : 'Area weighted mean downdraft w for '// numcha//' updraft segment')
690 : call addfld('qla_d'//numcha//'_SP', (/ 'ilev' /), 'A', 'kg/kg',&
691 : 'Area weighted mean downdraft ql for '// numcha//' updraft segment')
692 : call addfld('qia_d'//numcha//'_SP', (/ 'ilev' /), 'A', 'kg/kg',&
693 : 'Area weighted mean downdraft qi for '// numcha//' updraft segment')
694 : call addfld('cmf_d'//numcha//'_SP', (/ 'ilev' /), 'A', 'kg/s/m^2',&
695 : 'Net downdraft cmf for '// numcha//' updraft segment')
696 : call addfld('a_d'//numcha//'_SP' , (/ 'ilev' /), 'A', 'fraction',&
697 : 'Net downdraft a for '// numcha//' updraft segment')
698 : call addfld('nl_d'//numcha//'_SP' , (/ 'ilev' /), 'A', '#/kg',&
699 : 'Mass-flux weighted mean downdraft nl for '// numcha//' updraft segment')
700 : call addfld('ni_d'//numcha//'_SP' , (/ 'ilev' /), 'A', '#/kg',&
701 : 'Mass-flux weighted mean downdraft ni for '// numcha//' updraft segment')
702 :
703 : enddo
704 :
705 : ! Nov.16.2012. Additional detailed diagnostic output
706 :
707 : call addfld('thl_orgfce_SP', horiz_only, 'A', 'K/s', &
708 : 'Total organization forcing generating thl difference between non-wake and grid-mean areas')
709 : call addfld('qt_orgfce_SP', horiz_only, 'A', 'kg/kg/s', &
710 : 'Total organization forcing generating qt difference between non-wake and grid-mean areas')
711 : call addfld('u_orgfce_SP', horiz_only, 'A', 'm/s/s', &
712 : 'Total organization forcing generating u difference between non-wake and grid-mean areas')
713 : call addfld('v_orgfce_SP', horiz_only, 'A', 'm/s/s', &
714 : 'Total organization forcing generating v difference between non-wake and grid-mean areas')
715 : call addfld('awk_orgfce_SP', horiz_only, 'A', '1/s', &
716 : 'Total organization forcing generating wake area')
717 :
718 : call addfld('thl_orgfce_f_SP', horiz_only, 'A', 'K/s', &
719 : 'PBL top flux-related forcing generating thl difference between non-wake and grid-mean areas')
720 : call addfld('qt_orgfce_f_SP', horiz_only, 'A', 'kg/kg/s', &
721 : 'PBL top flux-related forcing generating qt difference between non-wake and grid-mean areas')
722 : call addfld('u_orgfce_f_SP', horiz_only, 'A', 'm/s/s', &
723 : 'PBL top flux-related forcing generating u difference between non-wake and grid-mean areas')
724 : call addfld('v_orgfce_f_SP', horiz_only, 'A', 'm/s/s', &
725 : 'PBL top flux-related forcing generating v difference between non-wake and grid-mean areas')
726 : call addfld('awk_orgfce_f_SP', horiz_only, 'A', '1/s', &
727 : 'PBL top flux-related forcing generating wake area')
728 :
729 : call addfld('thl_orgfce_u_SP', horiz_only, 'A', 'K/s', &
730 : 'Up-and-Down diabatic forcing generating thl difference between non-wake and grid-mean areas')
731 : call addfld('qt_orgfce_u_SP', horiz_only, 'A', 'kg/kg/s', &
732 : 'Up-and-Down diabatic forcing generating qt difference between non-wake and grid-mean areas')
733 : call addfld('u_orgfce_u_SP', horiz_only, 'A', 'm/s/s', &
734 : 'Up-and-Down diabatic forcing generating u difference between non-wake and grid-mean areas')
735 : call addfld('v_orgfce_u_SP', horiz_only, 'A', 'm/s/s', &
736 : 'Up-and-Down diabatic forcing generating v difference between non-wake and grid-mean areas')
737 : call addfld('awk_orgfce_m_SP', horiz_only, 'A', '1/s', &
738 : 'Lateral-Mixing forcing for wake area')
739 :
740 : call addfld('thl_orgfce_e_SP', horiz_only, 'A', 'K/s', &
741 : 'Environment diabatic forcing generating thl difference between non-wake and grid-mean areas')
742 : call addfld('qt_orgfce_e_SP', horiz_only, 'A', 'kg/kg/s', &
743 : 'Environment diabatic forcing generating qt difference between non-wake and grid-mean areas')
744 : call addfld('u_orgfce_e_SP', horiz_only, 'A', 'm/s/s', &
745 : 'Environment diabatic forcing generating u difference between non-wake and grid-mean areas')
746 : call addfld('v_orgfce_e_SP', horiz_only, 'A', 'm/s/s', &
747 : 'Environment diabatic forcing generating v difference between non-wake and grid-mean areas')
748 : call addfld('cmf_d_orgh_SP', horiz_only, 'A', 'kg/m^2/s', &
749 : 'Organization-inducing downdraft mass flux at the PBL top interface')
750 :
751 : call addfld('taui_thl_SP', horiz_only, 'A', '1/s', &
752 : 'Inverse of damping time scale of the difference between off-wake and grid-mean thl')
753 : call addfld('taui_qt_SP', horiz_only, 'A', '1/s', &
754 : 'Inverse of damping time scale of the difference between off-wake and grid-mean qt')
755 : call addfld('taui_u_SP', horiz_only, 'A', '1/s', &
756 : 'Inverse of damping time scale of the difference between off-wake and grid-mean u')
757 : call addfld('taui_v_SP', horiz_only, 'A', '1/s', &
758 : 'Inverse of damping time scale of the difference between off-wake and grid-mean v')
759 : call addfld('taui_awk_SP', horiz_only, 'A', '1/s', &
760 : 'Inverse of damping time scale of the wake area')
761 :
762 : call addfld('del_org_SP', horiz_only, 'A', '1/s', &
763 : 'Detrainment rate of the cold pool from UNICON')
764 : call addfld('del0_org_SP', horiz_only, 'A', '1/s', &
765 : 'Effective detrainment rate of the cold pool from UNICON')
766 :
767 : #endif
768 :
769 0 : end subroutine unicon_cam_init
770 :
771 : !==================================================================================================
772 :
773 0 : function unicon_implements_cnst(name)
774 :
775 : ! Return true if specified constituent is implemented by this package
776 :
777 : character(len=*), intent(in) :: name ! constituent name
778 : logical :: unicon_implements_cnst ! return value
779 :
780 : integer :: m
781 : !-----------------------------------------------------------------------
782 :
783 0 : unicon_implements_cnst = .false.
784 :
785 : #ifdef USE_UNICON
786 :
787 : do m = 1, n_org
788 : if (name == cnst_names(m)) then
789 : unicon_implements_cnst = .true.
790 : return
791 : end if
792 : end do
793 :
794 : #endif
795 :
796 0 : end function unicon_implements_cnst
797 :
798 : !==================================================================================================
799 :
800 0 : subroutine unicon_init_cnst(name, latvals, lonvals, mask, q)
801 :
802 : ! Initialize constituents if they are not read from the initial file
803 :
804 : character(len=*), intent(in) :: name ! constituent name
805 : real(r8), intent(in) :: latvals(:) ! lat in degrees (ncol)
806 : real(r8), intent(in) :: lonvals(:) ! lon in degrees (ncol)
807 : logical, intent(in) :: mask(:) ! Only initialize where .true.
808 : real(r8), intent(out) :: q(:,:) ! kg tracer/kg dry air (gcol, plev
809 : !-----------------------------------------------------------------------
810 : integer :: k
811 :
812 : #ifdef USE_UNICON
813 :
814 : do k = 1, size(q, 2)
815 : if ( name == 'ORGawk' ) then
816 : where(mask)
817 : q(:,k) = 0.0_r8
818 : end where
819 : else if ( name == 'ORGthl' ) then
820 : where(mask)
821 : q(:,k) = 100.0_r8
822 : end where
823 : else if ( name == 'ORGqto' ) then
824 : where(mask)
825 : q(:,k) = 100.0_r8
826 : end where
827 : else if ( name == 'ORGuoo' ) then
828 : where(mask)
829 : q(:,k) = 100.0_r8
830 : end where
831 : else if ( name == 'ORGvoo' ) then
832 : where(mask)
833 : q(:,k) = 100.0_r8
834 : end where
835 : end if
836 : end do
837 :
838 : #endif
839 :
840 0 : end subroutine unicon_init_cnst
841 :
842 : !==================================================================================================
843 :
844 0 : subroutine unicon_cam_tend(dt, state, cam_in, &
845 : pbuf, ptend, out)
846 :
847 :
848 : ! ---------------------- !
849 : ! Input-output Arguments !
850 : ! ---------------------- !
851 :
852 : real(r8), intent(in) :: dt ! Time step [s]
853 : type(physics_state), intent(in) :: state ! Physics state variables
854 : type(cam_in_t), intent(in) :: cam_in ! import state
855 : type(physics_buffer_desc), pointer :: pbuf(:) ! physics buffer
856 : type(physics_ptend), intent(out) :: ptend ! parameterization tendencies
857 : type(unicon_out_t), intent(out) :: out ! parameterization outputs
858 :
859 : ! -------------------------------------------------------- !
860 : ! Internal output and local variables for positive tracers !
861 : ! -------------------------------------------------------- !
862 :
863 : real(r8) :: sten_ori(mix,mkx) ! Tendency of dry static energy [ J / kg / s ]
864 : real(r8) :: qvten_ori(mix,mkx) ! Tendency of water vapor specific humidity [ kg / kg / s ]
865 : real(r8) :: qlten_ori(mix,mkx) ! Tendency of liquid water mixing ratio [ kg / kg / s ]
866 : real(r8) :: qiten_ori(mix,mkx) ! Tendency of ice mixing ratio [ kg / kg / s ]
867 : real(r8) :: trten_ori(mix,mkx,ncnst) ! Tendency of tracers [ # / kg / s, kg / kg / s ]
868 :
869 : real(r8) :: slten_pos_inv(mix,mkx) !
870 : real(r8) :: qtten_pos_inv(mix,mkx) !
871 : real(r8) :: uten_pos_inv(mix,mkx) !
872 : real(r8) :: vten_pos_inv(mix,mkx) !
873 : real(r8) :: sten_pos_inv(mix,mkx) !
874 : real(r8) :: qvten_pos_inv(mix,mkx) !
875 : real(r8) :: qlten_pos_inv(mix,mkx) !
876 : real(r8) :: qiten_pos_inv(mix,mkx) !
877 : real(r8) :: trten_pos_inv(mix,mkx,ncnst)
878 :
879 : ! --------------- !
880 : ! Local variables !
881 : ! --------------- !
882 :
883 : integer :: iend
884 : integer :: lchnk
885 : integer :: itim
886 :
887 : ! fields in physics buffer
888 : real(r8), pointer, dimension(:) :: & ! (mix)
889 : sgh30 ! Standard deviation of topography at 30s [m]
890 :
891 : real(r8), pointer, dimension(:,:) :: & ! (mix,mkx)
892 : ast0_inv ! Physical stratus fraction [ fraction ]
893 :
894 : real(r8), pointer, dimension(:,:) :: & ! (mix,mkx+1)
895 : tke0_inv, &! TKE at the interface [ m2/s2 ]
896 : bprod0_inv ! Buoyancy production at the interface [ m2/s3 ]
897 :
898 : integer(i4), pointer, dimension(:) :: & ! (mix)
899 : kpblh_inv ! Layer index with PBL top in it or at the base interface
900 :
901 : real(r8), pointer, dimension(:) :: & ! (mix)
902 : pblh, &! PBL top height [m]
903 : went, &! Entrainment rate at the PBL top interface directly from UW PBL scheme [ m / s ]. went = 0 for STL.
904 : cush, &! Cumulus top height [ m ]
905 : cushavg, &! Mean cumulus top height weighted by updraft mass flux at surface [ m ]
906 : cuorg, &! Convective organization parameter [ 0-1 ]
907 :
908 : awk_PBL, &! Wake area within PBL [ 0 - 1 ]
909 : delta_thl_PBL, &! Diff of thl between off-wake region and grid-mean value averaged over the PBL [ K ]
910 : delta_qt_PBL, &! Diff of qt between off-wake region and grid-mean value averaged over the PBL [ kg/kg ]
911 : delta_u_PBL, &! Diff of u between off-wake region and grid-mean value averaged over the PBL [ m/s ]
912 : delta_v_PBL ! Diff of v between off-wake region and grid-mean value averaged over the PBL [ m/s ]
913 :
914 : real(r8), pointer, dimension(:,:) :: & ! (mix,ncnst)
915 : delta_tr_PBL ! Diff of tr between off-wake region and grid-mean value avg over the PBL [ kg/kg, #/kg ]
916 :
917 : real(r8), dimension(mix,mkx) :: & ! (mix,mkx)
918 : cu_cmfum ! The mass involved in the updraft buoyancy sorting at the previous time step [ kg/s/m2 ]
919 :
920 : real(r8), pointer, dimension(:,:) :: & ! (mix,mkx)
921 : cu_cmfr, &! The detrained mass from convective up and downdraft at the previous time step [ kg/s/m2 ]
922 : cu_thlr, &! Mass-flux wghted mean 'thl' of detrained mass from conv up and downdraft at prev step [ K ]
923 : cu_qtr, &! Mass-flux wghted mean 'qt' of detrained mass from conv up and downdraft at prev step [ kg/kg ]
924 : cu_ur, &! Mass-flux wghted mean 'u' of detrained mass from conv up and downdraft at prev step [ m/s ]
925 : cu_vr, &! Mass-flux wghted mean 'v' of detrained mass from conv up and downdraft at prev step [ m/s ]
926 : cu_qlr, &! Mass-flux wghted mean 'ql' of detrained mass from conv up and downdraft at prev step [ kg/kg ]
927 : cu_qir, &! Mass-flux wghted mean 'qi' of detrained mass from conv up and downdraft at prev step [ kg/kg ]
928 : cmfr_det,&! The detrained mass from convective up and downdraft at the previous time step [ kg/s/m2 ]
929 : qlr_det, &! Mass-flux wghted mean 'ql' of detrained mass from conv up and downdraft at prev step [ kg/kg ]
930 : qir_det ! Mass-flux wghted mean 'qi' of detrained mass from conv up and downdraft at prev step [ kg/kg ]
931 :
932 : real(r8), pointer, dimension(:,:,:) :: & ! (mix,mkx,ncnst)
933 : cu_trr ! Mass-flux wghted mean 'tr' of detrained mass from conv up and downdraft at prev step [ kg/kg ]
934 :
935 : real(r8), dimension(mix,mkx) :: & ! (mix,mkx)
936 : cu_cmfrd,&! The amount of detrained mass from convective downdraft at the previous time step [ kg/s/m2 ]
937 : cu_thlrd,&! Mass-flux wghted mean 'thl' of detrained mass from conv downdraft at prev step [ K ]
938 : cu_qtrd, &! Mass-flux wghted mean 'qt' of detrained mass from conv downdraft at prev step [ kg/kg ]
939 : cu_urd, &! Mass-flux wghted mean 'u' of detrained mass from conv downdraft at prev step [ m/s ]
940 : cu_vrd, &! Mass-flux wghted mean 'v' of detrained mass from conv downdraft at prev step [ m/s ]
941 : cu_qlrd, &! Mass-flux wghted mean 'ql' of detrained mass from conv downdraft at prev step [ kg/kg ]
942 : cu_qird ! Mass-flux wghted mean 'qi' of detrained mass from conv downdraft at prev step [ kg/kg ]
943 :
944 : real(r8), dimension(mix,mkx,ncnst) :: & ! (mix,mkx,ncnst)
945 : cu_trrd ! Mass-flux wghted mean 'tr' of detrained mass from conv downdraft at prev step [ kg/kg ]
946 :
947 : real(r8), pointer, dimension(:,:) :: & ! (mix,mkx)
948 : shfrc, &! Convective updraft fractional area
949 : icwmrsh, &! In-cloud LWC+IWC within convective updraft
950 : rprdsh, &! Prod of rain+snow by lateral expels of cumulus condensate [ kg / kg / s ]
951 : evapc_inv, &! Evaporation rate of convective precipitation within environment [ kg/kg/s ]
952 : am_evp_st_inv, &! Evaporation area of stratiform precipitation [fraction]
953 : evprain_st_inv, &! Grid-mean evaporation rate of stratiform rain [kg/kg/s] >= 0.
954 : evpsnow_st_inv ! Grid-mean evaporation rate of stratiform snow [kg/kg/s] >= 0.
955 :
956 : real(r8), pointer, dimension(:) :: & ! (mix)
957 : precip, &! Precipitation flux at surface in flux unit [ m / s ]
958 : snow ! Snow flux at surface in flux unit [ m / s ]
959 :
960 : real(r8) :: ps0(mix,0:mkx) ! Environmental pressure at full sigma levels
961 : real(r8) :: zs0(mix,0:mkx) ! Environmental height at full sigma levels
962 : real(r8) :: p0(mix,mkx) ! Environmental pressure at half sigma levels
963 : real(r8) :: z0(mix,mkx) ! Environmental height at half sigma levels
964 : real(r8) :: dp0(mix,mkx) ! Environmental layer pressure thickness
965 : real(r8) :: dpdry0(mix,mkx) ! Environmental layer dry pressure thickness
966 : real(r8) :: u0(mix,mkx) ! Environmental zonal wind
967 : real(r8) :: v0(mix,mkx) ! Environmental meridional wind
968 : real(r8) :: qv0(mix,mkx) ! Environmental specific humidity
969 : real(r8) :: ql0(mix,mkx) ! Environmental liquid water mixing ratio
970 : real(r8) :: qi0(mix,mkx) ! Environmental ice mixing ratio
971 : real(r8) :: tr0(mix,mkx,ncnst) ! Environmental tracers [ #/kg, kg/kg ]
972 : real(r8) :: t0(mix,mkx) ! Environmental temperature
973 : real(r8) :: s0(mix,mkx) ! Environmental dry static energy
974 : real(r8) :: ast0(mix,mkx) ! Physical stratiform cloud fraction [ fraction ]
975 : real(r8) :: tke0(mix,0:mkx) ! TKE [ m2/s2 ]
976 : real(r8) :: bprod0(mix,0:mkx) ! Buoyancy production [ m2/s3 ]
977 : real(r8) :: am_evp_st(mix,mkx) ! Evaporation area of stratiform precipitation [fraction]
978 : real(r8) :: evprain_st(mix,mkx) ! Grid-mean evaporation rate of stratiform rain [kg/kg/s] >= 0.
979 : real(r8) :: evpsnow_st(mix,mkx) ! Grid-mean evaporation rate of stratiform snow [kg/kg/s] >= 0.
980 :
981 : integer(i4) :: kpblh(mix) ! Layer index with PBL top in it or at the base interface
982 :
983 :
984 : real(r8) :: am_u(mix,mkx) ! Convective updraft fractional area
985 : real(r8) :: qlm_u(mix,mkx) ! In-cloud LWC within convective updraft [ kg / kg ]
986 : real(r8) :: qim_u(mix,mkx) ! In-cloud IWC within convective updraft [ kg / kg ]
987 :
988 : real(r8) :: am_d(mix,mkx) ! Convective downdraft fractional area
989 : real(r8) :: qlm_d(mix,mkx) ! In-cloud LWC within downdraft updraft [ kg / kg ]
990 : real(r8) :: qim_d(mix,mkx) ! In-cloud IWC within downdraft updraft [ kg / kg ]
991 :
992 : real(r8) :: cmf_u(mix,0:mkx) ! Upward convective mass flux at the interface [ kg / s / m2 ]
993 : real(r8) :: slflx(mix,0:mkx) ! Net upward convective flux of liquid static energy [ J / s / m2 ]
994 : real(r8) :: qtflx(mix,0:mkx) ! Net upward convective flux of total specific humidity [ kg / s / m2 ]
995 :
996 : real(r8) :: qvten(mix,mkx) ! Tendency of water vapor specific humidity [ kg / kg / s ]
997 : real(r8) :: qlten(mix,mkx) ! Tendency of liquid water mixing ratio [ kg / kg / s ]
998 : real(r8) :: qiten(mix,mkx) ! Tendency of ice mixing ratio [ kg / kg / s ]
999 : real(r8) :: trten(mix,mkx,ncnst) ! Tendency of tracers [ # / kg / s, kg / kg / s ]
1000 :
1001 : real(r8) :: sten(mix,mkx) ! Tendency of dry static energy [ J / kg / s ]
1002 : real(r8) :: uten(mix,mkx) ! Tendency of zonal wind [ m / s / s ]
1003 : real(r8) :: vten(mix,mkx) ! Tendency of meridional wind [ m / s / s ]
1004 :
1005 : real(r8) :: qrten(mix,mkx) ! Production rate of rain by lateral expels of cumulus condensate [ kg / kg / s ]
1006 : real(r8) :: qsten(mix,mkx) ! Production rate of snow by lateral expels of cumulus condensate [ kg / kg / s ]
1007 :
1008 : real(r8) :: evapc(mix,mkx) ! Evaporation rate of convective precipitation within environment [ kg/kg/s ]
1009 :
1010 : real(r8) :: rqc(mix,mkx) ! Production rate of detrained LWC+IWC [ kg / kg / s ] > 0
1011 : real(r8) :: rqc_l(mix,mkx) ! Production rate of detrained LWC [ kg / kg / s ] > 0
1012 : real(r8) :: rqc_i(mix,mkx) ! Production rate of detrained IWC [ kg / kg / s ] > 0
1013 : real(r8) :: rnc_l(mix,mkx) ! Production rate of detrained droplet number of cloud liquid droplets [ # / kg / s ] > 0
1014 : real(r8) :: rnc_i(mix,mkx) ! Production rate of detrained droplet number of cloud ice droplets [ # / kg / s ] > 0
1015 :
1016 : real(r8) :: cnt(mix) ! Cloud top interface index ( ki = kpen )
1017 : real(r8) :: cnb(mix) ! Cloud base interface index ( ki = krel-1 )
1018 :
1019 : real(r8) :: pdel0(mix,mkx) ! Environmental pressure thickness ( either dry or moist ) [ Pa ]
1020 : real(r8) :: trmin ! Minimum concentration of tracer [ # / kg ]
1021 :
1022 : real(r8) :: cmf_det(mix,mkx) ! Det mass flux from convective updraft (not from environment) and downdraft [kg/s/m^2]
1023 : real(r8) :: ql_det(mix,mkx) ! Det conv LWC from convective updraft (not from environment) and downdraft [kg/s/m^2]
1024 : real(r8) :: qi_det(mix,mkx) ! Det conv IWC from convective updraft (not from environment) and downdraft [kg/s/m^2]
1025 :
1026 : ! For prognostically updated state variables
1027 :
1028 : real(r8) :: qv0_c(mix,mkx) ! Environmental specific humidity
1029 : real(r8) :: ql0_c(mix,mkx) ! Environmental liquid water mixing ratio
1030 : real(r8) :: qi0_c(mix,mkx) ! Environmental ice mixing ratio
1031 : real(r8) :: t0_c(mix,mkx) ! Environmental temperature
1032 : real(r8) :: s0_c(mix,mkx) ! Environmental dry static energy
1033 : real(r8) :: tr0_c(mix,mkx,ncnst) ! Environmental tracers [ # / kg, kg / kg ]
1034 :
1035 : ! Layer index variables
1036 : integer :: k ! Vertical index for local fields
1037 : integer :: k_inv ! Vertical index for incoming fields
1038 : integer :: mt ! Tracer index [ no ]
1039 : integer :: m
1040 :
1041 : ! For aerosol tendency output
1042 : character(len=30) :: varname
1043 :
1044 : logical :: lq(ncnst)
1045 :
1046 : ! --------- !
1047 : ! Main body !
1048 : ! --------- !
1049 :
1050 : #ifdef USE_UNICON
1051 :
1052 : iend = state%ncol
1053 : lchnk = state%lchnk
1054 :
1055 : ! Associate pointers with physics buffer fields
1056 :
1057 : itim = pbuf_old_tim_idx()
1058 : call pbuf_get_field(pbuf, sgh30_idx, sgh30)
1059 : call pbuf_get_field(pbuf, ast_idx, ast0_inv, start=(/1,1,itim/), kount=(/pcols,pver,1/))
1060 : call pbuf_get_field(pbuf, tke_idx, tke0_inv)
1061 : call pbuf_get_field(pbuf, bprod_idx, bprod0_inv)
1062 : call pbuf_get_field(pbuf, kpblh_idx, kpblh_inv)
1063 : call pbuf_get_field(pbuf, pblh_idx, pblh)
1064 : call pbuf_get_field(pbuf, went_idx, went)
1065 : call pbuf_get_field(pbuf, cush_idx, cush, start=(/1,itim/), kount=(/pcols,1/))
1066 : call pbuf_get_field(pbuf, cushavg_idx, cushavg, start=(/1,itim/), kount=(/pcols,1/))
1067 : call pbuf_get_field(pbuf, cuorg_idx, cuorg, start=(/1,itim/), kount=(/pcols,1/))
1068 :
1069 : call pbuf_get_field(pbuf, awk_PBL_idx, awk_PBL, start=(/1,itim/), kount=(/pcols,1/))
1070 : call pbuf_get_field(pbuf, delta_thl_PBL_idx, delta_thl_PBL, start=(/1,itim/), kount=(/pcols,1/))
1071 : call pbuf_get_field(pbuf, delta_qt_PBL_idx, delta_qt_PBL, start=(/1,itim/), kount=(/pcols,1/))
1072 : call pbuf_get_field(pbuf, delta_u_PBL_idx, delta_u_PBL, start=(/1,itim/), kount=(/pcols,1/))
1073 : call pbuf_get_field(pbuf, delta_v_PBL_idx, delta_v_PBL, start=(/1,itim/), kount=(/pcols,1/))
1074 : call pbuf_get_field(pbuf, delta_tr_PBL_idx, delta_tr_PBL, start=(/1,1,itim/), kount=(/pcols,pcnst,1/))
1075 :
1076 : call pbuf_get_field(pbuf, cu_cmfr_idx, cu_cmfr, start=(/1,1,itim/), kount=(/pcols,pver,1/))
1077 : call pbuf_get_field(pbuf, cu_thlr_idx, cu_thlr, start=(/1,1,itim/), kount=(/pcols,pver,1/))
1078 : call pbuf_get_field(pbuf, cu_qtr_idx, cu_qtr, start=(/1,1,itim/), kount=(/pcols,pver,1/))
1079 : call pbuf_get_field(pbuf, cu_ur_idx, cu_ur, start=(/1,1,itim/), kount=(/pcols,pver,1/))
1080 : call pbuf_get_field(pbuf, cu_vr_idx, cu_vr, start=(/1,1,itim/), kount=(/pcols,pver,1/))
1081 : call pbuf_get_field(pbuf, cu_qlr_idx, cu_qlr, start=(/1,1,itim/), kount=(/pcols,pver,1/))
1082 : call pbuf_get_field(pbuf, cu_qir_idx, cu_qir, start=(/1,1,itim/), kount=(/pcols,pver,1/))
1083 : call pbuf_get_field(pbuf, cu_trr_idx, cu_trr, start=(/1,1,1,itim/), kount=(/pcols,pver,pcnst,1/))
1084 :
1085 : call pbuf_get_field(pbuf, shfrc_idx, shfrc)
1086 : call pbuf_get_field(pbuf, icwmrsh_idx, icwmrsh)
1087 : call pbuf_get_field(pbuf, rprdsh_idx, rprdsh)
1088 : call pbuf_get_field(pbuf, prec_sh_idx, precip)
1089 : call pbuf_get_field(pbuf, snow_sh_idx, snow)
1090 : call pbuf_get_field(pbuf, nevapr_shcu_idx, evapc_inv)
1091 : call pbuf_get_field(pbuf, am_evp_st_idx, am_evp_st_inv)
1092 : call pbuf_get_field(pbuf, evprain_st_idx, evprain_st_inv)
1093 : call pbuf_get_field(pbuf, evpsnow_st_idx, evpsnow_st_inv)
1094 :
1095 : call pbuf_get_field(pbuf, cmfr_det_idx, cmfr_det)
1096 : call pbuf_get_field(pbuf, qlr_det_idx, qlr_det)
1097 : call pbuf_get_field(pbuf, qir_det_idx, qir_det)
1098 :
1099 : ! Reverse variables defined at the layer mid-point
1100 :
1101 : do k = 1, mkx
1102 : k_inv = mkx + 1 - k
1103 : p0(:iend,k) = state%pmid(:iend,k_inv)
1104 : u0(:iend,k) = state%u(:iend,k_inv)
1105 : v0(:iend,k) = state%v(:iend,k_inv)
1106 : z0(:iend,k) = state%zm(:iend,k_inv)
1107 : dp0(:iend,k) = state%pdel(:iend,k_inv)
1108 : dpdry0(:iend,k) = state%pdeldry(:iend,k_inv)
1109 : qv0(:iend,k) = state%q(:iend,k_inv,1)
1110 : ql0(:iend,k) = state%q(:iend,k_inv,ixcldliq)
1111 : qi0(:iend,k) = state%q(:iend,k_inv,ixcldice)
1112 : t0(:iend,k) = state%t(:iend,k_inv)
1113 : s0(:iend,k) = state%s(:iend,k_inv)
1114 : ast0(:iend,k) = ast0_inv(:iend,k_inv)
1115 : am_evp_st(:iend,k) = am_evp_st_inv(:iend,k_inv)
1116 : evprain_st(:iend,k) = evprain_st_inv(:iend,k_inv)
1117 : evpsnow_st(:iend,k) = evpsnow_st_inv(:iend,k_inv)
1118 : do mt = 1, ncnst
1119 : tr0(:iend,k,mt) = state%q(:iend,k_inv,mt)
1120 : end do
1121 : end do
1122 :
1123 : ! Reverse variables defined at the interfaces
1124 :
1125 : do k = 0, mkx
1126 : k_inv = mkx + 1 - k
1127 : ps0(:iend,k) = state%pint(:iend,k_inv)
1128 : zs0(:iend,k) = state%zi(:iend,k_inv)
1129 : tke0(:iend,k) = tke0_inv(:iend,k_inv)
1130 : bprod0(:iend,k) = bprod0_inv(:iend,k_inv)
1131 : end do
1132 :
1133 : ! Reverse the index of ambiguous layer
1134 :
1135 : kpblh(:iend) = mkx + 1 - kpblh_inv(:iend)
1136 :
1137 :
1138 : call compute_unicon( mix , mkx , iend , ncnst , dt , &
1139 : ps0 , zs0 , p0 , z0 , dp0 , dpdry0 , &
1140 : t0 , qv0 , ql0 , qi0 , tr0 , &
1141 : u0 , v0 , ast0 , tke0 , bprod0 , &
1142 : kpblh , pblh , went , &
1143 : cam_in%cflx(:,1), cam_in%shf, cam_in%wsx, cam_in%wsy, cam_in%cflx , &
1144 : cam_in%landfrac, sgh30 , &
1145 : am_evp_st , evprain_st , evpsnow_st, &
1146 : cush , cushavg , cuorg , &
1147 : awk_PBL , delta_thl_PBL , delta_qt_PBL , &
1148 : delta_u_PBL , delta_v_PBL , delta_tr_PBL , &
1149 : cu_cmfum , cu_cmfr , cu_thlr , cu_qtr , cu_ur , cu_vr , &
1150 : cu_qlr , cu_qir , cu_trr , &
1151 : cu_cmfrd , cu_thlrd , cu_qtrd , cu_urd , cu_vrd , &
1152 : cu_qlrd , cu_qird , cu_trrd , &
1153 : am_u , qlm_u , qim_u , &
1154 : am_d , qlm_d , qim_d , &
1155 : cmf_u , slflx , qtflx , &
1156 : qvten , qlten , qiten , trten , &
1157 : sten , uten , vten , &
1158 : qrten , qsten , &
1159 : rqc_l , rqc_i , rqc , rnc_l , rnc_i , &
1160 : out%rliq , precip , snow , evapc , &
1161 : cnt , cnb , cmf_det , ql_det , qi_det , &
1162 : lchnk )
1163 :
1164 :
1165 : ! Initialize output ptend
1166 : lq(:) = .true.
1167 : call physics_ptend_init(ptend, state%psetcols, 'unicon', ls=.true., lu=.true., lv=.true., lq=lq)
1168 :
1169 : ! ----------------------------------------------------- !
1170 : ! Treatment of reserved liquid-ice water for CAM !
1171 : ! All the reserved condensate is converted into liquid. !
1172 : ! Jan.04.2012. Relocated from below to here to prevent !
1173 : ! energy and water conservation error. !
1174 : ! Also add corresponding tendency of cloud !
1175 : ! droplet number concentration. !
1176 : ! Note that since 'positive_moisture' will convert !
1177 : ! negative 'ql,qi,nl,ni' into positive, below may cause !
1178 : ! larger 'ql,qi' and smaller 'qv' than it should be if !
1179 : ! 'ql0,qi0' further below becomes negative. !
1180 : ! This feature is inevitable at the current stage. !
1181 : ! Note that in future, I can impose consistency between !
1182 : ! positive_moisture and positive_tracer for nl,ni. !
1183 : ! ----------------------------------------------------- !
1184 :
1185 : qlten(:iend,:mkx) = qlten(:iend,:mkx) - rqc_l(:iend,:mkx)
1186 : qiten(:iend,:mkx) = qiten(:iend,:mkx) - rqc_i(:iend,:mkx)
1187 : sten(:iend,:mkx) = sten(:iend,:mkx) - ( xls - xlv ) * rqc_i(:iend,:mkx)
1188 : trten(:iend,:mkx,ixnumliq) = trten(:iend,:mkx,ixnumliq) - rnc_l(:iend,:mkx)
1189 : trten(:iend,:mkx,ixnumice) = trten(:iend,:mkx,ixnumice) - rnc_i(:iend,:mkx)
1190 :
1191 : ! --------------------------------------------- !
1192 : ! Prevent negative cloud condensate and tracers !
1193 : ! --------------------------------------------- !
1194 :
1195 : qv0_c(:iend,:mkx) = qv0(:iend,:mkx) + qvten(:iend,:mkx)*dt
1196 : ql0_c(:iend,:mkx) = ql0(:iend,:mkx) + qlten(:iend,:mkx)*dt
1197 : qi0_c(:iend,:mkx) = qi0(:iend,:mkx) + qiten(:iend,:mkx)*dt
1198 : t0_c(:iend,:mkx) = t0(:iend,:mkx) + (1._r8/cp)*sten(:iend,:mkx)*dt
1199 : s0_c(:iend,:mkx) = s0(:iend,:mkx) + sten(:iend,:mkx)*dt
1200 : do mt = 1, ncnst
1201 : tr0_c(:iend,:mkx,mt) = tr0(:iend,:mkx,mt) + trten(:iend,:mkx,mt)*dt
1202 : enddo
1203 :
1204 : ! Note : Since 'positive_moisture' will only perform condensation not the evaporation,
1205 : ! we don't need to modify corresponding 'nl,ni'.
1206 : ! Thus, current version is completely OK.
1207 :
1208 : sten_ori(:iend,:mkx) = sten(:iend,:mkx)
1209 : qvten_ori(:iend,:mkx) = qvten(:iend,:mkx)
1210 : qlten_ori(:iend,:mkx) = qlten(:iend,:mkx)
1211 : qiten_ori(:iend,:mkx) = qiten(:iend,:mkx)
1212 : do mt = 1, ncnst
1213 : trten_ori(:iend,:mkx,mt) = trten(:iend,:mkx,mt)
1214 : enddo
1215 :
1216 : call positive_moisture( &
1217 : cp, xlv, xls, mix, iend, &
1218 : mkx, dt, qmin(1), qmin(ixcldliq), qmin(ixcldice), &
1219 : dp0, qv0_c, ql0_c, qi0_c, t0_c, &
1220 : s0_c, qvten, qlten, qiten, sten )
1221 :
1222 : do mt = 1, ncnst
1223 :
1224 : if( cnst_get_type_byind(mt) .eq. 'wet' ) then
1225 : pdel0(:iend,:mkx) = dp0(:iend,:mkx)
1226 : else
1227 : pdel0(:iend,:mkx) = dpdry0(:iend,:mkx)
1228 : endif
1229 :
1230 : if (cnst_is_mam_num(mt)) then
1231 : trmin = 1.e-5_r8
1232 : else
1233 : trmin = qmin(mt)
1234 : end if
1235 :
1236 : call positive_tracer( mix, iend, mkx, dt, trmin, pdel0, tr0_c(:,:,mt), trten(:,:,mt) )
1237 :
1238 : enddo
1239 :
1240 : ! ----------------------------------------------------- !
1241 : ! Treatment of reserved liquid-ice water for CAM !
1242 : ! All the reserved condensate is converted into liquid. !
1243 : ! ----------------------------------------------------- !
1244 :
1245 : ! Important : While below treatment looks different from what is being done in the uwshcu.F90,
1246 : ! below is actually the same as what is being done in the uwshcu.F90. The process
1247 : ! used below is (1) convert detrained ice into liquid, which involves melting cooling,
1248 : ! and then (2) subtract this total detrained liquid + ice from the grid-mean qlten.
1249 : ! Sep.08.2010. In the new scheme, it will be 'rqc_l = rqc_i = 0'. Thus, below block does nothing.
1250 : ! But it is no harm to keep below block.
1251 : ! Jan.04.2012. Below block will be used again for many reasons : (1) inctease TGCLDLWP, (2) reduce PREH20,
1252 : ! (2) reduce too strong SWCF over the far eastern equatorial Pacific.
1253 : ! Currently, below produces conservation error of both energy and water.
1254 : ! Jan.04.2012. Conservation errors of energy and moisture dissappears if I locate below block
1255 : ! before applying 'positive_moisture'. Thus, I relocated below block above.
1256 :
1257 : ! qlten(:iend,:mkx) = qlten(:iend,:mkx) + rqc_i(:iend,:mkx)
1258 : ! qiten(:iend,:mkx) = qiten(:iend,:mkx) - rqc_i(:iend,:mkx)
1259 : ! sten(:iend,:mkx) = sten(:iend,:mkx) - ( xls - xlv ) * rqc_i(:iend,:mkx)
1260 : ! qlten(:iend,:mkx) = qlten(:iend,:mkx) - rqc_l(:iend,:mkx) - rqc_i(:iend,:mkx)
1261 : ! trten(:iend,:mkx,ixnumliq) = trten(:iend,:mkx,ixnumliq) - rnc_l(:iend,:mkx)
1262 : ! trten(:iend,:mkx,ixnumice) = trten(:iend,:mkx,ixnumice) - rnc_i(:iend,:mkx)
1263 :
1264 : ! --------------------------- !
1265 : ! Reverse vertical coordinate !
1266 : ! --------------------------- !
1267 :
1268 : ! Reverse cloud top/base interface indices
1269 :
1270 : out%cnt(:iend) = real(mkx,r8) + 1._r8 - cnt(:iend)
1271 : out%cnb(:iend) = real(mkx,r8) + 1._r8 - cnb(:iend)
1272 :
1273 : ! Reverse variables defined at the interfaces
1274 :
1275 : do k = 0, mkx
1276 : k_inv = mkx + 1 - k
1277 : out%cmfmc(:iend,k_inv) = cmf_u(:iend,k) ! Updraft mass flux at top of layer
1278 : out%slflx(:iend,k_inv) = slflx(:iend,k) ! Net liquid static energy flux
1279 : out%qtflx(:iend,k_inv) = qtflx(:iend,k) ! Net total water flux
1280 : end do
1281 :
1282 : ! Reverse variables defined at the layer mid-point
1283 :
1284 : do k = 1, mkx
1285 : k_inv = mkx + 1 - k
1286 : ptend%q(:iend,k_inv,1) = qvten(:iend,k) ! Convective tendency of specific humidity
1287 : ptend%q(:iend,k_inv,ixcldliq) = qlten(:iend,k) ! Convective tendency of liquid water mixing ratio
1288 : ptend%q(:iend,k_inv,ixcldice) = qiten(:iend,k) ! Convective tendency of ice mixing ratio
1289 : ptend%s(:iend,k_inv) = sten(:iend,k) ! Convective tendency of static energy
1290 : ptend%u(:iend,k_inv) = uten(:iend,k) ! Convective tendency of zonal wind
1291 : ptend%v(:iend,k_inv) = vten(:iend,k) ! Convective tendency of meridional wind
1292 :
1293 : out%rqc(:iend,k_inv) = rqc(:iend,k) ! Convective tendency of reserved (suspended) liquid + ice water
1294 :
1295 : ! These quantities are being put into physics buffer fields that were meant for shallow convection.
1296 : ! This should be fixed.
1297 : shfrc(:iend,k_inv) = am_u(:iend,k)
1298 : icwmrsh(:iend,k_inv) = qlm_u(:iend,k) + qim_u(:iend,k)
1299 : rprdsh(:iend,k_inv) = qrten(:iend,k) + qsten(:iend,k)
1300 : evapc_inv(:iend,k_inv) = evapc(:iend,k) ! Evaporation rate of convective precipitation within environment
1301 :
1302 : do mt = 2, ncnst
1303 : if (mt /= ixcldliq .and. mt /= ixcldice) then
1304 : ptend%q(:iend,k_inv,mt) = trten(:iend,k,mt)
1305 : end if
1306 : enddo
1307 :
1308 : ! Additional diagnostic output associated with positive tracer
1309 :
1310 : sten_pos_inv(:iend,k_inv) = sten(:iend,k) - sten_ori(:iend,k)
1311 : qvten_pos_inv(:iend,k_inv) = qvten(:iend,k) - qvten_ori(:iend,k)
1312 : qlten_pos_inv(:iend,k_inv) = qlten(:iend,k) - qlten_ori(:iend,k)
1313 : qiten_pos_inv(:iend,k_inv) = qiten(:iend,k) - qiten_ori(:iend,k)
1314 : qtten_pos_inv(:iend,k_inv) = qvten_pos_inv(:iend,k_inv) + qlten_pos_inv(:iend,k_inv) &
1315 : + qiten_pos_inv(:iend,k_inv)
1316 : slten_pos_inv(:iend,k_inv) = sten_pos_inv(:iend,k_inv) - xlv * qlten_pos_inv(:iend,k_inv) &
1317 : - xls * qiten_pos_inv(:iend,k_inv)
1318 : uten_pos_inv(:iend,k_inv) = 0._r8
1319 : vten_pos_inv(:iend,k_inv) = 0._r8
1320 : do mt = 1, ncnst
1321 : trten_pos_inv(:iend,k_inv,mt) = trten(:iend,k,mt) - trten_ori(:iend,k,mt)
1322 : enddo
1323 :
1324 : ! Save detrainment variables to pbuf for use in the cloud macrophysics
1325 : cmfr_det(:iend,k_inv) = cmf_det(:iend,k)
1326 : qlr_det(:iend,k_inv) = ql_det(:iend,k)
1327 : qir_det(:iend,k_inv) = qi_det(:iend,k)
1328 :
1329 : end do
1330 :
1331 : call outfld('slten_pos_SP' , slten_pos_inv, mix, lchnk)
1332 : call outfld('qtten_pos_SP' , qtten_pos_inv, mix, lchnk)
1333 : call outfld('uten_pos_SP' , uten_pos_inv, mix, lchnk)
1334 : call outfld('vten_pos_SP' , vten_pos_inv, mix, lchnk)
1335 : call outfld('sten_pos_SP' , sten_pos_inv, mix, lchnk)
1336 : call outfld('qvten_pos_SP' , qvten_pos_inv, mix, lchnk)
1337 : call outfld('qlten_pos_SP' , qlten_pos_inv, mix, lchnk)
1338 : call outfld('qiten_pos_SP' , qiten_pos_inv, mix, lchnk)
1339 : call outfld('nlten_pos_SP' , trten_pos_inv(:,:,ixnumliq), mix, lchnk)
1340 : call outfld('niten_pos_SP' , trten_pos_inv(:,:,ixnumice), mix, lchnk)
1341 : call outfld('CMFR_DET' , cmfr_det, pcols, lchnk)
1342 : call outfld('QLR_DET' , qlr_det, pcols, lchnk)
1343 : call outfld('QIR_DET' , qir_det, pcols, lchnk)
1344 :
1345 : do m = 1, ncnst
1346 : if (cnst_is_mam_num(m) .or. cnst_is_mam_mmr(m)) then
1347 : varname = trim(cnst_name(m))//'_pos_SP'
1348 : call outfld(trim(varname), trten_pos_inv(:,:,m), mix, lchnk)
1349 : end if
1350 : end do
1351 :
1352 : #endif
1353 :
1354 0 : end subroutine unicon_cam_tend
1355 :
1356 0 : subroutine unicon_cam_org_diags(state, pbuf)
1357 :
1358 : ! ------------------------------------------------------------------------
1359 : ! Insert the organization-related heterogeneities computed inside the
1360 : ! UNICON into the tracer arrays here before performing advection.
1361 : ! This is necessary to prevent any modifications of organization-related
1362 : ! heterogeneities by non convection-advection process, such as
1363 : ! dry and wet deposition of aerosols, MAM, etc.
1364 : ! Again, note that only UNICON and advection schemes are allowed to
1365 : ! changes to organization at this stage, although we can include the
1366 : ! effects of other physical processes in future.
1367 : ! ------------------------------------------------------------------------
1368 :
1369 : ! Arguments
1370 : type(physics_state), intent(inout) :: state
1371 : type(physics_buffer_desc), pointer :: pbuf(:)
1372 :
1373 : ! Local variables
1374 : real(r8), pointer, dimension(: ) :: awk_PBL
1375 : real(r8), pointer, dimension(: ) :: delta_thl_PBL
1376 : real(r8), pointer, dimension(: ) :: delta_qt_PBL
1377 : real(r8), pointer, dimension(: ) :: delta_u_PBL
1378 : real(r8), pointer, dimension(: ) :: delta_v_PBL
1379 :
1380 : integer :: i, itim, ncol
1381 : ! ------------------------------------------------------------------------
1382 :
1383 : #ifdef USE_UNICON
1384 :
1385 : ncol = state%ncol
1386 :
1387 : itim = pbuf_old_tim_idx()
1388 : call pbuf_get_field(pbuf, awk_PBL_idx, awk_PBL, start=(/1,itim/), kount=(/pcols,1/))
1389 : call pbuf_get_field(pbuf, delta_thl_PBL_idx, delta_thl_PBL, start=(/1,itim/), kount=(/pcols,1/))
1390 : call pbuf_get_field(pbuf, delta_qt_PBL_idx, delta_qt_PBL, start=(/1,itim/), kount=(/pcols,1/))
1391 : call pbuf_get_field(pbuf, delta_u_PBL_idx, delta_u_PBL, start=(/1,itim/), kount=(/pcols,1/))
1392 : call pbuf_get_field(pbuf, delta_v_PBL_idx, delta_v_PBL, start=(/1,itim/), kount=(/pcols,1/))
1393 :
1394 : do i = 1, ncol
1395 :
1396 : state%q(i,:,awk_cnst_ind) = awk_PBL(i)
1397 :
1398 : ! Add a constant offset of '100._r8' to 'delta_xxx' variables to be
1399 : ! consistent with the reading sentence of unicon.F90 and so to prevent
1400 : ! negative tracers.
1401 : state%q(i,:,thl_cnst_ind) = max( 0._r8, delta_thl_PBL(i) + 100._r8 )
1402 : state%q(i,:,qt_cnst_ind) = max( 0._r8, delta_qt_PBL(i) + 100._r8 )
1403 : state%q(i,:,u_cnst_ind) = max( 0._r8, delta_u_PBL(i) + 100._r8 )
1404 : state%q(i,:,v_cnst_ind) = max( 0._r8, delta_v_PBL(i) + 100._r8 )
1405 : end do
1406 :
1407 : #endif
1408 :
1409 0 : end subroutine unicon_cam_org_diags
1410 :
1411 : !==================================================================================================
1412 0 : end module unicon_cam
|