Line data Source code
1 : module zm_conv_intr
2 : !---------------------------------------------------------------------------------
3 : ! Purpose:
4 : !
5 : ! CAM interface to the Zhang-McFarlane deep convection scheme
6 : !
7 : ! Author: D.B. Coleman
8 : ! January 2010 modified by J. Kay to add COSP simulator fields to physics buffer
9 : !---------------------------------------------------------------------------------
10 : use shr_kind_mod, only: r8=>shr_kind_r8
11 : use physconst, only: cpair, epsilo, gravit, latvap, tmelt, rair
12 : use ppgrid, only: pver, pcols, pverp, begchunk, endchunk
13 : use zm_conv_evap, only: zm_conv_evap_run
14 : use zm_convr, only: zm_convr_init, zm_convr_run
15 : use zm_conv_convtran, only: zm_conv_convtran_run
16 : use zm_conv_momtran, only: zm_conv_momtran_run
17 :
18 : use rad_constituents, only: rad_cnst_get_info, rad_cnst_get_mode_num, rad_cnst_get_aer_mmr, &
19 : rad_cnst_get_aer_props, rad_cnst_get_mode_props !, &
20 : use cam_abortutils, only: endrun
21 : use physconst, only: pi
22 : use spmd_utils, only: masterproc
23 : use perf_mod
24 : use cam_logfile, only: iulog
25 : use constituents, only: cnst_add
26 :
27 : implicit none
28 : private
29 : save
30 :
31 : ! Public methods
32 :
33 : public ::&
34 : zm_conv_register, &! register fields in physics buffer
35 : zm_conv_readnl, &! read namelist
36 : zm_conv_init, &! initialize donner_deep module
37 : zm_conv_tend, &! return tendencies
38 : zm_conv_tend_2 ! return tendencies
39 :
40 : public zmconv_ke, zmconv_ke_lnd, zmconv_org ! needed by convect_shallow
41 :
42 : integer ::& ! indices for fields in the physics buffer
43 : zm_mu_idx, &
44 : zm_eu_idx, &
45 : zm_du_idx, &
46 : zm_md_idx, &
47 : zm_ed_idx, &
48 : zm_dp_idx, &
49 : zm_dsubcld_idx, &
50 : zm_jt_idx, &
51 : zm_maxg_idx, &
52 : zm_ideep_idx, &
53 : dp_flxprc_idx, &
54 : dp_flxsnw_idx, &
55 : dp_cldliq_idx, &
56 : ixorg, &
57 : dp_cldice_idx, &
58 : dlfzm_idx, & ! detrained convective cloud water mixing ratio.
59 : difzm_idx, & ! detrained convective cloud ice mixing ratio.
60 : dnlfzm_idx, & ! detrained convective cloud water num concen.
61 : dnifzm_idx, & ! detrained convective cloud ice num concen.
62 : prec_dp_idx, &
63 : snow_dp_idx, &
64 : mconzm_idx ! convective mass flux
65 :
66 : real(r8), parameter :: unset_r8 = huge(1.0_r8)
67 : real(r8) :: zmconv_c0_lnd = unset_r8
68 : real(r8) :: zmconv_c0_ocn = unset_r8
69 : real(r8) :: zmconv_ke = unset_r8
70 : real(r8) :: zmconv_ke_lnd = unset_r8
71 : real(r8) :: zmconv_momcu = unset_r8
72 : real(r8) :: zmconv_momcd = unset_r8
73 : integer :: zmconv_num_cin ! Number of negative buoyancy regions that are allowed
74 : ! before the convection top and CAPE calculations are completed.
75 : logical :: zmconv_org ! Parameterization for sub-grid scale convective organization for the ZM deep
76 : ! convective scheme based on Mapes and Neale (2011)
77 : real(r8) :: zmconv_dmpdz = unset_r8 ! Parcel fractional mass entrainment rate
78 : real(r8) :: zmconv_tiedke_add = unset_r8 ! Convective parcel temperature perturbation
79 : real(r8) :: zmconv_capelmt = unset_r8 ! Triggering thereshold for ZM convection
80 : logical :: zmconv_parcel_pbl = .false. ! switch for parcel pbl calculation
81 : real(r8) :: zmconv_tau = unset_r8 ! Timescale for convection
82 :
83 :
84 : ! indices for fields in the physics buffer
85 : integer :: cld_idx = 0
86 : integer :: icwmrdp_idx = 0
87 : integer :: rprddp_idx = 0
88 : integer :: fracis_idx = 0
89 : integer :: nevapr_dpcu_idx = 0
90 : integer :: dgnum_idx = 0
91 :
92 : integer :: nmodes
93 : integer :: nbulk
94 :
95 : !=========================================================================================
96 : contains
97 : !=========================================================================================
98 :
99 60360 : subroutine zm_conv_register
100 :
101 : !----------------------------------------
102 : ! Purpose: register fields with the physics buffer
103 : !----------------------------------------
104 :
105 : use physics_buffer, only : pbuf_add_field, dtype_r8, dtype_i4
106 :
107 : implicit none
108 :
109 : integer idx
110 :
111 1536 : call pbuf_add_field('ZM_MU', 'physpkg', dtype_r8, (/pcols,pver/), zm_mu_idx)
112 1536 : call pbuf_add_field('ZM_EU', 'physpkg', dtype_r8, (/pcols,pver/), zm_eu_idx)
113 1536 : call pbuf_add_field('ZM_DU', 'physpkg', dtype_r8, (/pcols,pver/), zm_du_idx)
114 1536 : call pbuf_add_field('ZM_MD', 'physpkg', dtype_r8, (/pcols,pver/), zm_md_idx)
115 1536 : call pbuf_add_field('ZM_ED', 'physpkg', dtype_r8, (/pcols,pver/), zm_ed_idx)
116 :
117 : ! wg layer thickness in mbs (between upper/lower interface).
118 1536 : call pbuf_add_field('ZM_DP', 'physpkg', dtype_r8, (/pcols,pver/), zm_dp_idx)
119 :
120 : ! wg layer thickness in mbs between lcl and maxi.
121 1536 : call pbuf_add_field('ZM_DSUBCLD', 'physpkg', dtype_r8, (/pcols/), zm_dsubcld_idx)
122 :
123 : ! wg top level index of deep cumulus convection.
124 1536 : call pbuf_add_field('ZM_JT', 'physpkg', dtype_i4, (/pcols/), zm_jt_idx)
125 :
126 : ! wg gathered values of maxi.
127 1536 : call pbuf_add_field('ZM_MAXG', 'physpkg', dtype_i4, (/pcols/), zm_maxg_idx)
128 :
129 : ! map gathered points to chunk index
130 1536 : call pbuf_add_field('ZM_IDEEP', 'physpkg', dtype_i4, (/pcols/), zm_ideep_idx)
131 :
132 : ! Flux of precipitation from deep convection (kg/m2/s)
133 1536 : call pbuf_add_field('DP_FLXPRC','global',dtype_r8,(/pcols,pverp/),dp_flxprc_idx)
134 :
135 : ! Flux of snow from deep convection (kg/m2/s)
136 1536 : call pbuf_add_field('DP_FLXSNW','global',dtype_r8,(/pcols,pverp/),dp_flxsnw_idx)
137 :
138 : ! deep gbm cloud liquid water (kg/kg)
139 1536 : call pbuf_add_field('DP_CLDLIQ','global',dtype_r8,(/pcols,pver/), dp_cldliq_idx)
140 :
141 : ! deep gbm cloud liquid water (kg/kg)
142 1536 : call pbuf_add_field('DP_CLDICE','global',dtype_r8,(/pcols,pver/), dp_cldice_idx)
143 :
144 1536 : call pbuf_add_field('ICWMRDP', 'physpkg',dtype_r8,(/pcols,pver/),icwmrdp_idx)
145 1536 : call pbuf_add_field('RPRDDP', 'physpkg',dtype_r8,(/pcols,pver/),rprddp_idx)
146 1536 : call pbuf_add_field('NEVAPR_DPCU','physpkg',dtype_r8,(/pcols,pver/),nevapr_dpcu_idx)
147 1536 : call pbuf_add_field('PREC_DP', 'physpkg',dtype_r8,(/pcols/), prec_dp_idx)
148 1536 : call pbuf_add_field('SNOW_DP', 'physpkg',dtype_r8,(/pcols/), snow_dp_idx)
149 :
150 : ! detrained convective cloud water mixing ratio.
151 1536 : call pbuf_add_field('DLFZM', 'physpkg', dtype_r8, (/pcols,pver/), dlfzm_idx)
152 : ! detrained convective cloud ice mixing ratio.
153 1536 : call pbuf_add_field('DIFZM', 'physpkg', dtype_r8, (/pcols,pver/), difzm_idx)
154 : ! convective mass fluxes
155 1536 : call pbuf_add_field('CMFMC_DP', 'physpkg', dtype_r8, (/pcols,pverp/), mconzm_idx)
156 :
157 : !CACNOTE - Is zm_org really a constituent or was it just a handy structure to use for an allocatable which persists in the run?
158 1536 : if (zmconv_org) then
159 0 : call cnst_add('ZM_ORG',0._r8,0._r8,0._r8,ixorg,longname='organization parameter')
160 : endif
161 :
162 1536 : end subroutine zm_conv_register
163 :
164 : !=========================================================================================
165 :
166 1536 : subroutine zm_conv_readnl(nlfile)
167 :
168 1536 : use spmd_utils, only: mpicom, masterproc, masterprocid, mpi_real8, mpi_integer, mpi_logical
169 : use namelist_utils, only: find_group_name
170 : use units, only: getunit, freeunit
171 :
172 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
173 :
174 : ! Local variables
175 : integer :: unitn, ierr
176 : character(len=*), parameter :: subname = 'zm_conv_readnl'
177 :
178 : namelist /zmconv_nl/ zmconv_c0_lnd, zmconv_c0_ocn, zmconv_num_cin, &
179 : zmconv_ke, zmconv_ke_lnd, zmconv_org, &
180 : zmconv_momcu, zmconv_momcd, &
181 : zmconv_dmpdz, zmconv_tiedke_add, zmconv_capelmt, &
182 : zmconv_parcel_pbl, zmconv_tau
183 : !-----------------------------------------------------------------------------
184 :
185 1536 : if (masterproc) then
186 2 : unitn = getunit()
187 2 : open( unitn, file=trim(nlfile), status='old' )
188 2 : call find_group_name(unitn, 'zmconv_nl', status=ierr)
189 2 : if (ierr == 0) then
190 2 : read(unitn, zmconv_nl, iostat=ierr)
191 2 : if (ierr /= 0) then
192 0 : call endrun(subname // ':: ERROR reading namelist')
193 : end if
194 : end if
195 2 : close(unitn)
196 2 : call freeunit(unitn)
197 :
198 : end if
199 :
200 : ! Broadcast namelist variables
201 1536 : call mpi_bcast(zmconv_num_cin, 1, mpi_integer, masterprocid, mpicom, ierr)
202 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_num_cin")
203 1536 : call mpi_bcast(zmconv_c0_lnd, 1, mpi_real8, masterprocid, mpicom, ierr)
204 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_c0_lnd")
205 1536 : call mpi_bcast(zmconv_c0_ocn, 1, mpi_real8, masterprocid, mpicom, ierr)
206 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_c0_ocn")
207 1536 : call mpi_bcast(zmconv_ke, 1, mpi_real8, masterprocid, mpicom, ierr)
208 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_ke")
209 1536 : call mpi_bcast(zmconv_ke_lnd, 1, mpi_real8, masterprocid, mpicom, ierr)
210 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_ke_lnd")
211 1536 : call mpi_bcast(zmconv_momcu, 1, mpi_real8, masterprocid, mpicom, ierr)
212 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_momcu")
213 1536 : call mpi_bcast(zmconv_momcd, 1, mpi_real8, masterprocid, mpicom, ierr)
214 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_momcd")
215 1536 : call mpi_bcast(zmconv_org, 1, mpi_logical, masterprocid, mpicom, ierr)
216 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_org")
217 1536 : call mpi_bcast(zmconv_dmpdz, 1, mpi_real8, masterprocid, mpicom, ierr)
218 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_dmpdz")
219 1536 : call mpi_bcast(zmconv_tiedke_add, 1, mpi_real8, masterprocid, mpicom, ierr)
220 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_tiedke_add")
221 1536 : call mpi_bcast(zmconv_capelmt, 1, mpi_real8, masterprocid, mpicom, ierr)
222 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_capelmt")
223 1536 : call mpi_bcast(zmconv_parcel_pbl, 1, mpi_logical, masterprocid, mpicom, ierr)
224 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_parcel_pbl")
225 1536 : call mpi_bcast(zmconv_tau, 1, mpi_real8, masterprocid, mpicom, ierr)
226 1536 : if (ierr /= 0) call endrun("zm_conv_readnl: FATAL: mpi_bcast: zmconv_tau")
227 :
228 1536 : end subroutine zm_conv_readnl
229 :
230 : !=========================================================================================
231 :
232 1536 : subroutine zm_conv_init(pref_edge)
233 :
234 : !----------------------------------------
235 : ! Purpose: declare output fields, initialize variables needed by convection
236 : !----------------------------------------
237 :
238 : use cam_history, only: addfld, add_default, horiz_only
239 : use ppgrid, only: pcols, pver
240 : use zm_convr, only: zm_convr_init
241 : use pmgrid, only: plev,plevp
242 : use spmd_utils, only: masterproc
243 : use phys_control, only: phys_deepconv_pbl, phys_getopts, cam_physpkg_is
244 : use physics_buffer, only: pbuf_get_index
245 :
246 : implicit none
247 :
248 : real(r8),intent(in) :: pref_edge(plevp) ! reference pressures at interfaces
249 :
250 : ! local variables
251 : real(r8), parameter :: scale_height = 7000._r8 ! std atm scale height (m)
252 : real(r8), parameter :: dz_min = 100._r8 ! minimum thickness for using
253 : ! zmconv_parcel_pbl=.false.
254 : real(r8) :: dz_bot_layer ! thickness of bottom layer (m)
255 :
256 : character(len=512) :: errmsg
257 : integer :: errflg
258 :
259 : logical :: no_deep_pbl ! if true, no deep convection in PBL
260 : integer limcnv ! top interface level limit for convection
261 : integer k, istat
262 : logical :: history_budget ! output tendencies and state variables for CAM4
263 : ! temperature, water vapor, cloud ice and cloud
264 : ! liquid budgets.
265 : integer :: history_budget_histfile_num ! output history file number for budget fields
266 :
267 : !
268 : ! Register fields with the output buffer
269 : !
270 :
271 1536 : if (zmconv_org) then
272 0 : call addfld ('ZM_ORG ', (/ 'lev' /), 'A', '- ','Organization parameter')
273 0 : call addfld ('ZM_ORG2D ', (/ 'lev' /), 'A', '- ','Organization parameter 2D')
274 : endif
275 1536 : call addfld ('PRECZ', horiz_only, 'A', 'm/s','total precipitation from ZM convection')
276 3072 : call addfld ('ZMDT', (/ 'lev' /), 'A', 'K/s','T tendency - Zhang-McFarlane moist convection')
277 3072 : call addfld ('ZMDQ', (/ 'lev' /), 'A', 'kg/kg/s','Q tendency - Zhang-McFarlane moist convection')
278 3072 : call addfld ('ZMDICE', (/ 'lev' /), 'A', 'kg/kg/s','Cloud ice tendency - Zhang-McFarlane convection')
279 3072 : call addfld ('ZMDLIQ', (/ 'lev' /), 'A', 'kg/kg/s','Cloud liq tendency - Zhang-McFarlane convection')
280 3072 : call addfld ('EVAPTZM', (/ 'lev' /), 'A', 'K/s','T tendency - Evaporation/snow prod from Zhang convection')
281 3072 : call addfld ('FZSNTZM', (/ 'lev' /), 'A', 'K/s','T tendency - Rain to snow conversion from Zhang convection')
282 3072 : call addfld ('EVSNTZM', (/ 'lev' /), 'A', 'K/s','T tendency - Snow to rain prod from Zhang convection')
283 3072 : call addfld ('EVAPQZM', (/ 'lev' /), 'A', 'kg/kg/s','Q tendency - Evaporation from Zhang-McFarlane moist convection')
284 :
285 3072 : call addfld ('ZMFLXPRC', (/ 'ilev' /), 'A', 'kg/m2/s','Flux of precipitation from ZM convection' )
286 3072 : call addfld ('ZMFLXSNW', (/ 'ilev' /), 'A', 'kg/m2/s','Flux of snow from ZM convection' )
287 3072 : call addfld ('ZMNTPRPD', (/ 'lev' /) , 'A', 'kg/kg/s','Net precipitation production from ZM convection')
288 3072 : call addfld ('ZMNTSNPD', (/ 'lev' /) , 'A', 'kg/kg/s','Net snow production from ZM convection' )
289 3072 : call addfld ('ZMEIHEAT', (/ 'lev' /) , 'A', 'W/kg' ,'Heating by ice and evaporation in ZM convection')
290 :
291 3072 : call addfld ('CMFMC_DP', (/ 'ilev' /), 'A', 'kg/m2/s','Convection mass flux from ZM deep ')
292 1536 : call addfld ('PRECCDZM', horiz_only, 'A', 'm/s','Convective precipitation rate from ZM deep')
293 :
294 1536 : call addfld ('PCONVB', horiz_only , 'A', 'Pa' ,'convection base pressure')
295 1536 : call addfld ('PCONVT', horiz_only , 'A', 'Pa' ,'convection top pressure')
296 :
297 1536 : call addfld ('CAPE', horiz_only, 'A', 'J/kg', 'Convectively available potential energy')
298 1536 : call addfld ('FREQZM', horiz_only , 'A', 'fraction', 'Fractional occurance of ZM convection')
299 :
300 3072 : call addfld ('ZMMTT', (/ 'lev' /), 'A', 'K/s', 'T tendency - ZM convective momentum transport')
301 3072 : call addfld ('ZMMTU', (/ 'lev' /), 'A', 'm/s2', 'U tendency - ZM convective momentum transport')
302 3072 : call addfld ('ZMMTV', (/ 'lev' /), 'A', 'm/s2', 'V tendency - ZM convective momentum transport')
303 :
304 3072 : call addfld ('ZMMU', (/ 'lev' /), 'A', 'kg/m2/s', 'ZM convection updraft mass flux')
305 3072 : call addfld ('ZMMD', (/ 'lev' /), 'A', 'kg/m2/s', 'ZM convection downdraft mass flux')
306 :
307 3072 : call addfld ('ZMUPGU', (/ 'lev' /), 'A', 'm/s2', 'zonal force from ZM updraft pressure gradient term')
308 3072 : call addfld ('ZMUPGD', (/ 'lev' /), 'A', 'm/s2', 'zonal force from ZM downdraft pressure gradient term')
309 3072 : call addfld ('ZMVPGU', (/ 'lev' /), 'A', 'm/s2', 'meridional force from ZM updraft pressure gradient term')
310 3072 : call addfld ('ZMVPGD', (/ 'lev' /), 'A', 'm/s2', 'merdional force from ZM downdraft pressure gradient term')
311 :
312 3072 : call addfld ('ZMICUU', (/ 'lev' /), 'A', 'm/s', 'ZM in-cloud U updrafts')
313 3072 : call addfld ('ZMICUD', (/ 'lev' /), 'A', 'm/s', 'ZM in-cloud U downdrafts')
314 3072 : call addfld ('ZMICVU', (/ 'lev' /), 'A', 'm/s', 'ZM in-cloud V updrafts')
315 3072 : call addfld ('ZMICVD', (/ 'lev' /), 'A', 'm/s', 'ZM in-cloud V downdrafts')
316 :
317 3072 : call addfld ('DIFZM' ,(/ 'lev' /), 'A','kg/kg/s ','Detrained ice water from ZM convection')
318 3072 : call addfld ('DLFZM' ,(/ 'lev' /), 'A','kg/kg/s ','Detrained liquid water from ZM convection')
319 :
320 : call phys_getopts( history_budget_out = history_budget, &
321 1536 : history_budget_histfile_num_out = history_budget_histfile_num)
322 :
323 1536 : if (zmconv_org) then
324 0 : call add_default('ZM_ORG', 1, ' ')
325 0 : call add_default('ZM_ORG2D', 1, ' ')
326 : endif
327 1536 : if ( history_budget ) then
328 0 : call add_default('EVAPTZM ', history_budget_histfile_num, ' ')
329 0 : call add_default('EVAPQZM ', history_budget_histfile_num, ' ')
330 0 : call add_default('ZMDT ', history_budget_histfile_num, ' ')
331 0 : call add_default('ZMDQ ', history_budget_histfile_num, ' ')
332 0 : call add_default('ZMDLIQ ', history_budget_histfile_num, ' ')
333 0 : call add_default('ZMDICE ', history_budget_histfile_num, ' ')
334 0 : call add_default('ZMMTT ', history_budget_histfile_num, ' ')
335 : end if
336 :
337 : !
338 : ! Limit deep convection to regions below 40 mb
339 : ! Note this calculation is repeated in the shallow convection interface
340 : !
341 1536 : limcnv = 0 ! null value to check against below
342 1536 : if (pref_edge(1) >= 4.e3_r8) then
343 0 : limcnv = 1
344 : else
345 52224 : do k=1,plev
346 52224 : if (pref_edge(k) < 4.e3_r8 .and. pref_edge(k+1) >= 4.e3_r8) then
347 1536 : limcnv = k
348 1536 : exit
349 : end if
350 : end do
351 1536 : if ( limcnv == 0 ) limcnv = plevp
352 : end if
353 :
354 1536 : if (masterproc) then
355 2 : write(iulog,*)'ZM_CONV_INIT: Deep convection will be capped at intfc ',limcnv, &
356 4 : ' which is ',pref_edge(limcnv),' pascals'
357 : end if
358 :
359 : ! If thickness of bottom layer is less than dz_min, and zmconv_parcel_pbl=.false.,
360 : ! then issue a warning.
361 1536 : dz_bot_layer = scale_height * log(pref_edge(pverp)/pref_edge(pver))
362 1536 : if (dz_bot_layer < dz_min .and. .not. zmconv_parcel_pbl) then
363 0 : if (masterproc) then
364 0 : write(iulog,*)'********** WARNING **********'
365 0 : write(iulog,*)' ZM_CONV_INIT: Bottom layer thickness (m) is ', dz_bot_layer
366 0 : write(iulog,*)' The namelist variable zmconv_parcel_pbl should be set to .true.'
367 0 : write(iulog,*)' when the bottom layer thickness is < ', dz_min
368 0 : write(iulog,*)'********** WARNING **********'
369 : end if
370 : end if
371 :
372 1536 : no_deep_pbl = phys_deepconv_pbl()
373 : !CACNOTE - Need to check errflg and report errors
374 : call zm_convr_init(cpair, epsilo, gravit, latvap, tmelt, rair, &
375 : limcnv,zmconv_c0_lnd, zmconv_c0_ocn, zmconv_ke, zmconv_ke_lnd, &
376 : zmconv_momcu, zmconv_momcd, zmconv_num_cin, zmconv_org, &
377 : no_deep_pbl, zmconv_tiedke_add, &
378 : zmconv_capelmt, zmconv_dmpdz,zmconv_parcel_pbl, zmconv_tau, &
379 1536 : masterproc, iulog, errmsg, errflg)
380 :
381 1536 : cld_idx = pbuf_get_index('CLD')
382 1536 : fracis_idx = pbuf_get_index('FRACIS')
383 :
384 1536 : end subroutine zm_conv_init
385 : !=========================================================================================
386 : !subroutine zm_conv_tend(state, ptend, tdt)
387 :
388 2730672 : subroutine zm_conv_tend(pblh ,mcon ,cme , &
389 : tpert ,zdu , &
390 : rliq ,rice ,ztodt , &
391 : jctop ,jcbot , &
392 : state ,ptend_all ,landfrac, pbuf)
393 :
394 :
395 1536 : use cam_history, only: outfld
396 : use physics_types, only: physics_state, physics_ptend
397 : use physics_types, only: physics_ptend_init, physics_update
398 : use physics_types, only: physics_state_copy, physics_state_dealloc
399 : use physics_types, only: physics_ptend_sum, physics_ptend_dealloc
400 :
401 : use time_manager, only: get_nstep, is_first_step
402 : use physics_buffer, only : pbuf_get_field, physics_buffer_desc, pbuf_old_tim_idx
403 : use constituents, only: pcnst, cnst_get_ind, cnst_is_convtran1
404 : use check_energy, only: check_energy_chng
405 : use physconst, only: gravit, latice, latvap, tmelt, cpwv, cpliq, rh2o
406 :
407 : use phys_control, only: cam_physpkg_is
408 :
409 : ! Arguments
410 :
411 : type(physics_state), intent(in),target :: state ! Physics state variables
412 : type(physics_ptend), intent(out) :: ptend_all ! individual parameterization tendencies
413 : type(physics_buffer_desc), pointer :: pbuf(:)
414 :
415 : real(r8), intent(in) :: ztodt ! 2 delta t (model time increment)
416 : real(r8), intent(in) :: pblh(pcols) ! Planetary boundary layer height
417 : real(r8), intent(in) :: tpert(pcols) ! Thermal temperature excess
418 : real(r8), intent(in) :: landfrac(pcols) ! RBN - Landfrac
419 :
420 : real(r8), intent(out) :: mcon(pcols,pverp) ! Convective mass flux--m sub c
421 : real(r8), intent(out) :: cme(pcols,pver) ! cmf condensation - evaporation
422 : real(r8), intent(out) :: zdu(pcols,pver) ! detraining mass flux
423 :
424 : real(r8), intent(out) :: rliq(pcols) ! reserved liquid (not yet in cldliq) for energy integrals
425 : real(r8), intent(out) :: rice(pcols) ! reserved ice (not yet in cldice) for energy integrals
426 :
427 :
428 : ! Local variables
429 : character(len=512) :: errmsg
430 : integer :: errflg
431 :
432 : integer :: i,k,l,m
433 : integer :: ilon ! global longitude index of a column
434 : integer :: ilat ! global latitude index of a column
435 : integer :: nstep
436 : integer :: ixcldice, ixcldliq ! constituent indices for cloud liquid and ice water.
437 : integer :: lchnk ! chunk identifier
438 : integer :: ncol ! number of atmospheric columns
439 : integer :: itim_old ! for physics buffer fields
440 :
441 : real(r8) :: ftem(pcols,pver) ! Temporary workspace for outfld variables
442 : real(r8) :: ntprprd(pcols,pver) ! evap outfld: net precip production in layer
443 : real(r8) :: ntsnprd(pcols,pver) ! evap outfld: net snow production in layer
444 : real(r8) :: tend_s_snwprd (pcols,pver) ! Heating rate of snow production
445 : real(r8) :: tend_s_snwevmlt(pcols,pver) ! Heating rate of evap/melting of snow
446 : real(r8) :: fake_dpdry(pcols,pver) ! used in convtran call
447 :
448 : ! physics types
449 65016 : type(physics_state) :: state1 ! locally modify for evaporation to use, not returned
450 2730672 : type(physics_ptend),target :: ptend_loc ! package tendencies
451 :
452 : ! physics buffer fields
453 65016 : real(r8), pointer, dimension(:) :: prec ! total precipitation
454 65016 : real(r8), pointer, dimension(:) :: snow ! snow from ZM convection
455 65016 : real(r8), pointer, dimension(:,:) :: cld
456 65016 : real(r8), pointer, dimension(:,:) :: ql ! wg grid slice of cloud liquid water.
457 65016 : real(r8), pointer, dimension(:,:) :: rprd ! rain production rate
458 65016 : real(r8), pointer, dimension(:,:,:) :: fracis ! fraction of transported species that are insoluble
459 65016 : real(r8), pointer, dimension(:,:) :: evapcdp ! Evaporation of deep convective precipitation
460 65016 : real(r8), pointer, dimension(:,:) :: flxprec ! Convective-scale flux of precip at interfaces (kg/m2/s)
461 65016 : real(r8), pointer, dimension(:,:) :: flxsnow ! Convective-scale flux of snow at interfaces (kg/m2/s)
462 65016 : real(r8), pointer, dimension(:,:) :: dp_cldliq
463 65016 : real(r8), pointer, dimension(:,:) :: dp_cldice
464 65016 : real(r8), pointer :: dlf(:,:) ! detrained convective cloud water mixing ratio.
465 65016 : real(r8), pointer :: dif(:,:) ! detrained convective cloud ice mixing ratio.
466 65016 : real(r8), pointer :: dnlf(:,:) ! detrained convective cloud water num concen.
467 65016 : real(r8), pointer :: dnif(:,:) ! detrained convective cloud ice num concen.
468 : real(r8), pointer :: lambdadpcu(:,:) ! slope of cloud liquid size distr
469 : real(r8), pointer :: mudpcu(:,:) ! width parameter of droplet size distr
470 65016 : real(r8), pointer :: mconzm(:,:) !convective mass fluxes
471 :
472 65016 : real(r8), pointer :: mu(:,:) ! (pcols,pver)
473 65016 : real(r8), pointer :: eu(:,:) ! (pcols,pver)
474 65016 : real(r8), pointer :: du(:,:) ! (pcols,pver)
475 65016 : real(r8), pointer :: md(:,:) ! (pcols,pver)
476 65016 : real(r8), pointer :: ed(:,:) ! (pcols,pver)
477 65016 : real(r8), pointer :: dp(:,:) ! (pcols,pver)
478 65016 : real(r8), pointer :: dsubcld(:) ! (pcols)
479 65016 : integer, pointer :: jt(:) ! (pcols)
480 65016 : integer, pointer :: maxg(:) ! (pcols)
481 65016 : integer, pointer :: ideep(:) ! (pcols)
482 : integer :: lengath
483 :
484 : real(r8) :: jctop(pcols) ! o row of top-of-deep-convection indices passed out.
485 : real(r8) :: jcbot(pcols) ! o row of base of cloud indices passed out.
486 :
487 : real(r8) :: pcont(pcols), pconb(pcols), freqzm(pcols)
488 :
489 : ! history output fields
490 : real(r8) :: cape(pcols) ! w convective available potential energy.
491 : real(r8) :: mu_out(pcols,pver)
492 : real(r8) :: md_out(pcols,pver)
493 :
494 : ! used in momentum transport calculation
495 : real(r8) :: pguallu(pcols, pver)
496 : real(r8) :: pguallv(pcols, pver)
497 : real(r8) :: pgdallu(pcols, pver)
498 : real(r8) :: pgdallv(pcols, pver)
499 : real(r8) :: icwuu(pcols,pver)
500 : real(r8) :: icwuv(pcols,pver)
501 : real(r8) :: icwdu(pcols,pver)
502 : real(r8) :: icwdv(pcols,pver)
503 : real(r8) :: seten(pcols, pver)
504 : logical :: l_windt(2)
505 : real(r8) :: tfinal1, tfinal2
506 : integer :: ii
507 :
508 65016 : real(r8),pointer :: zm_org2d(:,:)
509 65016 : real(r8),allocatable :: orgt_alloc(:,:), org_alloc(:,:)
510 :
511 130032 : real(r8) :: zm_org2d_ncol(state%ncol,pver)
512 195048 : real(r8) :: orgt_ncol(state%ncol,pver), org_ncol(state%ncol,pver)
513 :
514 : logical :: lq(pcnst)
515 :
516 : !----------------------------------------------------------------------
517 :
518 : ! initialize
519 65016 : lchnk = state%lchnk
520 65016 : ncol = state%ncol
521 65016 : nstep = get_nstep()
522 :
523 65016 : ftem = 0._r8
524 65016 : mu_out(:,:) = 0._r8
525 65016 : md_out(:,:) = 0._r8
526 :
527 65016 : call physics_state_copy(state,state1) ! copy state to local state1.
528 :
529 65016 : lq(:) = .FALSE.
530 65016 : lq(1) = .TRUE.
531 65016 : if (zmconv_org) then
532 0 : lq(ixorg) = .TRUE.
533 : endif
534 65016 : call physics_ptend_init(ptend_loc, state%psetcols, 'zm_convr_run', ls=.true., lq=lq)! initialize local ptend type
535 :
536 : !
537 : ! Associate pointers with physics buffer fields
538 : !
539 65016 : itim_old = pbuf_old_tim_idx()
540 260064 : call pbuf_get_field(pbuf, cld_idx, cld, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
541 :
542 65016 : call pbuf_get_field(pbuf, icwmrdp_idx, ql )
543 65016 : call pbuf_get_field(pbuf, rprddp_idx, rprd )
544 65016 : call pbuf_get_field(pbuf, fracis_idx, fracis, start=(/1,1,1/), kount=(/pcols, pver, pcnst/) )
545 65016 : call pbuf_get_field(pbuf, nevapr_dpcu_idx, evapcdp )
546 65016 : call pbuf_get_field(pbuf, prec_dp_idx, prec )
547 65016 : call pbuf_get_field(pbuf, snow_dp_idx, snow )
548 :
549 65016 : call pbuf_get_field(pbuf, zm_mu_idx, mu)
550 65016 : call pbuf_get_field(pbuf, zm_eu_idx, eu)
551 65016 : call pbuf_get_field(pbuf, zm_du_idx, du)
552 65016 : call pbuf_get_field(pbuf, zm_md_idx, md)
553 65016 : call pbuf_get_field(pbuf, zm_ed_idx, ed)
554 65016 : call pbuf_get_field(pbuf, zm_dp_idx, dp)
555 65016 : call pbuf_get_field(pbuf, zm_dsubcld_idx, dsubcld)
556 65016 : call pbuf_get_field(pbuf, zm_jt_idx, jt)
557 65016 : call pbuf_get_field(pbuf, zm_maxg_idx, maxg)
558 65016 : call pbuf_get_field(pbuf, zm_ideep_idx, ideep)
559 :
560 65016 : call pbuf_get_field(pbuf, dlfzm_idx, dlf)
561 65016 : call pbuf_get_field(pbuf, difzm_idx, dif)
562 65016 : call pbuf_get_field(pbuf, mconzm_idx, mconzm)
563 :
564 65016 : allocate(dnlf(pcols,pver), dnif(pcols,pver))
565 :
566 : !
567 : ! Begin with Zhang-McFarlane (1996) convection parameterization
568 : !
569 65016 : call t_startf ('zm_convr_run')
570 :
571 65016 : if (zmconv_org) then
572 0 : allocate(zm_org2d(pcols,pver))
573 0 : allocate(org_alloc(ncol,pver))
574 0 : allocate(orgt_alloc(ncol,pver))
575 0 : org_ncol(:ncol,:) = state%q(1:ncol,:,ixorg)
576 : endif
577 :
578 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
579 102855312 : ptend_loc%q(:,:,1) = 0._r8
580 102855312 : ptend_loc%s(:,:) = 0._r8
581 65016 : mcon(:,:) = 0._r8
582 102855312 : dlf(:,:) = 0._r8
583 65016 : cme(:,:) = 0._r8
584 65016 : cape(:) = 0._r8
585 65016 : zdu(:,:) = 0._r8
586 102855312 : rprd(:,:) = 0._r8
587 102855312 : dif(:,:) = 0._r8
588 102855312 : dnlf(:,:) = 0._r8
589 102855312 : dnif(:,:) = 0._r8
590 102855312 : mu(:,:) = 0._r8
591 102855312 : eu(:,:) = 0._r8
592 102855312 : du(:,:) = 0._r8
593 102855312 : md(:,:) = 0._r8
594 102855312 : ed(:,:) = 0._r8
595 102855312 : dp(:,:) = 0._r8
596 1105272 : dsubcld(:) = 0._r8
597 65016 : jctop(:) = 0._r8
598 65016 : jcbot(:) = 0._r8
599 1105272 : prec(:) = 0._r8
600 65016 : rliq(:) = 0._r8
601 65016 : rice(:) = 0._r8
602 1105272 : ideep(:) = 0._r8
603 : !REMOVECAM_END
604 :
605 : !CACNOTE - Need to check errflg and report errors
606 : call zm_convr_run(ncol, pver, &
607 : pverp, gravit, latice, cpwv, cpliq, rh2o, &
608 0 : state%t(:ncol,:), state%q(:ncol,:,1), prec(:ncol), &
609 0 : pblh(:ncol), state%zm(:ncol,:), state%phis(:ncol), state%zi(:ncol,:), ptend_loc%q(:ncol,:,1), &
610 0 : ptend_loc%s(:ncol,:), state%pmid(:ncol,:), state%pint(:ncol,:), state%pdel(:ncol,:), &
611 : .5_r8*ztodt, mcon(:ncol,:), cme(:ncol,:), cape(:ncol), &
612 0 : tpert(:ncol), dlf(:ncol,:), zdu(:ncol,:), rprd(:ncol,:), &
613 0 : mu(:ncol,:), md(:ncol,:), du(:ncol,:), eu(:ncol,:), ed(:ncol,:), &
614 0 : dp(:ncol,:), dsubcld(:ncol), jt(:ncol), maxg(:ncol), ideep(:ncol), &
615 0 : ql(:ncol,:), rliq(:ncol), landfrac(:ncol), &
616 130032 : org_ncol(:ncol,:), orgt_ncol(:ncol,:), zm_org2d_ncol(:ncol,:), &
617 0 : dif(:ncol,:), dnlf(:ncol,:), dnif(:ncol,:), &
618 195048 : rice(:ncol), errmsg, errflg)
619 :
620 :
621 65016 : if (zmconv_org) then
622 0 : ptend_loc%q(:,:,ixorg)=orgt_ncol(:ncol,:)
623 0 : zm_org2d(:ncol,:) = zm_org2d_ncol(:ncol,:)
624 : endif
625 :
626 1105272 : lengath = count(ideep > 0)
627 65016 : if (lengath > ncol) lengath = ncol ! should not happen, but force it to not be larger than ncol for safety sake
628 :
629 1105272 : jctop(:) = real(pver,r8)
630 1105272 : jcbot(:) = 1._r8
631 247119 : do i = 1,lengath
632 182103 : jctop(ideep(i)) = real(jt(i), r8)
633 247119 : jcbot(ideep(i)) = real(maxg(i), r8)
634 : end do
635 :
636 65016 : call outfld('CAPE', cape, pcols, lchnk) ! RBN - CAPE output
637 : !
638 : ! Output fractional occurance of ZM convection
639 : !
640 65016 : freqzm(:) = 0._r8
641 247119 : do i = 1,lengath
642 247119 : freqzm(ideep(i)) = 1.0_r8
643 : end do
644 65016 : call outfld('FREQZM ',freqzm ,pcols ,lchnk )
645 : !
646 : ! Convert mass flux from reported mb/s to kg/m^2/s
647 : !
648 102112920 : mcon(:ncol,:pverp) = mcon(:ncol,:pverp) * 100._r8/gravit
649 102112920 : mconzm(:ncol,:pverp) = mcon(:ncol,:pverp)
650 :
651 65016 : call outfld('CMFMC_DP', mconzm, pcols, lchnk)
652 :
653 : ! Store upward and downward mass fluxes in un-gathered arrays
654 : ! + convert from mb/s to kg/m^2/s
655 247119 : do i=1,lengath
656 17182698 : do k=1,pver
657 16935579 : ii = ideep(i)
658 16935579 : mu_out(ii,k) = mu(i,k) * 100._r8/gravit
659 17117682 : md_out(ii,k) = md(i,k) * 100._r8/gravit
660 : end do
661 : end do
662 :
663 65016 : call outfld('ZMMU', mu_out, pcols, lchnk)
664 65016 : call outfld('ZMMD', md_out, pcols, lchnk)
665 :
666 101027304 : ftem(:ncol,:pver) = ptend_loc%s(:ncol,:pver)/cpair
667 65016 : call outfld('ZMDT ',ftem ,pcols ,lchnk )
668 65016 : call outfld('ZMDQ ',ptend_loc%q(1,1,1) ,pcols ,lchnk )
669 65016 : call t_stopf ('zm_convr_run')
670 :
671 65016 : call outfld('DIFZM' ,dif ,pcols, lchnk)
672 65016 : call outfld('DLFZM' ,dlf ,pcols, lchnk)
673 :
674 1085616 : pcont(:ncol) = state%ps(:ncol)
675 1085616 : pconb(:ncol) = state%ps(:ncol)
676 247119 : do i = 1,lengath
677 247119 : if (maxg(i).gt.jt(i)) then
678 182103 : pcont(ideep(i)) = state%pmid(ideep(i),jt(i)) ! gathered array (or jctop ungathered)
679 182103 : pconb(ideep(i)) = state%pmid(ideep(i),maxg(i))! gathered array
680 : endif
681 : ! write(iulog,*) ' pcont, pconb ', pcont(i), pconb(i), cnt(i), cnb(i)
682 : end do
683 65016 : call outfld('PCONVT ',pcont ,pcols ,lchnk )
684 65016 : call outfld('PCONVB ',pconb ,pcols ,lchnk )
685 :
686 65016 : call physics_ptend_init(ptend_all, state%psetcols, 'zm_conv_tend')
687 :
688 : ! add tendency from this process to tendencies from other processes
689 65016 : call physics_ptend_sum(ptend_loc,ptend_all, ncol)
690 :
691 : ! update physics state type state1 with ptend_loc
692 65016 : call physics_update(state1, ptend_loc, ztodt)
693 :
694 : ! initialize ptend for next process
695 65016 : lq(:) = .FALSE.
696 65016 : lq(1) = .TRUE.
697 65016 : if (zmconv_org) then
698 0 : lq(ixorg) = .TRUE.
699 : endif
700 65016 : call physics_ptend_init(ptend_loc, state1%psetcols, 'zm_conv_evap_run', ls=.true., lq=lq)
701 :
702 65016 : call t_startf ('zm_conv_evap_run')
703 : !
704 : ! Determine the phase of the precipitation produced and add latent heat of fusion
705 : ! Evaporate some of the precip directly into the environment (Sundqvist)
706 : ! Allow this to use the updated state1 and the fresh ptend_loc type
707 : ! heating and specific humidity tendencies produced
708 : !
709 :
710 65016 : call pbuf_get_field(pbuf, dp_flxprc_idx, flxprec )
711 65016 : call pbuf_get_field(pbuf, dp_flxsnw_idx, flxsnow )
712 65016 : call pbuf_get_field(pbuf, dp_cldliq_idx, dp_cldliq )
713 65016 : call pbuf_get_field(pbuf, dp_cldice_idx, dp_cldice )
714 101027304 : dp_cldliq(:ncol,:) = 0._r8
715 101027304 : dp_cldice(:ncol,:) = 0._r8
716 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
717 103960584 : flxprec(:,:) = 0._r8
718 103960584 : flxsnow(:,:) = 0._r8
719 1105272 : snow(:) = 0._r8
720 : !REMOVECAM_END
721 :
722 : call zm_conv_evap_run(state1%ncol, pver, pverp, &
723 : gravit, latice, latvap, tmelt, &
724 : cpair, zmconv_ke, zmconv_ke_lnd, zmconv_org, &
725 0 : state1%t(:ncol,:),state1%pmid(:ncol,:),state1%pdel(:ncol,:),state1%q(:ncol,:pver,1), &
726 0 : landfrac(:ncol), &
727 0 : ptend_loc%s(:ncol,:), tend_s_snwprd(:ncol,:), tend_s_snwevmlt(:ncol,:), ptend_loc%q(:ncol,:pver,1), &
728 0 : rprd(:ncol,:), cld(:ncol,:), ztodt, &
729 65016 : prec(:ncol), snow(:ncol), ntprprd(:ncol,:), ntsnprd(:ncol,:), flxprec(:ncol,:), flxsnow(:ncol,:) )
730 :
731 101027304 : evapcdp(:ncol,:pver) = ptend_loc%q(:ncol,:pver,1)
732 :
733 65016 : if (zmconv_org) then
734 0 : ptend_loc%q(:ncol,:pver,ixorg) = min(1._r8,max(0._r8,(50._r8*1000._r8*1000._r8*abs(evapcdp(:ncol,:pver))) &
735 0 : -(state%q(:ncol,:pver,ixorg)/10800._r8)))
736 0 : ptend_loc%q(:ncol,:pver,ixorg) = (ptend_loc%q(:ncol,:pver,ixorg) - state%q(:ncol,:pver,ixorg))/ztodt
737 : endif
738 :
739 : !
740 : ! Write out variables from zm_conv_evap_run
741 : !
742 101027304 : ftem(:ncol,:pver) = ptend_loc%s(:ncol,:pver)/cpair
743 65016 : call outfld('EVAPTZM ',ftem ,pcols ,lchnk )
744 101027304 : ftem(:ncol,:pver) = tend_s_snwprd (:ncol,:pver)/cpair
745 65016 : call outfld('FZSNTZM ',ftem ,pcols ,lchnk )
746 101027304 : ftem(:ncol,:pver) = tend_s_snwevmlt(:ncol,:pver)/cpair
747 65016 : call outfld('EVSNTZM ',ftem ,pcols ,lchnk )
748 65016 : call outfld('EVAPQZM ',ptend_loc%q(1,1,1) ,pcols ,lchnk )
749 65016 : call outfld('ZMFLXPRC', flxprec, pcols, lchnk)
750 65016 : call outfld('ZMFLXSNW', flxsnow, pcols, lchnk)
751 65016 : call outfld('ZMNTPRPD', ntprprd, pcols, lchnk)
752 65016 : call outfld('ZMNTSNPD', ntsnprd, pcols, lchnk)
753 65016 : call outfld('ZMEIHEAT', ptend_loc%s, pcols, lchnk)
754 65016 : call outfld('CMFMC_DP ',mcon , pcols ,lchnk )
755 65016 : call outfld('PRECCDZM ',prec, pcols ,lchnk )
756 :
757 :
758 65016 : call t_stopf ('zm_conv_evap_run')
759 :
760 65016 : call outfld('PRECZ ', prec , pcols, lchnk)
761 :
762 : ! add tendency from this process to tend from other processes here
763 65016 : call physics_ptend_sum(ptend_loc,ptend_all, ncol)
764 :
765 : ! update physics state type state1 with ptend_loc
766 65016 : call physics_update(state1, ptend_loc, ztodt)
767 :
768 :
769 : ! Momentum Transport (non-cam3 physics)
770 :
771 65016 : if ( .not. cam_physpkg_is('cam3')) then
772 :
773 65016 : call physics_ptend_init(ptend_loc, state1%psetcols, 'zm_conv_momtran_run', ls=.true., lu=.true., lv=.true.)
774 :
775 65016 : l_windt(1) = .true.
776 65016 : l_windt(2) = .true.
777 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
778 102855312 : ptend_loc%s(:,:) = 0._r8
779 102855312 : ptend_loc%u(:,:) = 0._r8
780 102855312 : ptend_loc%v(:,:) = 0._r8
781 : !REMOVECAM_END
782 :
783 65016 : call t_startf ('zm_conv_momtran_run')
784 :
785 : call zm_conv_momtran_run (ncol, pver, pverp, &
786 0 : l_windt,state1%u(:ncol,:), state1%v(:ncol,:), 2, mu(:ncol,:), md(:ncol,:), &
787 : zmconv_momcu, zmconv_momcd, &
788 0 : du(:ncol,:), eu(:ncol,:), ed(:ncol,:), dp(:ncol,:), dsubcld(:ncol), &
789 0 : jt(:ncol), maxg(:ncol), ideep(:ncol), 1, lengath, &
790 0 : nstep, ptend_loc%u(:ncol,:), ptend_loc%v(:ncol,:),&
791 0 : pguallu(:ncol,:), pguallv(:ncol,:), pgdallu(:ncol,:), pgdallv(:ncol,:), &
792 65016 : icwuu(:ncol,:), icwuv(:ncol,:), icwdu(:ncol,:), icwdv(:ncol,:), ztodt, seten(:ncol,:) )
793 65016 : call t_stopf ('zm_conv_momtran_run')
794 :
795 101027304 : ptend_loc%s(:ncol,:pver) = seten(:ncol,:pver)
796 :
797 65016 : call physics_ptend_sum(ptend_loc,ptend_all, ncol)
798 :
799 : ! Output ptend variables before they are set to zero with physics_update
800 65016 : call outfld('ZMMTU', ptend_loc%u, pcols, lchnk)
801 65016 : call outfld('ZMMTV', ptend_loc%v, pcols, lchnk)
802 :
803 : ! update physics state type state1 with ptend_loc
804 65016 : call physics_update(state1, ptend_loc, ztodt)
805 :
806 101027304 : ftem(:ncol,:pver) = seten(:ncol,:pver)/cpair
807 65016 : if (zmconv_org) then
808 0 : call outfld('ZM_ORG', state%q(:,:,ixorg), pcols, lchnk)
809 0 : call outfld('ZM_ORG2D', zm_org2d, pcols, lchnk)
810 : endif
811 65016 : call outfld('ZMMTT', ftem , pcols, lchnk)
812 :
813 : ! Output apparent force from pressure gradient
814 65016 : call outfld('ZMUPGU', pguallu, pcols, lchnk)
815 65016 : call outfld('ZMUPGD', pgdallu, pcols, lchnk)
816 65016 : call outfld('ZMVPGU', pguallv, pcols, lchnk)
817 65016 : call outfld('ZMVPGD', pgdallv, pcols, lchnk)
818 :
819 : ! Output in-cloud winds
820 65016 : call outfld('ZMICUU', icwuu, pcols, lchnk)
821 65016 : call outfld('ZMICUD', icwdu, pcols, lchnk)
822 65016 : call outfld('ZMICVU', icwuv, pcols, lchnk)
823 65016 : call outfld('ZMICVD', icwdv, pcols, lchnk)
824 :
825 : end if
826 :
827 : ! Transport cloud water and ice only
828 65016 : call cnst_get_ind('CLDLIQ', ixcldliq)
829 65016 : call cnst_get_ind('CLDICE', ixcldice)
830 :
831 65016 : lq(:) = .FALSE.
832 2665656 : lq(2:) = cnst_is_convtran1(2:)
833 65016 : call physics_ptend_init(ptend_loc, state1%psetcols, 'convtran1', lq=lq)
834 :
835 :
836 : ! dpdry is not used in this call to convtran since the cloud liquid and ice mixing
837 : ! ratios are moist
838 65016 : fake_dpdry(:,:) = 0._r8
839 :
840 65016 : call t_startf ('convtran1')
841 :
842 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
843 4217132808 : ptend_loc%q(:,:,:) = 0._r8
844 : !REMOVECAM_END
845 :
846 : call zm_conv_convtran_run (ncol, pver, &
847 0 : ptend_loc%lq,state1%q(:ncol,:,:), pcnst, mu(:ncol,:), md(:ncol,:), &
848 0 : du(:ncol,:), eu(:ncol,:), ed(:ncol,:), dp(:ncol,:), dsubcld(:ncol), &
849 0 : jt(:ncol), maxg(:ncol), ideep(:ncol), 1, lengath, &
850 65016 : nstep, fracis(:ncol,:,:), ptend_loc%q(:ncol,:,:), fake_dpdry(:ncol,:), ztodt)
851 65016 : call t_stopf ('convtran1')
852 :
853 65016 : call outfld('ZMDICE ',ptend_loc%q(1,1,ixcldice) ,pcols ,lchnk )
854 65016 : call outfld('ZMDLIQ ',ptend_loc%q(1,1,ixcldliq) ,pcols ,lchnk )
855 :
856 : ! add tendency from this process to tend from other processes here
857 65016 : call physics_ptend_sum(ptend_loc,ptend_all, ncol)
858 :
859 65016 : call physics_state_dealloc(state1)
860 65016 : call physics_ptend_dealloc(ptend_loc)
861 :
862 65016 : if (zmconv_org) then
863 0 : deallocate(zm_org2d)
864 : end if
865 :
866 65016 : deallocate(dnlf, dnif)
867 :
868 260064 : end subroutine zm_conv_tend
869 : !=========================================================================================
870 :
871 :
872 2529432 : subroutine zm_conv_tend_2( state, ptend, ztodt, pbuf)
873 :
874 65016 : use physics_types, only: physics_state, physics_ptend, physics_ptend_init
875 : use time_manager, only: get_nstep
876 : use physics_buffer, only: pbuf_get_index, pbuf_get_field, physics_buffer_desc
877 : use constituents, only: pcnst, cnst_is_convtran2
878 :
879 : ! Arguments
880 : type(physics_state), intent(in ) :: state ! Physics state variables
881 : type(physics_ptend), intent(out) :: ptend ! indivdual parameterization tendencies
882 :
883 : type(physics_buffer_desc), pointer :: pbuf(:)
884 :
885 : real(r8), intent(in) :: ztodt ! 2 delta t (model time increment)
886 :
887 : ! Local variables
888 : integer :: i, lchnk, istat
889 : integer :: lengath ! number of columns with deep convection
890 : integer :: nstep
891 : integer :: ncol
892 :
893 : real(r8), dimension(pcols,pver) :: dpdry
894 :
895 : ! physics buffer fields
896 58824 : real(r8), pointer :: fracis(:,:,:) ! fraction of transported species that are insoluble
897 58824 : real(r8), pointer :: mu(:,:) ! (pcols,pver)
898 58824 : real(r8), pointer :: eu(:,:) ! (pcols,pver)
899 58824 : real(r8), pointer :: du(:,:) ! (pcols,pver)
900 58824 : real(r8), pointer :: md(:,:) ! (pcols,pver)
901 58824 : real(r8), pointer :: ed(:,:) ! (pcols,pver)
902 58824 : real(r8), pointer :: dp(:,:) ! (pcols,pver)
903 58824 : real(r8), pointer :: dsubcld(:) ! (pcols)
904 58824 : integer, pointer :: jt(:) ! (pcols)
905 58824 : integer, pointer :: maxg(:) ! (pcols)
906 58824 : integer, pointer :: ideep(:) ! (pcols)
907 : !-----------------------------------------------------------------------------------
908 :
909 :
910 58824 : call physics_ptend_init(ptend, state%psetcols, 'convtran2', lq=cnst_is_convtran2 )
911 :
912 58824 : call pbuf_get_field(pbuf, fracis_idx, fracis)
913 58824 : call pbuf_get_field(pbuf, zm_mu_idx, mu)
914 58824 : call pbuf_get_field(pbuf, zm_eu_idx, eu)
915 58824 : call pbuf_get_field(pbuf, zm_du_idx, du)
916 58824 : call pbuf_get_field(pbuf, zm_md_idx, md)
917 58824 : call pbuf_get_field(pbuf, zm_ed_idx, ed)
918 58824 : call pbuf_get_field(pbuf, zm_dp_idx, dp)
919 58824 : call pbuf_get_field(pbuf, zm_dsubcld_idx, dsubcld)
920 58824 : call pbuf_get_field(pbuf, zm_jt_idx, jt)
921 58824 : call pbuf_get_field(pbuf, zm_maxg_idx, maxg)
922 58824 : call pbuf_get_field(pbuf, zm_ideep_idx, ideep)
923 :
924 :
925 58824 : lchnk = state%lchnk
926 58824 : ncol = state%ncol
927 58824 : nstep = get_nstep()
928 :
929 1000008 : lengath = count(ideep > 0)
930 58824 : if (lengath > ncol) lengath = ncol ! should not happen, but force it to not be larger than ncol for safety sake
931 :
932 823536 : if (any(ptend%lq(:))) then
933 : ! initialize dpdry for call to convtran
934 : ! it is used for tracers of dry mixing ratio type
935 58824 : dpdry = 0._r8
936 223698 : do i = 1, lengath
937 15556980 : dpdry(i,:) = state%pdeldry(ideep(i),:)/100._r8
938 : end do
939 :
940 58824 : call t_startf ('convtran2')
941 :
942 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
943 3815501112 : ptend%q(:,:,:) = 0._r8
944 : !REMOVECAM_END
945 :
946 : call zm_conv_convtran_run (ncol, pver, &
947 0 : ptend%lq,state%q(:ncol,:,:), pcnst, mu(:ncol,:), md(:ncol,:), &
948 0 : du(:ncol,:), eu(:ncol,:), ed(:ncol,:), dp(:ncol,:), dsubcld(:ncol), &
949 0 : jt(:ncol), maxg(:ncol), ideep(:ncol), 1, lengath, &
950 58824 : nstep, fracis(:ncol,:,:), ptend%q(:ncol,:,:), dpdry(:ncol,:), ztodt)
951 58824 : call t_stopf ('convtran2')
952 : end if
953 :
954 117648 : end subroutine zm_conv_tend_2
955 :
956 : !=========================================================================================
957 :
958 :
959 : end module zm_conv_intr
|