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