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