Line data Source code
1 : module rk_stratiform_cam
2 :
3 : !-------------------------------------------------------------------------------------------------------
4 : !
5 : ! Provides the CAM interface to the Rasch and Kristjansson (RK)
6 : ! prognostic cloud microphysics, and the cam4 macrophysics.
7 : !
8 : !-------------------------------------------------------------------------------------------------------
9 :
10 : use shr_kind_mod, only: r8=>shr_kind_r8
11 : use ppgrid, only: pcols, pver, pverp
12 : use physconst, only: gravit, latvap, latice
13 : use phys_control, only: phys_getopts
14 : use constituents, only: pcnst
15 : use spmd_utils, only: masterproc
16 :
17 : use cam_logfile, only: iulog
18 : use cam_abortutils, only: endrun
19 : use perf_mod
20 :
21 : implicit none
22 : private
23 : save
24 :
25 : public :: rk_stratiform_cam_register, rk_stratiform_cam_init_cnst, rk_stratiform_cam_implements_cnst
26 : public :: rk_stratiform_cam_init
27 : public :: rk_stratiform_cam_tend
28 : public :: rk_stratiform_cam_readnl
29 :
30 : ! Physics buffer indices
31 : integer :: landm_idx = 0
32 :
33 : integer :: qcwat_idx = 0
34 : integer :: lcwat_idx = 0
35 : integer :: tcwat_idx = 0
36 :
37 : integer :: cld_idx = 0
38 : integer :: ast_idx = 0
39 : integer :: concld_idx = 0
40 : integer :: fice_idx = 0
41 :
42 : integer :: qme_idx = 0
43 : integer :: prain_idx = 0
44 : integer :: nevapr_idx = 0
45 :
46 : integer :: wsedl_idx = 0
47 :
48 : integer :: rei_idx = 0
49 : integer :: rel_idx = 0
50 :
51 : integer :: shfrc_idx = 0
52 : integer :: cmfmc_sh_idx = 0
53 :
54 : integer :: prec_str_idx = 0
55 : integer :: snow_str_idx = 0
56 : integer :: prec_sed_idx = 0
57 : integer :: snow_sed_idx = 0
58 : integer :: prec_pcw_idx = 0
59 : integer :: snow_pcw_idx = 0
60 :
61 : integer :: ls_flxprc_idx = 0
62 : integer :: ls_flxsnw_idx = 0
63 :
64 : ! Physics buffer indices for convective_cloud_cover
65 : integer :: sh_frac_idx = 0
66 : integer :: dp_frac_idx = 0
67 :
68 : integer, parameter :: ncnst = 2 ! Number of constituents
69 : character(len=8), dimension(ncnst), parameter & ! Constituent names
70 : :: cnst_names = (/'CLDLIQ', 'CLDICE'/)
71 : logical :: use_shfrc ! Local copy of flag from convect_shallow_use_shfrc
72 :
73 : logical :: do_cnst = .false. ! True when this module has registered constituents.
74 :
75 : integer :: &
76 : ixcldliq, &! cloud liquid amount index
77 : ixcldice ! cloud ice amount index
78 :
79 : real(r8), parameter :: unset_r8 = huge(1.0_r8)
80 : logical :: do_psrhmin
81 :
82 : !===============================================================================
83 : contains
84 : !===============================================================================
85 1024 : subroutine rk_stratiform_cam_readnl(nlfile)
86 :
87 : use namelist_utils, only: find_group_name
88 : use units, only: getunit, freeunit
89 :
90 : use prognostic_cloud_water, only: prognostic_cloud_water_init
91 : use physconst, only: tmelt, rhodair, pi
92 :
93 : use mpishorthand
94 :
95 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
96 :
97 : ! Local variables
98 : integer :: unitn, ierr
99 : character(len=*), parameter :: subname = 'rk_stratiform_readnl'
100 :
101 : character(len=512) :: errmsg
102 : integer :: errflg
103 :
104 : ! Namelist variables
105 : real(r8) :: rk_strat_icritw = unset_r8 ! icritw = threshold for autoconversion of warm ice
106 : real(r8) :: rk_strat_icritc = unset_r8 ! icritc = threshold for autoconversion of cold ice
107 : real(r8) :: rk_strat_conke = unset_r8 ! conke = tunable constant for evaporation of precip
108 : real(r8) :: rk_strat_r3lcrit = unset_r8 ! r3lcrit = critical radius where liq conversion begins
109 : real(r8) :: rk_strat_polstrat_rhmin = unset_r8 ! condensation threadhold in polar stratosphere
110 :
111 : namelist /rk_stratiform_nl/ rk_strat_icritw, rk_strat_icritc, rk_strat_conke, rk_strat_r3lcrit
112 : namelist /rk_stratiform_nl/ rk_strat_polstrat_rhmin
113 :
114 : !-----------------------------------------------------------------------------
115 :
116 1026 : if (masterproc) then
117 2 : unitn = getunit()
118 2 : open( unitn, file=trim(nlfile), status='old' )
119 2 : call find_group_name(unitn, 'rk_stratiform_nl', status=ierr)
120 2 : if (ierr == 0) then
121 2 : read(unitn, rk_stratiform_nl, iostat=ierr)
122 2 : if (ierr /= 0) then
123 0 : call endrun(subname // ':: ERROR reading namelist')
124 : end if
125 : end if
126 2 : close(unitn)
127 2 : call freeunit(unitn)
128 : end if
129 :
130 : #ifdef SPMD
131 : ! Broadcast namelist variables
132 1024 : call mpibcast(rk_strat_icritw, 1, mpir8, 0, mpicom)
133 1024 : call mpibcast(rk_strat_icritc, 1, mpir8, 0, mpicom)
134 1024 : call mpibcast(rk_strat_conke, 1, mpir8, 0, mpicom)
135 1024 : call mpibcast(rk_strat_r3lcrit, 1, mpir8, 0, mpicom)
136 1024 : call mpibcast(rk_strat_polstrat_rhmin, 1, mpir8, 0, mpicom)
137 : #endif
138 :
139 1024 : do_psrhmin = rk_strat_polstrat_rhmin .ne. unset_r8
140 :
141 : call prognostic_cloud_water_init(&
142 : amIRoot = masterproc, &
143 : iulog = iulog, &
144 : tmelt = tmelt, &
145 : rhodair = rhodair, &
146 : pi = pi, &
147 : icritc_in = rk_strat_icritc, &
148 : icritw_in = rk_strat_icritw, &
149 : conke_in = rk_strat_conke, &
150 : r3lcrit_in = rk_strat_r3lcrit, &
151 : do_psrhmin_in = do_psrhmin, &
152 : psrhmin_in = rk_strat_polstrat_rhmin, &
153 : errmsg = errmsg, &
154 1024 : errflg = errflg)
155 :
156 1024 : if(errflg /= 0) then
157 0 : call endrun(subname // errmsg)
158 : endif
159 :
160 : ! cloud_particle_sedimentation is initialized in pkg_cld_sediment,
161 : ! which has its own namelist read.
162 :
163 1024 : end subroutine rk_stratiform_cam_readnl
164 :
165 1024 : subroutine rk_stratiform_cam_register
166 :
167 : !---------------------------------------------------------------------- !
168 : ! !
169 : ! Register the constituents (cloud liquid and cloud ice) and the fields !
170 : ! in the physics buffer. !
171 : ! !
172 : !---------------------------------------------------------------------- !
173 :
174 : use constituents, only: cnst_add, pcnst
175 : use physconst, only: mwh2o, cpair
176 :
177 : use physics_buffer, only : pbuf_add_field, dtype_r8, dyn_time_lvls
178 :
179 : !-----------------------------------------------------------------------
180 :
181 : ! Take note of the fact that we are registering constituents.
182 1024 : do_cnst = .true.
183 :
184 : ! Register cloud water and save indices.
185 : call cnst_add(cnst_names(1), mwh2o, cpair, 0._r8, ixcldliq, &
186 1024 : longname='Grid box averaged cloud liquid amount', is_convtran1=.true.)
187 : call cnst_add(cnst_names(2), mwh2o, cpair, 0._r8, ixcldice, &
188 1024 : longname='Grid box averaged cloud ice amount', is_convtran1=.true.)
189 :
190 4096 : call pbuf_add_field('QCWAT', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), qcwat_idx)
191 4096 : call pbuf_add_field('LCWAT', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), lcwat_idx)
192 4096 : call pbuf_add_field('TCWAT', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), tcwat_idx)
193 :
194 4096 : call pbuf_add_field('CLD', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), cld_idx) ! cloud_area_fraction
195 4096 : call pbuf_add_field('AST', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), ast_idx) ! stratiform_cloud_area_fraction
196 4096 : call pbuf_add_field('CONCLD', 'global', dtype_r8, (/pcols,pver,dyn_time_lvls/), concld_idx) ! convective_cloud_area_fraction
197 :
198 1024 : call pbuf_add_field('FICE', 'physpkg', dtype_r8, (/pcols,pver/), fice_idx) ! mass_fraction_of_ice_content_within_stratiform_cloud
199 :
200 1024 : call pbuf_add_field('QME', 'physpkg', dtype_r8, (/pcols,pver/), qme_idx) ! net_condensation_rate_due_to_microphysics
201 1024 : call pbuf_add_field('PRAIN', 'physpkg', dtype_r8, (/pcols,pver/), prain_idx) ! precipitation_production_due_to_microphysics
202 1024 : call pbuf_add_field('NEVAPR', 'physpkg', dtype_r8, (/pcols,pver/), nevapr_idx) ! precipitation_evaporation_due_to_microphysics
203 :
204 1024 : call pbuf_add_field('WSEDL', 'physpkg', dtype_r8, (/pcols,pver/), wsedl_idx) ! sedimentation velocity of liquid stratus cloud droplet [m s-1]
205 :
206 1024 : call pbuf_add_field('REI', 'physpkg', dtype_r8, (/pcols,pver/), rei_idx) ! effective_radius_of_stratiform_cloud_ice_particle
207 1024 : call pbuf_add_field('REL', 'physpkg', dtype_r8, (/pcols,pver/), rel_idx) ! effective_radius_of_stratiform_cloud_liquid_water_particle
208 :
209 1024 : call pbuf_add_field('LS_FLXPRC', 'physpkg', dtype_r8, (/pcols,pverp/), ls_flxprc_idx) ! effective_radius_of_stratiform_cloud_liquid_water_particle
210 1024 : call pbuf_add_field('LS_FLXSNW', 'physpkg', dtype_r8, (/pcols,pverp/), ls_flxsnw_idx) ! stratiform_snow_flux_at_interface
211 :
212 1024 : end subroutine rk_stratiform_cam_register
213 :
214 : !===============================================================================
215 :
216 10800 : function rk_stratiform_cam_implements_cnst(name)
217 :
218 : !----------------------------------------------------------------------------- !
219 : ! !
220 : ! Return true if specified constituent is implemented by this package !
221 : ! !
222 : !----------------------------------------------------------------------------- !
223 :
224 : character(len=*), intent(in) :: name ! constituent name
225 : logical :: rk_stratiform_cam_implements_cnst ! return value
226 :
227 : !-----------------------------------------------------------------------
228 :
229 16200 : rk_stratiform_cam_implements_cnst = (do_cnst .and. any(name == cnst_names))
230 :
231 22624 : end function rk_stratiform_cam_implements_cnst
232 :
233 : !===============================================================================
234 :
235 10800 : subroutine rk_stratiform_cam_init_cnst(name, latvals, lonvals, mask, q)
236 :
237 : !----------------------------------------------------------------------- !
238 : ! !
239 : ! Initialize the cloud water mixing ratios (liquid and ice), if they are !
240 : ! not read from the initial file !
241 : ! !
242 : !----------------------------------------------------------------------- !
243 :
244 : character(len=*), intent(in) :: name ! constituent name
245 : real(r8), intent(in) :: latvals(:) ! lat in degrees (ncol)
246 : real(r8), intent(in) :: lonvals(:) ! lon in degrees (ncol)
247 : logical, intent(in) :: mask(:) ! Only initialize where .true.
248 : real(r8), intent(out) :: q(:,:) ! kg tracer/kg dry air (gcol, plev
249 : !-----------------------------------------------------------------------
250 : integer :: k
251 :
252 16200 : if (any(name == cnst_names)) then
253 291600 : do k = 1, size(q, 2)
254 5065200 : where(mask)
255 561600 : q(:, k) = 0.0_r8
256 : end where
257 : end do
258 : end if
259 :
260 10800 : end subroutine rk_stratiform_cam_init_cnst
261 :
262 : !===============================================================================
263 :
264 1024 : subroutine rk_stratiform_cam_init()
265 :
266 : !-------------------------------------------- !
267 : ! !
268 : ! Initialize the cloud water parameterization !
269 : ! !
270 : !-------------------------------------------- !
271 :
272 : use physics_buffer, only: physics_buffer_desc, pbuf_get_index
273 : use constituents, only: cnst_get_ind, cnst_name, cnst_longname, sflxnam, apcnst, bpcnst
274 : use cam_history, only: addfld, add_default, horiz_only
275 : use convect_shallow, only: convect_shallow_use_shfrc
276 : use phys_control, only: cam_physpkg_is
277 : use physconst, only: tmelt, rhodair, rh2o
278 :
279 : integer :: m, mm
280 : logical :: history_amwg ! output the variables used by the AMWG diag package
281 : logical :: history_aerosol ! Output the MAM aerosol tendencies
282 : logical :: history_budget ! Output tendencies and state variables for CAM4
283 : ! temperature, water vapor, cloud ice and cloud
284 : ! liquid budgets.
285 : integer :: history_budget_histfile_num ! output history file number for budget fields
286 : !-----------------------------------------------------------------------
287 :
288 : call phys_getopts( history_aerosol_out = history_aerosol , &
289 : history_amwg_out = history_amwg , &
290 : history_budget_out = history_budget , &
291 1024 : history_budget_histfile_num_out = history_budget_histfile_num)
292 :
293 1024 : landm_idx = pbuf_get_index("LANDM")
294 :
295 : ! Find out whether shfrc from convect_shallow will be used in cldfrc
296 1024 : if( convect_shallow_use_shfrc() ) then
297 0 : use_shfrc = .true.
298 0 : shfrc_idx = pbuf_get_index('shfrc')
299 : else
300 1024 : use_shfrc = .false.
301 : endif
302 :
303 : ! Register history variables
304 :
305 3072 : do m = 1, ncnst
306 2048 : call cnst_get_ind( cnst_names(m), mm )
307 4096 : call addfld( cnst_name(mm), (/ 'lev' /), 'A', 'kg/kg', cnst_longname(mm) )
308 2048 : call addfld( sflxnam (mm), horiz_only, 'A', 'kg/m2/s', trim(cnst_name(mm))//' surface flux' )
309 5120 : if (history_amwg) then
310 2048 : call add_default( cnst_name(mm), 1, ' ' )
311 2048 : call add_default( sflxnam (mm), 1, ' ' )
312 : endif
313 : enddo
314 :
315 2048 : call addfld (apcnst(ixcldliq), (/ 'lev' /), 'A', 'kg/kg', trim(cnst_name(ixcldliq))//' after physics' )
316 2048 : call addfld (apcnst(ixcldice), (/ 'lev' /), 'A', 'kg/kg', trim(cnst_name(ixcldice))//' after physics' )
317 2048 : call addfld (bpcnst(ixcldliq), (/ 'lev' /), 'A', 'kg/kg', trim(cnst_name(ixcldliq))//' before physics' )
318 2048 : call addfld (bpcnst(ixcldice), (/ 'lev' /), 'A', 'kg/kg', trim(cnst_name(ixcldice))//' before physics' )
319 :
320 1024 : if( history_budget) then
321 0 : call add_default (cnst_name(ixcldliq), history_budget_histfile_num, ' ')
322 0 : call add_default (cnst_name(ixcldice), history_budget_histfile_num, ' ')
323 0 : call add_default (apcnst (ixcldliq), history_budget_histfile_num, ' ')
324 0 : call add_default (apcnst (ixcldice), history_budget_histfile_num, ' ')
325 0 : call add_default (bpcnst (ixcldliq), history_budget_histfile_num, ' ')
326 0 : call add_default (bpcnst (ixcldice), history_budget_histfile_num, ' ')
327 : end if
328 :
329 2048 : call addfld ('FWAUT', (/ 'lev' /), 'A', 'fraction', 'Relative importance of liquid autoconversion' )
330 2048 : call addfld ('FSAUT', (/ 'lev' /), 'A', 'fraction', 'Relative importance of ice autoconversion' )
331 2048 : call addfld ('FRACW', (/ 'lev' /), 'A', 'fraction', 'Relative importance of rain accreting liquid' )
332 2048 : call addfld ('FSACW', (/ 'lev' /), 'A', 'fraction', 'Relative importance of snow accreting liquid' )
333 2048 : call addfld ('FSACI', (/ 'lev' /), 'A', 'fraction', 'Relative importance of snow accreting ice' )
334 2048 : call addfld ('CME', (/ 'lev' /), 'A', 'kg/kg/s' , 'Rate of cond-evap within the cloud' )
335 2048 : call addfld ('CMEICE', (/ 'lev' /), 'A', 'kg/kg/s' , 'Rate of cond-evap of ice within the cloud' )
336 2048 : call addfld ('CMELIQ', (/ 'lev' /), 'A', 'kg/kg/s' , 'Rate of cond-evap of liq within the cloud' )
337 2048 : call addfld ('ICE2PR', (/ 'lev' /), 'A', 'kg/kg/s' , 'Rate of conversion of ice to precip' )
338 2048 : call addfld ('LIQ2PR', (/ 'lev' /), 'A', 'kg/kg/s' , 'Rate of conversion of liq to precip' )
339 2048 : call addfld ('ZMDLF', (/ 'lev' /), 'A', 'kg/kg/s' , 'Detrained liquid water from ZM convection' )
340 2048 : call addfld ('SHDLF', (/ 'lev' /), 'A', 'kg/kg/s' , 'Detrained liquid water from shallow convection' )
341 :
342 2048 : call addfld ('PRODPREC', (/ 'lev' /), 'A', 'kg/kg/s' , 'Rate of conversion of condensate to precip' )
343 2048 : call addfld ('EVAPPREC', (/ 'lev' /), 'A', 'kg/kg/s' , 'Rate of evaporation of falling precip' )
344 2048 : call addfld ('EVAPSNOW', (/ 'lev' /), 'A', 'kg/kg/s' , 'Rate of evaporation of falling snow' )
345 2048 : call addfld ('HPROGCLD', (/ 'lev' /), 'A', 'W/kg' , 'Heating from prognostic clouds' )
346 2048 : call addfld ('HCME', (/ 'lev' /), 'A', 'W/kg' , 'Heating from cond-evap within the cloud' )
347 2048 : call addfld ('HEVAP', (/ 'lev' /), 'A', 'W/kg' , 'Heating from evaporation of falling precip' )
348 2048 : call addfld ('HFREEZ', (/ 'lev' /), 'A', 'W/kg' , 'Heating rate due to freezing of precip' )
349 2048 : call addfld ('HMELT', (/ 'lev' /), 'A', 'W/kg' , 'Heating from snow melt' )
350 2048 : call addfld ('HREPART', (/ 'lev' /), 'A', 'W/kg' , 'Heating from cloud ice/liquid repartitioning' )
351 2048 : call addfld ('REPARTICE', (/ 'lev' /), 'A', 'kg/kg/s' , 'Cloud ice tendency from cloud ice/liquid repartitioning' )
352 2048 : call addfld ('REPARTLIQ', (/ 'lev' /), 'A', 'kg/kg/s' , 'Cloud liq tendency from cloud ice/liquid repartitioning' )
353 2048 : call addfld ('FICE', (/ 'lev' /), 'A', 'fraction', 'Fractional ice content within cloud' )
354 2048 : call addfld ('ICWMR', (/ 'lev' /), 'A', 'kg/kg' , 'Prognostic in-cloud water mixing ratio' )
355 2048 : call addfld ('ICIMR', (/ 'lev' /), 'A', 'kg/kg' , 'Prognostic in-cloud ice mixing ratio' )
356 1024 : call addfld ('PCSNOW', horiz_only , 'A', 'm/s' , 'Snow fall from prognostic clouds' )
357 :
358 2048 : call addfld ('DQSED', (/ 'lev' /), 'A', 'kg/kg/s' , 'Water vapor tendency from cloud sedimentation' )
359 2048 : call addfld ('DLSED', (/ 'lev' /), 'A', 'kg/kg/s' , 'Cloud liquid tendency from sedimentation' )
360 2048 : call addfld ('DISED', (/ 'lev' /), 'A', 'kg/kg/s' , 'Cloud ice tendency from sedimentation' )
361 2048 : call addfld ('HSED', (/ 'lev' /), 'A', 'W/kg' , 'Heating from cloud sediment evaporation' )
362 1024 : call addfld ('SNOWSED', horiz_only, 'A', 'm/s' , 'Snow from cloud ice sedimentation' )
363 1024 : call addfld ('RAINSED', horiz_only, 'A', 'm/s' , 'Rain from cloud liquid sedimentation' )
364 1024 : call addfld ('PRECSED', horiz_only, 'A', 'm/s' , 'Precipitation from cloud sedimentation' )
365 :
366 2048 : call addfld ('CLDST', (/ 'lev' /), 'A', 'fraction', 'Stratus cloud fraction' )
367 2048 : call addfld ('CONCLD', (/ 'lev' /), 'A', 'fraction', 'Convective cloud cover' )
368 :
369 2048 : call addfld ('AST', (/ 'lev' /), 'A', 'fraction', 'Stratus cloud fraction' )
370 2048 : call addfld ('LIQCLDF', (/ 'lev' /), 'A', 'fraction', 'Stratus Liquid cloud fraction' )
371 2048 : call addfld ('ICECLDF', (/ 'lev' /), 'A', 'fraction', 'Stratus ICE cloud fraction' )
372 2048 : call addfld ('IWC', (/ 'lev' /), 'A', 'kg/m3' , 'Grid box average ice water content' )
373 2048 : call addfld ('LWC', (/ 'lev' /), 'A', 'kg/m3' , 'Grid box average liquid water content' )
374 2048 : call addfld ('REL', (/ 'lev' /), 'A', 'micron' , 'effective liquid drop radius' )
375 2048 : call addfld ('REI', (/ 'lev' /), 'A', 'micron' , 'effective ice particle radius' )
376 :
377 1024 : if ( history_budget ) then
378 :
379 0 : call add_default ('EVAPSNOW ', history_budget_histfile_num, ' ')
380 0 : call add_default ('EVAPPREC ', history_budget_histfile_num, ' ')
381 0 : call add_default ('CMELIQ ', history_budget_histfile_num, ' ')
382 :
383 0 : if( cam_physpkg_is('cam4') ) then
384 :
385 0 : call add_default ('ZMDLF ', history_budget_histfile_num, ' ')
386 0 : call add_default ('CME ', history_budget_histfile_num, ' ')
387 0 : call add_default ('DQSED ', history_budget_histfile_num, ' ')
388 0 : call add_default ('DISED ', history_budget_histfile_num, ' ')
389 0 : call add_default ('DLSED ', history_budget_histfile_num, ' ')
390 0 : call add_default ('HSED ', history_budget_histfile_num, ' ')
391 0 : call add_default ('CMEICE ', history_budget_histfile_num, ' ')
392 0 : call add_default ('LIQ2PR ', history_budget_histfile_num, ' ')
393 0 : call add_default ('ICE2PR ', history_budget_histfile_num, ' ')
394 0 : call add_default ('HCME ', history_budget_histfile_num, ' ')
395 0 : call add_default ('HEVAP ', history_budget_histfile_num, ' ')
396 0 : call add_default ('HFREEZ ', history_budget_histfile_num, ' ')
397 0 : call add_default ('HMELT ', history_budget_histfile_num, ' ')
398 0 : call add_default ('HREPART ', history_budget_histfile_num, ' ')
399 0 : call add_default ('HPROGCLD ', history_budget_histfile_num, ' ')
400 0 : call add_default ('REPARTLIQ', history_budget_histfile_num, ' ')
401 0 : call add_default ('REPARTICE', history_budget_histfile_num, ' ')
402 :
403 : end if
404 :
405 : end if
406 :
407 1024 : if (history_amwg) then
408 1024 : call add_default ('ICWMR', 1, ' ')
409 1024 : call add_default ('ICIMR', 1, ' ')
410 1024 : call add_default ('CONCLD ', 1, ' ')
411 1024 : call add_default ('FICE ', 1, ' ')
412 : endif
413 :
414 : ! History Variables for COSP/CFMIP
415 2048 : call addfld ('LS_FLXPRC', (/ 'ilev' /), 'A', 'kg/m2/s', 'ls stratiform gbm interface rain+snow flux')
416 2048 : call addfld ('LS_FLXSNW', (/ 'ilev' /), 'A', 'kg/m2/s', 'ls stratiform gbm interface snow flux')
417 2048 : call addfld ('PRACWO', (/ 'lev' /), 'A', '1/s', 'Accretion of cloud water by rain')
418 2048 : call addfld ('PSACWO', (/ 'lev' /), 'A', '1/s', 'Accretion of cloud water by snow')
419 2048 : call addfld ('PSACIO', (/ 'lev' /), 'A', '1/s', 'Accretion of cloud ice by snow')
420 :
421 2048 : call addfld ('CLDLIQSTR', (/ 'lev' /), 'A', 'kg/kg', 'Stratiform CLDLIQ')
422 2048 : call addfld ('CLDICESTR', (/ 'lev' /), 'A', 'kg/kg', 'Stratiform CLDICE')
423 2048 : call addfld ('CLDLIQCON', (/ 'lev' /), 'A', 'kg/kg', 'Convective CLDLIQ')
424 2048 : call addfld ('CLDICECON', (/ 'lev' /), 'A', 'kg/kg', 'Convective CLDICE')
425 :
426 1024 : cmfmc_sh_idx = pbuf_get_index('CMFMC_SH') ! atmosphere_convective_mass_flux_due_to_shallow_convection
427 1024 : prec_str_idx = pbuf_get_index('PREC_STR') ! stratiform_rain_and_snow_surface_flux_due_to_sedimentation
428 1024 : snow_str_idx = pbuf_get_index('SNOW_STR') ! lwe_snow_and_cloud_ice_precipitation_rate_at_surface_due_to_microphysics
429 1024 : prec_pcw_idx = pbuf_get_index('PREC_PCW') ! lwe_stratiform_precipitation_rate_at_surface
430 1024 : snow_pcw_idx = pbuf_get_index('SNOW_PCW') ! lwe_snow_precipitation_rate_at_surface_due_to_microphysics
431 1024 : prec_sed_idx = pbuf_get_index('PREC_SED') ! stratiform_cloud_water_surface_flux_due_to_sedimentation
432 1024 : snow_sed_idx = pbuf_get_index('SNOW_SED') ! lwe_cloud_ice_sedimentation_rate_at_surface_due_to_microphysics
433 :
434 1024 : sh_frac_idx = pbuf_get_index('SH_FRAC')
435 1024 : dp_frac_idx = pbuf_get_index('DP_FRAC')
436 :
437 1024 : end subroutine rk_stratiform_cam_init
438 :
439 : !===============================================================================
440 :
441 281568 : subroutine rk_stratiform_cam_tend( &
442 : state, ptend_all, pbuf, dtime, icefrac, &
443 : landfrac, ocnfrac, snowh, dlf, &
444 : dlf2, rliq, cmfmc, ts, &
445 : sst, zdu)
446 :
447 : !-------------------------------------------------------- !
448 : ! !
449 : ! Interface to sedimentation, detrain, cloud fraction and !
450 : ! cloud macro - microphysics subroutines !
451 : ! !
452 : !-------------------------------------------------------- !
453 :
454 1024 : use cloud_fraction_fice, only: cloud_fraction_fice_run
455 : use physics_types, only: physics_state, physics_ptend
456 : use physics_types, only: physics_ptend_init, physics_update
457 : use physics_types, only: physics_ptend_sum, physics_state_copy
458 : use physics_types, only: physics_state_dealloc
459 : use cam_history, only: outfld
460 : use physics_buffer, only: physics_buffer_desc, pbuf_get_field, pbuf_old_tim_idx
461 :
462 : ! ccppized version of the RK scheme
463 : use rk_stratiform, only: rk_stratiform_sedimentation_run
464 : use cloud_particle_sedimentation, only: cloud_particle_sedimentation_run
465 : use rk_stratiform, only: rk_stratiform_detrain_convective_condensate_run
466 : use compute_cloud_fraction, only: compute_cloud_fraction_run
467 : use rk_stratiform, only: rk_stratiform_cloud_fraction_perturbation_run
468 : use rk_stratiform, only: rk_stratiform_condensate_repartioning_run
469 : use rk_stratiform, only: rk_stratiform_external_forcings_run
470 : use prognostic_cloud_water, only: prognostic_cloud_water_run
471 : use rk_stratiform, only: rk_stratiform_prognostic_cloud_water_tendencies_run
472 : use convective_cloud_cover, only: convective_cloud_cover_run
473 : use rk_stratiform, only: rk_stratiform_save_qtlcwat_run
474 : use rk_stratiform, only: rk_stratiform_cloud_optical_properties_run
475 :
476 : use prognostic_cloud_water, only: icritc !REMOVECAM no longer need to be public after CAM is retired
477 : use physconst, only: gravit, rhoh2o, epsilo, latvap, latice, cpair, tmelt, cappa
478 : use physconst, only: rair, rh2o, pi, pref, lapse_rate
479 : use ref_pres, only: trop_cloud_top_lev
480 :
481 : use cloud_optical_properties, only: cldefr
482 : use phys_control, only: cam_physpkg_is
483 : use tropopause, only: tropopause_find_cam
484 : use phys_grid, only: get_rlat_all_p
485 : use constituents, only: cnst_get_ind
486 :
487 : ! Arguments
488 : type(physics_state), intent(in) :: state ! State variables
489 : type(physics_ptend), intent(out) :: ptend_all ! Package tendencies
490 : type(physics_buffer_desc), pointer :: pbuf(:)
491 :
492 : real(r8), intent(in) :: dtime ! Timestep
493 : real(r8), intent(in) :: icefrac (pcols) ! Sea ice fraction (fraction)
494 : real(r8), intent(in) :: landfrac(pcols) ! Land fraction (fraction)
495 : real(r8), intent(in) :: ocnfrac (pcols) ! Ocean fraction (fraction)
496 : real(r8), intent(in) :: snowh(pcols) ! Snow depth over land, water equivalent (m)
497 :
498 : real(r8), intent(in) :: dlf(pcols,pver) ! Detrained water from convection schemes
499 : real(r8), intent(in) :: dlf2(pcols,pver) ! Detrained water from shallow convection scheme
500 : real(r8), intent(in) :: rliq(pcols) ! Vertical integral of liquid not yet in q(ixcldliq)
501 : real(r8), intent(in) :: cmfmc(pcols,pverp) ! Deep + Shallow Convective mass flux [ kg /s/m^2 ]
502 :
503 : real(r8), intent(in) :: ts(pcols) ! Surface temperature
504 : real(r8), intent(in) :: sst(pcols) ! Sea surface temperature
505 : real(r8), intent(in) :: zdu(pcols,pver) ! Detrainment rate from deep convection
506 :
507 : ! Local variables
508 :
509 70392 : type(physics_state) :: state1 ! Local copy of the state variable
510 281568 : type(physics_ptend) :: ptend_loc ! Package tendencies
511 :
512 : integer :: i, k, m
513 : integer :: lchnk ! Chunk identifier
514 : integer :: ncol ! Number of atmospheric columns
515 : integer :: itim_old
516 :
517 : ! Physics buffer fields
518 70392 : real(r8), pointer :: landm(:) ! Land fraction ramped over water
519 :
520 70392 : real(r8), pointer :: prec_str(:) ! [Total] Sfc flux of precip from stratiform [ m/s ]
521 70392 : real(r8), pointer :: snow_str(:) ! [Total] Sfc flux of snow from stratiform [ m/s ]
522 70392 : real(r8), pointer :: prec_sed(:) ! Surface flux of total cloud water from sedimentation
523 70392 : real(r8), pointer :: snow_sed(:) ! Surface flux of cloud ice from sedimentation
524 70392 : real(r8), pointer :: prec_pcw(:) ! Sfc flux of precip from microphysics [ m/s ]
525 70392 : real(r8), pointer :: snow_pcw(:) ! Sfc flux of snow from microphysics [ m/s ]
526 70392 : real(r8), pointer, dimension(:,:) :: deepcu ! deep convection cloud fraction
527 70392 : real(r8), pointer, dimension(:,:) :: shallowcu ! shallow convection cloud fraction
528 :
529 70392 : real(r8), pointer, dimension(:,:) :: qcwat ! Cloud water old q
530 70392 : real(r8), pointer, dimension(:,:) :: tcwat ! Cloud water old temperature
531 70392 : real(r8), pointer, dimension(:,:) :: lcwat ! Cloud liquid water old q
532 70392 : real(r8), pointer, dimension(:,:) :: cld ! Total cloud fraction
533 70392 : real(r8), pointer, dimension(:,:) :: fice ! Cloud ice/water partitioning ratio.
534 70392 : real(r8), pointer, dimension(:,:) :: ast ! Relative humidity cloud fraction
535 70392 : real(r8), pointer, dimension(:,:) :: concld ! Convective cloud fraction
536 70392 : real(r8), pointer, dimension(:,:) :: qme ! rate of cond-evap of condensate (1/s)
537 70392 : real(r8), pointer, dimension(:,:) :: prain ! Total precipitation (rain + snow)
538 70392 : real(r8), pointer, dimension(:,:) :: nevapr ! Evaporation of total precipitation (rain + snow)
539 70392 : real(r8), pointer, dimension(:,:) :: rel ! Liquid effective drop radius (microns)
540 70392 : real(r8), pointer, dimension(:,:) :: rei ! Ice effective drop size (microns)
541 70392 : real(r8), pointer, dimension(:,:) :: wsedl ! Sedimentation velocity of liquid stratus cloud droplet [ m/s ]
542 70392 : real(r8), pointer, dimension(:,:) :: shfrc ! Cloud fraction from shallow convection scheme
543 70392 : real(r8), pointer, dimension(:,:) :: cmfmc_sh ! Shallow convective mass flux (pcols,pverp) [ kg/s/m^2 ]
544 :
545 : real(r8), target :: shfrc_local(pcols,pver)
546 :
547 : ! physics buffer fields for COSP simulator (RK only)
548 70392 : real(r8), pointer, dimension(:,:) :: rkflxprc ! RK grid-box mean flux_large_scale_cloud_rain+snow at interfaces (kg/m2/s)
549 70392 : real(r8), pointer, dimension(:,:) :: rkflxsnw ! RK grid-box mean flux_large_scale_cloud_snow at interfaces (kg/m2/s)
550 :
551 : ! Local variables for stratiform_sediment
552 : real(r8) :: rain(pcols) ! Surface flux of cloud liquid
553 : real(r8) :: pvliq(pcols,pverp) ! Vertical velocity of cloud liquid drops (Pa/s)
554 : real(r8) :: pvice(pcols,pverp) ! Vertical velocity of cloud ice particles (Pa/s)
555 :
556 : ! Local variables for cldfrc
557 :
558 : real(r8) :: cldst(pcols,pver) ! Stratus cloud fraction
559 : real(r8) :: rhcloud(pcols,pver) ! Relative humidity cloud (last timestep)
560 : real(r8) :: relhum(pcols,pver) ! RH, output to determine drh/da
561 : real(r8) :: rhu00(pcols,pver)
562 : real(r8) :: rhdfda(pcols,pver)
563 : real(r8) :: icecldf(pcols,pver) ! Ice cloud fraction
564 : real(r8) :: liqcldf(pcols,pver) ! Liquid cloud fraction (combined into cloud)
565 :
566 : ! Local variables for microphysics
567 :
568 : real(r8) :: qtend(pcols,pver) ! Moisture tendencies
569 : real(r8) :: ttend(pcols,pver) ! Temperature tendencies
570 : real(r8) :: ltend(pcols,pver) ! Cloud liquid water tendencies
571 : real(r8) :: evapheat(pcols,pver) ! Heating rate due to evaporation of precip
572 : real(r8) :: evapsnow(pcols,pver) ! Local evaporation of snow
573 : real(r8) :: prfzheat(pcols,pver) ! Heating rate due to freezing of precip (W/kg)
574 : real(r8) :: meltheat(pcols,pver) ! Heating rate due to phase change of precip
575 : real(r8) :: cmeheat (pcols,pver) ! Heating rate due to phase change of precip
576 : real(r8) :: prodsnow(pcols,pver) ! Local production of snow
577 : real(r8) :: totcw(pcols,pver) ! Total cloud water mixing ratio
578 : real(r8) :: fsnow(pcols,pver) ! Fractional snow production
579 : real(r8) :: repartht(pcols,pver) ! Heating rate due to phase repartition of input precip
580 : real(r8) :: icimr(pcols,pver) ! In cloud ice mixing ratio
581 : real(r8) :: icwmr(pcols,pver) ! In cloud water mixing ratio
582 : real(r8) :: fwaut(pcols,pver)
583 : real(r8) :: fsaut(pcols,pver)
584 : real(r8) :: fracw(pcols,pver)
585 : real(r8) :: fsacw(pcols,pver)
586 : real(r8) :: fsaci(pcols,pver)
587 : real(r8) :: cmeice(pcols,pver) ! Rate of cond-evap of ice within the cloud
588 : real(r8) :: cmeliq(pcols,pver) ! Rate of cond-evap of liq within the cloud
589 : real(r8) :: ice2pr(pcols,pver) ! Rate of conversion of ice to precip
590 : real(r8) :: liq2pr(pcols,pver) ! Rate of conversion of liquid to precip
591 : real(r8) :: liq2snow(pcols,pver) ! Rate of conversion of liquid to snow
592 :
593 : ! Local variables for CFMIP calculations
594 : real(r8) :: mr_lsliq(pcols,pver) ! mixing_ratio_large_scale_cloud_liquid (kg/kg)
595 : real(r8) :: mr_lsice(pcols,pver) ! mixing_ratio_large_scale_cloud_ice (kg/kg)
596 : real(r8) :: mr_ccliq(pcols,pver) ! mixing_ratio_convective_cloud_liquid (kg/kg)
597 : real(r8) :: mr_ccice(pcols,pver) ! mixing_ratio_convective_cloud_ice (kg/kg)
598 :
599 : real(r8) :: pracwo(pcols,pver) ! RK accretion of cloud water by rain (1/s)
600 : real(r8) :: psacwo(pcols,pver) ! RK accretion of cloud water by snow (1/s)
601 : real(r8) :: psacio(pcols,pver) ! RK accretion of cloud ice by snow (1/s)
602 :
603 : real(r8) :: iwc(pcols,pver) ! Grid box average ice water content
604 : real(r8) :: lwc(pcols,pver) ! Grid box average liquid water content
605 :
606 : logical :: lq(pcnst)
607 : integer :: troplev(pcols)
608 : real(r8) :: rlat(pcols)
609 : real(r8) :: dlat(pcols)
610 : real(r8), parameter :: rad2deg = 180._r8/pi
611 :
612 : integer :: top_lev
613 : integer :: ixq
614 :
615 : character(len=512) :: errmsg
616 : integer :: errflg
617 :
618 :
619 : ! ======================================================================
620 :
621 70392 : lchnk = state%lchnk
622 70392 : ncol = state%ncol
623 :
624 70392 : call cnst_get_ind('Q', ixq) ! water vapor index in const array.
625 :
626 70392 : call physics_state_copy(state,state1) ! Copy state to local state1.
627 :
628 : ! Associate pointers with physics buffer fields
629 :
630 70392 : call pbuf_get_field(pbuf, landm_idx, landm) ! smoothed_land_area_fraction
631 :
632 70392 : itim_old = pbuf_old_tim_idx()
633 281568 : call pbuf_get_field(pbuf, qcwat_idx, qcwat, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
634 281568 : call pbuf_get_field(pbuf, tcwat_idx, tcwat, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
635 281568 : call pbuf_get_field(pbuf, lcwat_idx, lcwat, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
636 :
637 281568 : call pbuf_get_field(pbuf, cld_idx, cld, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
638 281568 : call pbuf_get_field(pbuf, concld_idx, concld, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
639 281568 : call pbuf_get_field(pbuf, ast_idx, ast, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
640 :
641 70392 : call pbuf_get_field(pbuf, fice_idx, fice)
642 :
643 70392 : call pbuf_get_field(pbuf, cmfmc_sh_idx, cmfmc_sh)
644 :
645 70392 : call pbuf_get_field(pbuf, prec_str_idx, prec_str)
646 70392 : call pbuf_get_field(pbuf, snow_str_idx, snow_str)
647 70392 : call pbuf_get_field(pbuf, prec_sed_idx, prec_sed)
648 70392 : call pbuf_get_field(pbuf, snow_sed_idx, snow_sed)
649 70392 : call pbuf_get_field(pbuf, prec_pcw_idx, prec_pcw)
650 70392 : call pbuf_get_field(pbuf, snow_pcw_idx, snow_pcw)
651 :
652 70392 : call pbuf_get_field(pbuf, qme_idx, qme )
653 70392 : call pbuf_get_field(pbuf, prain_idx, prain)
654 70392 : call pbuf_get_field(pbuf, nevapr_idx, nevapr)
655 :
656 70392 : call pbuf_get_field(pbuf, rel_idx, rel)
657 70392 : call pbuf_get_field(pbuf, rei_idx, rei)
658 :
659 70392 : call pbuf_get_field(pbuf, wsedl_idx, wsedl)
660 :
661 : ! check that qcwat and tcwat were initialized; if not then do it now.
662 70392 : if (qcwat(1,1) == huge(1._r8)) then
663 1354104 : qcwat(:ncol,:) = state%q(:ncol,:,1)
664 : end if
665 70392 : if (tcwat(1,1) == huge(1._r8)) then
666 1354104 : tcwat(:ncol,:) = state%t(:ncol,:)
667 : end if
668 :
669 70392 : if ( do_psrhmin ) then
670 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
671 0 : troplev(:) = 0
672 : !REMOVECAM_END
673 0 : call tropopause_find_cam(state, troplev)
674 : endif
675 :
676 : ! dlat required for CCPPized scheme
677 : ! have to use get_rlat_all_p here because dlat is expected to be r8
678 : ! and phys_grid does not implement dlat in r8, only int
679 70392 : call get_rlat_all_p(lchnk,ncol,rlat)
680 1090992 : dlat(:ncol) = rlat(:ncol)*rad2deg
681 :
682 : ! ------------- !
683 : ! Sedimentation !
684 : ! ------------- !
685 :
686 : ! Allow the cloud liquid drops and ice particles to sediment.
687 : ! This is done before adding convectively detrained cloud water,
688 : ! because the phase of the detrained water is unknown.
689 :
690 70392 : call t_startf('stratiform_sediment')
691 :
692 70392 : lq(:) = .FALSE.
693 70392 : lq(ixq) = .TRUE.
694 70392 : lq(ixcldice) = .TRUE.
695 70392 : lq(ixcldliq) = .TRUE.
696 70392 : call physics_ptend_init(ptend_loc, state%psetcols, 'pcwsediment', ls=.true., lq=lq)! Initialize local ptend type
697 :
698 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
699 70392 : pvliq(:,:) = 0._r8
700 70392 : pvice(:,:) = 0._r8
701 70392 : rain(:) = 0._r8
702 1196664 : snow_sed(:) = 0._r8
703 : !REMOVECAM_END
704 : ! Call CCPP-ized subroutine for cld_sediment_vel and cld_sediment_tend combined
705 : call cloud_particle_sedimentation_run( &
706 : ncol = ncol, &
707 : pver = pver, &
708 : pverp = pverp, &
709 : dtime = dtime, &
710 : tmelt = tmelt, &
711 : gravit = gravit, &
712 : latvap = latvap, &
713 : latice = latice, &
714 : rair = rair, &
715 : rhoh2o = rhoh2o, &
716 : icritc = icritc, &
717 0 : pint = state1%pint(:ncol,:), &
718 0 : pmid = state1%pmid(:ncol,:), &
719 0 : pdel = state1%pdel(:ncol,:), &
720 0 : t = state1%t(:ncol,:), &
721 0 : cloud = cld(:ncol,:), &
722 0 : icefrac = icefrac(:ncol), &
723 0 : landfrac = landfrac(:ncol), &
724 0 : ocnfrac = ocnfrac(:ncol), &
725 0 : cldliq = state1%q(:ncol,:,ixcldliq), &
726 0 : cldice = state1%q(:ncol,:,ixcldice), &
727 0 : snowh = snowh(:ncol), &
728 0 : landm = landm(:ncol), &
729 : ! Outputs
730 0 : pvliq = pvliq(:ncol,:), &
731 0 : pvice = pvice(:ncol,:), &
732 0 : liqtend = ptend_loc%q(:ncol,:,ixcldliq), &
733 0 : icetend = ptend_loc%q(:ncol,:,ixcldice), &
734 0 : wvtend = ptend_loc%q(:ncol,:,ixq), &
735 0 : htend = ptend_loc%s(:ncol,:), &
736 0 : sfliq = rain(:ncol), & ! mass units (not precip units)
737 0 : sfice = snow_sed(:ncol), & ! precip units (as stored in pbuf)
738 : errmsg = errmsg, &
739 70392 : errflg = errflg)
740 :
741 70392 : if(errflg /= 0) then
742 0 : call endrun('rk_stratiform_tend:' // errmsg)
743 : endif
744 :
745 28436184 : wsedl(:ncol,:pver) = pvliq(:ncol,:pver)/gravit/(state1%pmid(:ncol,:pver)/(rair*state1%t(:ncol,:pver)))
746 :
747 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
748 1196664 : prec_sed(:) = 0._r8
749 1196664 : prec_str(:) = 0._r8
750 1196664 : snow_str(:) = 0._r8
751 : !REMOVECAM_END
752 : ! Call the CCPPized scheme for accumulation of total stratiform fluxes [m s-1]
753 : call rk_stratiform_sedimentation_run( &
754 : ncol = ncol, &
755 0 : sfliq = rain(:ncol), & ! in from sedimentation
756 0 : snow_sed = snow_sed(:ncol), & ! in from sedimentation
757 0 : prec_sed = prec_sed(:ncol), & ! out -- converted to precip units
758 0 : prec_str = prec_str(:ncol), & ! out -- initial assignment for accumulation
759 0 : snow_str = snow_str(:ncol), & ! out -- initial assignment for accumulation
760 : errmsg = errmsg, &
761 70392 : errflg = errflg)
762 :
763 70392 : if(errflg /= 0) then
764 0 : call endrun('rk_stratiform_tend:' // errmsg)
765 : endif
766 :
767 : ! Record history variables
768 70392 : call outfld( 'DQSED' ,ptend_loc%q(:,:,ixq) , pcols,lchnk )
769 70392 : call outfld( 'DISED' ,ptend_loc%q(:,:,ixcldice), pcols,lchnk )
770 70392 : call outfld( 'DLSED' ,ptend_loc%q(:,:,ixcldliq), pcols,lchnk )
771 70392 : call outfld( 'HSED' ,ptend_loc%s , pcols,lchnk )
772 70392 : call outfld( 'PRECSED' ,prec_sed , pcols,lchnk )
773 70392 : call outfld( 'SNOWSED' ,snow_sed , pcols,lchnk )
774 1196664 : call outfld( 'RAINSED' ,rain/1000._r8 , pcols,lchnk ) ! convert from kg m-2 s-1 to m s-1 (precip units) for output
775 :
776 : ! Add tendency from this process to tend from other processes here
777 70392 : call physics_ptend_init(ptend_all, state%psetcols, 'stratiform')
778 70392 : call physics_ptend_sum( ptend_loc, ptend_all, ncol )
779 :
780 : ! Update physics state type state1 with ptend_loc
781 70392 : call physics_update( state1, ptend_loc, dtime )
782 :
783 70392 : call t_stopf('stratiform_sediment')
784 :
785 : ! ----------------------------------------------------------------------------- !
786 : ! Detrainment of convective condensate into the environment or stratiform cloud !
787 : ! ----------------------------------------------------------------------------- !
788 :
789 : ! Put all of the detraining cloud water from convection into the large scale cloud.
790 : ! It all goes in liquid for the moment.
791 : ! Strictly speaking, this approach is detraining all the cconvective water into
792 : ! the environment, not the large-scale cloud.
793 :
794 70392 : lq(:) = .FALSE.
795 70392 : lq(ixcldliq) = .TRUE.
796 70392 : call physics_ptend_init( ptend_loc, state1%psetcols, 'pcwdetrain', lq=lq)
797 :
798 : call rk_stratiform_detrain_convective_condensate_run( &
799 : ncol = ncol, &
800 0 : dlf = dlf(:ncol,:), &
801 0 : rliq = rliq(:ncol), &
802 0 : prec_str = prec_str(:ncol), &
803 0 : tend_cldliq = ptend_loc%q(:ncol,:,ixcldliq), &
804 : errmsg = errmsg, &
805 70392 : errflg = errflg)
806 :
807 70392 : if(errflg /= 0) then
808 0 : call endrun('rk_stratiform_tend:' // errmsg)
809 : endif
810 :
811 70392 : call outfld( 'ZMDLF', dlf, pcols, lchnk )
812 70392 : call outfld( 'SHDLF', dlf2, pcols, lchnk )
813 :
814 : ! Add hie detrainment tendency to tend from the other prior processes
815 :
816 70392 : call physics_ptend_sum( ptend_loc, ptend_all, ncol )
817 70392 : call physics_update( state1, ptend_loc, dtime )
818 :
819 : ! -------------------------------------- !
820 : ! Computation of Various Cloud Fractions !
821 : ! -------------------------------------- !
822 :
823 : ! ----------------------------------------------------------------------------- !
824 : ! Treatment of cloud fraction in CAM4 and CAM5 differs !
825 : ! (1) CAM4 !
826 : ! . Cumulus AMT = Deep Cumulus AMT ( empirical fcn of mass flux ) + !
827 : ! Shallow Cumulus AMT ( empirical fcn of mass flux ) !
828 : ! . Stratus AMT = max( RH stratus AMT, Stability Stratus AMT ) !
829 : ! . Cumulus and Stratus are 'minimally' overlapped without hierarchy. !
830 : ! . Cumulus LWC,IWC is assumed to be the same as Stratus LWC,IWC !
831 : ! (2) CAM5 !
832 : ! . Cumulus AMT = Deep Cumulus AMT ( empirical fcn of mass flux ) + !
833 : ! Shallow Cumulus AMT ( internally fcn of mass flux and w ) !
834 : ! . Stratus AMT = fcn of environmental-mean RH ( no Stability Stratus ) !
835 : ! . Cumulus and Stratus are non-overlapped with higher priority on Cumulus !
836 : ! . Cumulus ( both Deep and Shallow ) has its own LWC and IWC. !
837 : ! ----------------------------------------------------------------------------- !
838 :
839 70392 : if( use_shfrc ) then
840 0 : call pbuf_get_field(pbuf, shfrc_idx, shfrc )
841 : else
842 70392 : shfrc=>shfrc_local
843 31183656 : shfrc(:,:) = 0._r8
844 : endif
845 :
846 70392 : call pbuf_get_field(pbuf, sh_frac_idx, shallowcu )
847 70392 : call pbuf_get_field(pbuf, dp_frac_idx, deepcu )
848 :
849 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
850 31183656 : shallowcu(:,:) = 0._r8
851 31183656 : deepcu(:,:) = 0._r8
852 31183656 : concld(:,:) = 0._r8
853 : !REMOVECAM_END
854 :
855 : ! compute convective cloud fraction using CCPP-ized subroutine
856 : call convective_cloud_cover_run( &
857 : ncol = ncol, &
858 : pver = pver, &
859 : top_lev_cloudphys = 1, & ! CAM4 macrophysics.
860 : use_shfrc = use_shfrc, &
861 0 : shfrc = shfrc(:ncol,:), &
862 0 : cmfmc_total = cmfmc(:ncol,:), &
863 0 : cmfmc_sh = cmfmc_sh(:ncol,:), &
864 0 : shallowcu = shallowcu(:ncol,:), &
865 0 : deepcu = deepcu(:ncol,:), &
866 0 : concld = concld(:ncol,:), &
867 : errmsg = errmsg, &
868 70392 : errflg = errflg)
869 :
870 : ! write out convective cloud fraction diagnostic.
871 70392 : call outfld( 'SH_CLD ', shallowcu , pcols, lchnk )
872 70392 : call outfld( 'DP_CLD ', deepcu , pcols, lchnk )
873 70392 : call outfld( 'CONCLD ', concld , pcols, lchnk )
874 :
875 : ! Stratus ('ast' = max(alst,aist)) and total cloud fraction ('cld = ast + concld')
876 : ! will be computed using this updated 'concld' in the stratiform macrophysics
877 : ! scheme (mmacro_pcond) later below.
878 :
879 70392 : call t_startf("cldfrc")
880 :
881 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
882 31183656 : cld(:,:) = 0._r8
883 70392 : rhcloud(:,:) = 0._r8
884 70392 : cldst(:,:) = 0._r8
885 70392 : rhu00(:,:) = 0._r8
886 70392 : icecldf(:,:) = 0._r8
887 70392 : liqcldf(:,:) = 0._r8
888 70392 : relhum(:,:) = 0._r8
889 : !REMOVECAM_END
890 :
891 : ! call CCPPized cloud fraction scheme
892 : call compute_cloud_fraction_run( &
893 : ncol = ncol, &
894 : pver = pver, &
895 : cappa = cappa, &
896 : gravit = gravit, &
897 : rair = rair, &
898 : tmelt = tmelt, &
899 : pref = pref, &
900 : lapse_rate = lapse_rate, &
901 : top_lev_cloudphys = 1, & ! CAM4 macrophysics.
902 0 : pmid = state1%pmid(:ncol,:), &
903 0 : ps = state1%pint(:ncol,pverp), &
904 0 : temp = state1%t(:ncol,:), &
905 0 : sst = sst(:ncol), &
906 0 : q = state1%q(:ncol,:,ixq), &
907 0 : cldice = state1%q(:ncol,:,ixcldice), &
908 0 : phis = state1%phis(:ncol), &
909 0 : shallowcu = shallowcu(:ncol,:), &
910 0 : deepcu = deepcu(:ncol,:), &
911 0 : concld = concld(:ncol,:), &
912 0 : landfrac = landfrac(:ncol), &
913 0 : ocnfrac = ocnfrac(:ncol), &
914 0 : snowh = snowh(:ncol), &
915 : rhpert_flag = .false., & ! below output:
916 0 : cloud = cld(:ncol, :), &
917 0 : rhcloud = rhcloud(:ncol, :), &
918 0 : cldst = cldst(:ncol,:), &
919 0 : rhu00 = rhu00(:ncol,:), &
920 0 : icecldf = icecldf(:ncol,:), &
921 0 : liqcldf = liqcldf(:ncol,:), &
922 0 : relhum = relhum(:ncol,:), &
923 : errmsg = errmsg, &
924 70392 : errflg = errflg)
925 :
926 : ! The AST history output is from cld from this call of compute_cloud_fraction,
927 : ! and not in the subsequent calls. (hplin, 3/6/25)
928 70392 : call outfld( 'AST', cld, pcols, lchnk )
929 :
930 : ! Re-calculate cloud with perturbed rh add call cldfrc to estimate rhdfda.
931 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
932 70392 : rhdfda(:,:) = 0._r8
933 : !REMOVECAM_END
934 :
935 : ! Re-calculate cloud with perturbed rh add call cldfrc to estimate rhdfda.
936 : call rk_stratiform_cloud_fraction_perturbation_run( &
937 : ncol = ncol, &
938 : pver = pver, &
939 : cappa = cappa, &
940 : gravit = gravit, &
941 : rair = rair, &
942 : tmelt = tmelt, &
943 : pref = pref, &
944 : lapse_rate = lapse_rate, &
945 : top_lev_cloudphys = 1, & ! CAM4 macrophysics.
946 0 : pmid = state1%pmid(:ncol,:), &
947 0 : ps = state1%pint(:ncol,pverp), &
948 0 : temp = state1%t(:ncol,:), &
949 0 : sst = sst(:ncol), &
950 0 : q_wv = state1%q(:ncol,:,ixq), &
951 0 : cldice = state1%q(:ncol,:,ixcldice), &
952 0 : phis = state1%phis(:ncol), &
953 0 : shallowcu = shallowcu(:ncol,:), &
954 0 : deepcu = deepcu(:ncol,:), &
955 0 : concld = concld(:ncol,:), &
956 0 : landfrac = landfrac(:ncol), &
957 0 : ocnfrac = ocnfrac(:ncol), &
958 0 : snowh = snowh(:ncol), &
959 0 : cloud = cld(:ncol,:), & ! unperturbed input
960 0 : relhum = relhum(:ncol,:), & ! unperturbed input
961 0 : rhu00 = rhu00(:ncol,:), & ! unperturbed input - output
962 0 : rhdfda = rhdfda(:ncol,:), & ! output for prognostic_cloud_water
963 : errmsg = errmsg, &
964 70392 : errflg = errflg)
965 :
966 70392 : call t_stopf("cldfrc")
967 :
968 : ! ---------------------------------------------- !
969 : ! Stratiform Cloud Macrophysics and Microphysics !
970 : ! ---------------------------------------------- !
971 :
972 70392 : call t_startf('stratiform_microphys')
973 :
974 : !-------------------------------------------------------
975 : ! Non-micro and non-macrophysical external advective forcings to compute net condensation rate.
976 : ! Note that advective forcing of condensate is aggregated into liquid phase.
977 : !-------------------------------------------------------
978 :
979 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
980 70392 : qtend(:,:) = 0._r8
981 70392 : ttend(:,:) = 0._r8
982 70392 : ltend(:,:) = 0._r8
983 : !REMOVECAM_END
984 :
985 : call rk_stratiform_external_forcings_run( &
986 : ncol = ncol, &
987 : pver = pver, &
988 : dtime = dtime, &
989 0 : t = state1%t(:ncol,:), &
990 0 : q_wv = state1%q(:ncol,:pver,ixq), &
991 0 : cldice = state1%q(:ncol,:pver,ixcldice), & ! using before-repartitioning values for bit-for-bit.
992 0 : cldliq = state1%q(:ncol,:pver,ixcldliq), & ! using before-repartitioning values for bit-for-bit.
993 0 : qcwat = qcwat(:ncol,:), &
994 0 : tcwat = tcwat(:ncol,:), &
995 0 : lcwat = lcwat(:ncol,:), &
996 0 : qtend = qtend(:ncol,:), &
997 0 : ttend = ttend(:ncol,:), &
998 0 : ltend = ltend(:ncol,:), &
999 : errmsg = errmsg, &
1000 70392 : errflg = errflg)
1001 :
1002 : !-------------------------------------------------------
1003 : ! Define fractional amount of stratus condensate and precipitation in ice phase.
1004 : ! This uses a ramp ( -30 ~ -10 for fice, -5 ~ 0 for fsnow ).
1005 : ! The ramp within convective cloud may be different
1006 : !
1007 : ! This ice fraction is used both in the repartitioning, and also in prognostic_cloud_water.
1008 : !-------------------------------------------------------
1009 :
1010 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
1011 31183656 : fice(:,:) = 0._r8
1012 70392 : fsnow(:,:) = 0._r8
1013 : !REMOVECAM_END
1014 70392 : top_lev = 1
1015 70392 : call cloud_fraction_fice_run(ncol, state1%t(:ncol,:), tmelt, top_lev, pver, fice(:ncol,:), fsnow(:ncol,:), errmsg, errflg)
1016 :
1017 : !-------------------------------------------------------
1018 : ! Compute Stratiform Macro-Microphysical Tendencies
1019 : !-------------------------------------------------------
1020 :
1021 : ! Add rain and snow fluxes as output variables from pcond, and into physics buffer
1022 70392 : call pbuf_get_field(pbuf, ls_flxprc_idx, rkflxprc)
1023 70392 : call pbuf_get_field(pbuf, ls_flxsnw_idx, rkflxsnw)
1024 :
1025 70392 : call t_startf('pcond')
1026 :
1027 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
1028 31183656 : qme(:,:) = 0._r8
1029 31183656 : prain(:,:) = 0._r8
1030 70392 : prodsnow(:,:) = 0._r8
1031 31183656 : nevapr(:,:) = 0._r8
1032 70392 : evapsnow(:,:) = 0._r8
1033 70392 : evapheat(:,:) = 0._r8
1034 70392 : prfzheat(:,:) = 0._r8
1035 70392 : meltheat(:,:) = 0._r8
1036 1196664 : prec_pcw(:) = 0._r8
1037 1196664 : snow_pcw(:) = 0._r8
1038 70392 : ice2pr(:,:) = 0._r8
1039 70392 : liq2pr(:,:) = 0._r8
1040 70392 : liq2snow(:,:) = 0._r8
1041 32380320 : rkflxprc(:,:) = 0._r8
1042 32380320 : rkflxsnw(:,:) = 0._r8
1043 70392 : pracwo(:,:) = 0._r8
1044 70392 : psacwo(:,:) = 0._r8
1045 70392 : psacio(:,:) = 0._r8
1046 : !REMOVECAM_END
1047 : ! call CCPP-ized subroutine
1048 : call prognostic_cloud_water_run( &
1049 : ncol = ncol, &
1050 : pver = pver, &
1051 : top_lev = trop_cloud_top_lev, &
1052 : deltat = dtime, &
1053 : iulog = iulog, &
1054 : pi = pi, &
1055 : gravit = gravit, &
1056 : rh2o = rh2o, &
1057 : epsilo = epsilo, &
1058 : latvap = latvap, &
1059 : latice = latice, &
1060 : cpair = cpair, &
1061 0 : dlat = dlat(:ncol), &
1062 0 : pmid = state1%pmid(:ncol,:), &
1063 0 : pdel = state1%pdel(:ncol,:), &
1064 0 : zi = state1%zi(:ncol,:), &
1065 0 : troplev = troplev(:ncol), &
1066 0 : ttend = ttend(:ncol,:), &
1067 0 : tn = state1%t(:ncol,:), &
1068 0 : qtend = qtend(:ncol,:), &
1069 0 : qn = state1%q(:ncol,:,ixq), &
1070 0 : ltend = ltend(:ncol,:), &
1071 0 : cldice = state1%q(:ncol,:pver,ixcldice), & ! for bit-for-bit, use values pre-repartitioning.
1072 0 : cldliq = state1%q(:ncol,:pver,ixcldliq), & ! for bit-for-bit, use values pre-repartitioning.
1073 : !cwat = totcw(:ncol,:), &
1074 0 : omega = state1%omega(:ncol,:), &
1075 0 : cldn = cld(:ncol,:), &
1076 0 : fice = fice(:ncol,:), &
1077 0 : fsnow = fsnow(:ncol,:), &
1078 0 : rhdfda = rhdfda(:ncol,:), &
1079 0 : rhu00 = rhu00(:ncol,:), &
1080 0 : landm = landm(:ncol), &
1081 0 : seaicef = icefrac(:ncol), &
1082 0 : snowh = snowh(:ncol), &
1083 : ! below output
1084 0 : qme = qme(:ncol,:), &
1085 0 : prodprec = prain(:ncol,:), &
1086 0 : prodsnow = prodsnow(:ncol,:), & ! ignored in rk.
1087 0 : evapprec = nevapr(:ncol,:), &
1088 0 : evapsnow = evapsnow(:ncol,:), &
1089 0 : evapheat = evapheat(:ncol,:), &
1090 0 : prfzheat = prfzheat(:ncol,:), &
1091 0 : meltheat = meltheat(:ncol,:), &
1092 0 : precip = prec_pcw(:ncol), &
1093 0 : snowab = snow_pcw(:ncol), &
1094 0 : ice2pr = ice2pr(:ncol,:), &
1095 0 : liq2pr = liq2pr(:ncol,:), &
1096 0 : liq2snow = liq2snow(:ncol,:), &
1097 0 : lsflxprc = rkflxprc(:ncol,:), &
1098 0 : lsflxsnw = rkflxsnw(:ncol,:), &
1099 0 : pracwo = pracwo(:ncol,:), &
1100 0 : psacwo = psacwo(:ncol,:), &
1101 0 : psacio = psacio(:ncol,:), &
1102 0 : fwaut = fwaut(:ncol,:), &
1103 0 : fsaut = fsaut(:ncol,:), &
1104 0 : fracw = fracw(:ncol,:), &
1105 0 : fsacw = fsacw(:ncol,:), &
1106 0 : fsaci = fsaci(:ncol,:), &
1107 : errmsg = errmsg, &
1108 70392 : errflg = errflg)
1109 70392 : if(errflg /= 0) then
1110 0 : call endrun('rk_stratiform_tend:' // errmsg)
1111 : endif
1112 :
1113 70392 : call t_stopf('pcond')
1114 :
1115 : !-------------------------------------------------------
1116 : ! Perform repartitioning of stratiform condensate.
1117 : ! Corresponding heating tendency will be added later.
1118 : !-------------------------------------------------------
1119 :
1120 70392 : lq(:) = .FALSE.
1121 70392 : lq(ixcldice) = .true.
1122 70392 : lq(ixcldliq) = .true.
1123 70392 : call physics_ptend_init( ptend_loc, state1%psetcols, 'cldwat-repartition', lq=lq )
1124 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
1125 70392 : repartht(:,:) = 0._r8
1126 : !REMOVECAM_END
1127 :
1128 : call rk_stratiform_condensate_repartioning_run( &
1129 : ncol = ncol, &
1130 : pver = pver, &
1131 : dtime = dtime, &
1132 : latice = latice, &
1133 0 : cldice = state1%q(:ncol,:,ixcldice), &
1134 0 : cldliq = state1%q(:ncol,:,ixcldliq), &
1135 0 : fice = fice(:ncol,:), & ! below output
1136 0 : tend_cldice = ptend_loc%q(:ncol,:,ixcldice), &
1137 0 : tend_cldliq = ptend_loc%q(:ncol,:,ixcldliq), &
1138 0 : repartht = repartht(:ncol,:), &
1139 : errmsg = errmsg, &
1140 70392 : errflg = errflg)
1141 :
1142 70392 : call outfld( 'REPARTICE', ptend_loc%q(:,:,ixcldice), pcols, lchnk )
1143 70392 : call outfld( 'REPARTLIQ', ptend_loc%q(:,:,ixcldliq), pcols, lchnk )
1144 :
1145 : ! note: to restore old way of calculating repartht
1146 : ! because apparently it causes errors greater than roundoff when accumulated
1147 : ! in certain compsets of the CAM regression tests.
1148 : ! in any case, I do not believe this is a good representation. It sets repartht
1149 : ! to an intermediate value with the wrong units...
1150 : ! repartht(:ncol,:pver) = state1%q(:ncol,:pver,ixcldice)
1151 : ! still bit differences:
1152 : ! repartht(:ncol,:pver) = (latice/dtime) * state1%q(:ncol,:pver,ixcldice)
1153 : ! // note
1154 :
1155 70392 : call physics_ptend_sum( ptend_loc, ptend_all, ncol )
1156 70392 : call physics_update( state1, ptend_loc, dtime )
1157 :
1158 : ! note: second part of restoring old way of calculating repartht
1159 : ! ...then does the difference here to determine repartition heating
1160 : ! repartht(:ncol,:pver) = (latice/dtime) * ( state1%q(:ncol,:pver,ixcldice) - repartht(:ncol,:pver) )
1161 : ! still bit differences if we switch to maintain proper units:
1162 : ! repartht(:ncol,:pver) = (latice/dtime) * state1%q(:ncol,:pver,ixcldice) - repartht(:ncol,:pver)
1163 :
1164 :
1165 : !-------------------------------------------------------
1166 : ! APPLY Stratiform Macro-Microphysical Tendencies
1167 : !-------------------------------------------------------
1168 :
1169 70392 : lq(:) = .FALSE.
1170 70392 : lq(ixq) = .true.
1171 70392 : lq(ixcldice) = .true.
1172 70392 : lq(ixcldliq) = .true.
1173 70392 : call physics_ptend_init( ptend_loc, state1%psetcols, 'cldwat', ls=.true., lq=lq)
1174 :
1175 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
1176 70392 : cmeheat(:,:) = 0._r8
1177 70392 : cmeice(:,:) = 0._r8
1178 70392 : cmeliq(:,:) = 0._r8
1179 : !REMOVECAM_END
1180 : ! call CCPP-ized subroutine to compute tendencies from prognostic_cloud_water
1181 : ! (to heating rate, q, cldice, cldliq)
1182 : call rk_stratiform_prognostic_cloud_water_tendencies_run( &
1183 : ncol = ncol, &
1184 : pver = pver, &
1185 : dtime = dtime, &
1186 : latvap = latvap, &
1187 : latice = latice, &
1188 0 : qme = qme(:ncol,:), &
1189 0 : fice = fice(:ncol,:), &
1190 0 : evapheat = evapheat(:ncol,:), &
1191 0 : prfzheat = prfzheat(:ncol,:), &
1192 0 : meltheat = meltheat(:ncol,:), &
1193 0 : repartht = repartht(:ncol,:), &
1194 0 : evapprec = nevapr(:ncol,:), &
1195 0 : ice2pr = ice2pr(:ncol,:), &
1196 0 : liq2pr = liq2pr(:ncol,:), &
1197 0 : cmeheat = cmeheat(:ncol,:), &
1198 0 : cmeice = cmeice(:ncol,:), &
1199 0 : cmeliq = cmeliq(:ncol,:), &
1200 0 : prec_pcw = prec_pcw(:ncol), &
1201 0 : snow_pcw = snow_pcw(:ncol), &
1202 0 : prec_str = prec_str(:ncol), &
1203 0 : snow_str = snow_str(:ncol), &
1204 0 : tend_s = ptend_loc%s(:ncol,:), &
1205 0 : tend_q = ptend_loc%q(:ncol,:,ixq), &
1206 0 : tend_cldice = ptend_loc%q(:ncol,:,ixcldice), &
1207 0 : tend_cldliq = ptend_loc%q(:ncol,:,ixcldliq), &
1208 : errmsg = errmsg, &
1209 70392 : errflg = errflg)
1210 :
1211 : ! Record history variables
1212 70392 : call outfld( 'CME' , qme, pcols, lchnk )
1213 70392 : call outfld( 'PRODPREC' , prain, pcols, lchnk )
1214 70392 : call outfld( 'EVAPPREC' , nevapr, pcols, lchnk )
1215 70392 : call outfld( 'EVAPSNOW' , evapsnow, pcols, lchnk )
1216 :
1217 70392 : call outfld( 'FWAUT' , fwaut, pcols, lchnk )
1218 70392 : call outfld( 'FSAUT' , fsaut, pcols, lchnk )
1219 70392 : call outfld( 'FRACW' , fracw, pcols, lchnk )
1220 70392 : call outfld( 'FSACW' , fsacw, pcols, lchnk )
1221 70392 : call outfld( 'FSACI' , fsaci, pcols, lchnk )
1222 :
1223 70392 : call outfld( 'PCSNOW' , snow_pcw, pcols, lchnk )
1224 70392 : call outfld( 'FICE' , fice, pcols, lchnk )
1225 70392 : call outfld( 'CMEICE' , cmeice, pcols, lchnk )
1226 70392 : call outfld( 'CMELIQ' , cmeliq, pcols, lchnk )
1227 70392 : call outfld( 'ICE2PR' , ice2pr, pcols, lchnk )
1228 70392 : call outfld( 'LIQ2PR' , liq2pr, pcols, lchnk )
1229 70392 : call outfld( 'HPROGCLD', ptend_loc%s, pcols, lchnk )
1230 70392 : call outfld( 'HEVAP ', evapheat, pcols, lchnk )
1231 70392 : call outfld( 'HMELT' , meltheat, pcols, lchnk )
1232 70392 : call outfld( 'HCME' , cmeheat , pcols, lchnk )
1233 70392 : call outfld( 'HFREEZ' , prfzheat, pcols, lchnk )
1234 70392 : call outfld( 'HREPART' , repartht, pcols, lchnk )
1235 70392 : call outfld('LS_FLXPRC', rkflxprc, pcols, lchnk )
1236 70392 : call outfld('LS_FLXSNW', rkflxsnw, pcols, lchnk )
1237 70392 : call outfld('PRACWO' , pracwo, pcols, lchnk )
1238 70392 : call outfld('PSACWO' , psacwo, pcols, lchnk )
1239 70392 : call outfld('PSACIO' , psacio, pcols, lchnk )
1240 :
1241 : ! initialize local variables for CFMIP diagnostics
1242 28436184 : mr_ccliq(1:ncol,1:pver) = 0._r8
1243 28436184 : mr_ccice(1:ncol,1:pver) = 0._r8
1244 28436184 : mr_lsliq(1:ncol,1:pver) = 0._r8
1245 28436184 : mr_lsice(1:ncol,1:pver) = 0._r8
1246 :
1247 1900584 : do k=1,pver
1248 28436184 : do i=1,ncol
1249 28365792 : if (cld(i,k) .gt. 0._r8) then
1250 6204147 : mr_ccliq(i,k) = (state%q(i,k,ixcldliq)/cld(i,k))*concld(i,k)
1251 6204147 : mr_ccice(i,k) = (state%q(i,k,ixcldice)/cld(i,k))*concld(i,k)
1252 6204147 : mr_lsliq(i,k) = (state%q(i,k,ixcldliq)/cld(i,k))*(cld(i,k)-concld(i,k))
1253 6204147 : mr_lsice(i,k) = (state%q(i,k,ixcldice)/cld(i,k))*(cld(i,k)-concld(i,k))
1254 : else
1255 20331453 : mr_ccliq(i,k) = 0._r8
1256 20331453 : mr_ccice(i,k) = 0._r8
1257 20331453 : mr_lsliq(i,k) = 0._r8
1258 20331453 : mr_lsice(i,k) = 0._r8
1259 : end if
1260 : end do
1261 : end do
1262 :
1263 70392 : call outfld( 'CLDLIQSTR ', mr_lsliq, pcols, lchnk )
1264 70392 : call outfld( 'CLDICESTR ', mr_lsice, pcols, lchnk )
1265 70392 : call outfld( 'CLDLIQCON ', mr_ccliq, pcols, lchnk )
1266 70392 : call outfld( 'CLDICECON ', mr_ccice, pcols, lchnk )
1267 :
1268 : ! ------------------------------- !
1269 : ! Update microphysical tendencies !
1270 : ! ------------------------------- !
1271 :
1272 70392 : call physics_ptend_sum( ptend_loc, ptend_all, ncol )
1273 70392 : call physics_update( state1, ptend_loc, dtime )
1274 :
1275 70392 : call t_startf("cldfrc")
1276 :
1277 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
1278 31183656 : cld(:,:) = 0._r8
1279 70392 : rhcloud(:,:) = 0._r8
1280 70392 : cldst(:,:) = 0._r8
1281 70392 : rhu00(:,:) = 0._r8
1282 70392 : icecldf(:,:) = 0._r8
1283 70392 : liqcldf(:,:) = 0._r8
1284 70392 : relhum(:,:) = 0._r8
1285 : !REMOVECAM_END
1286 :
1287 : ! call CCPPized cloud fraction scheme
1288 : call compute_cloud_fraction_run( &
1289 : ncol = ncol, &
1290 : pver = pver, &
1291 : cappa = cappa, &
1292 : gravit = gravit, &
1293 : rair = rair, &
1294 : tmelt = tmelt, &
1295 : pref = pref, &
1296 : lapse_rate = lapse_rate, &
1297 : top_lev_cloudphys = 1, & ! CAM4 macrophysics.
1298 0 : pmid = state1%pmid(:ncol,:), &
1299 0 : ps = state1%pint(:ncol,pverp), &
1300 0 : temp = state1%t(:ncol,:), &
1301 0 : sst = sst(:ncol), &
1302 0 : q = state1%q(:ncol,:,ixq), &
1303 0 : cldice = state1%q(:ncol,:,ixcldice), &
1304 0 : phis = state1%phis(:ncol), &
1305 0 : shallowcu = shallowcu(:ncol,:), &
1306 0 : deepcu = deepcu(:ncol,:), &
1307 0 : concld = concld(:ncol,:), &
1308 0 : landfrac = landfrac(:ncol), &
1309 0 : ocnfrac = ocnfrac(:ncol), &
1310 0 : snowh = snowh(:ncol), &
1311 : rhpert_flag = .false., & ! below output:
1312 0 : cloud = cld(:ncol, :), &
1313 0 : rhcloud = rhcloud(:ncol, :), &
1314 0 : cldst = cldst(:ncol,:), &
1315 0 : rhu00 = rhu00(:ncol,:), &
1316 0 : icecldf = icecldf(:ncol,:), &
1317 0 : liqcldf = liqcldf(:ncol,:), &
1318 0 : relhum = relhum(:ncol,:), &
1319 : errmsg = errmsg, &
1320 70392 : errflg = errflg)
1321 70392 : call t_stopf("cldfrc")
1322 :
1323 70392 : call outfld( 'CLDST ', cldst, pcols, lchnk )
1324 :
1325 1900584 : do k = 1, pver
1326 28436184 : do i = 1, ncol
1327 26535600 : iwc(i,k) = state1%q(i,k,ixcldice)*state1%pmid(i,k)/(rair*state1%t(i,k))
1328 26535600 : lwc(i,k) = state1%q(i,k,ixcldliq)*state1%pmid(i,k)/(rair*state1%t(i,k))
1329 26535600 : icimr(i,k) = state1%q(i,k,ixcldice) / max(0.01_r8,rhcloud(i,k))
1330 28365792 : icwmr(i,k) = state1%q(i,k,ixcldliq) / max(0.01_r8,rhcloud(i,k))
1331 : end do
1332 : end do
1333 :
1334 70392 : call outfld( 'IWC' , iwc, pcols, lchnk )
1335 70392 : call outfld( 'LWC' , lwc, pcols, lchnk )
1336 70392 : call outfld( 'ICIMR' , icimr, pcols, lchnk )
1337 70392 : call outfld( 'ICWMR' , icwmr, pcols, lchnk )
1338 :
1339 70392 : call t_stopf('stratiform_microphys')
1340 :
1341 : ! Save variables for use in the macrophysics at the next time step
1342 : call rk_stratiform_save_qtlcwat_run( &
1343 : ncol = ncol, &
1344 : pver = pver, &
1345 0 : t = state1%t(:ncol,:), &
1346 0 : q_wv = state1%q(:ncol,:,ixq), &
1347 0 : cldice = state1%q(:ncol,:,ixcldice), &
1348 0 : cldliq = state1%q(:ncol,:,ixcldliq), &
1349 0 : qcwat = qcwat(:ncol,:), &
1350 0 : tcwat = tcwat(:ncol,:), &
1351 0 : lcwat = lcwat(:ncol,:), &
1352 : errmsg = errmsg, &
1353 70392 : errflg = errflg)
1354 :
1355 : ! Cloud water and ice particle sizes, saved in physics buffer for radiation
1356 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
1357 31183656 : rel(:,:) = 0._r8
1358 31183656 : rei(:,:) = 0._r8
1359 : !REMOVECAM_END
1360 :
1361 : call rk_stratiform_cloud_optical_properties_run( &
1362 : ncol = ncol, &
1363 : pver = pver, &
1364 : tmelt = tmelt, &
1365 0 : landfrac = landfrac(:ncol), &
1366 0 : icefrac = icefrac(:ncol), &
1367 0 : snowh = snowh(:ncol), &
1368 0 : landm = landm(:ncol), &
1369 0 : t = state1%t(:ncol,:), &
1370 0 : ps = state1%ps(:ncol), &
1371 0 : pmid = state1%pmid(:ncol,:), & ! below output:
1372 0 : rel = rel(:ncol,:), &
1373 0 : rei = rei(:ncol,:), &
1374 : errmsg = errmsg, &
1375 70392 : errflg = errflg)
1376 :
1377 70392 : call outfld('REL', rel, pcols, lchnk)
1378 70392 : call outfld('REI', rei, pcols, lchnk)
1379 :
1380 70392 : call physics_state_dealloc(state1)
1381 :
1382 211176 : end subroutine rk_stratiform_cam_tend
1383 :
1384 : !============================================================================ !
1385 : ! !
1386 : !============================================================================ !
1387 :
1388 0 : subroutine debug_microphys_1(state1,ptend,i,k, &
1389 : dtime,qme,fice,snow_pcw,prec_pcw, &
1390 : prain,nevapr,prodsnow, evapsnow, &
1391 : ice2pr,liq2pr,liq2snow)
1392 :
1393 70392 : use physics_types, only: physics_state, physics_ptend
1394 : use physconst, only: tmelt
1395 :
1396 : implicit none
1397 :
1398 : integer, intent(in) :: i,k
1399 : type(physics_state), intent(in) :: state1 ! local copy of the state variable
1400 : type(physics_ptend), intent(in) :: ptend ! local copy of the ptend variable
1401 : real(r8), intent(in) :: dtime ! timestep
1402 : real(r8), intent(in) :: qme(pcols,pver) ! local condensation - evaporation of cloud water
1403 :
1404 : real(r8), intent(in) :: prain(pcols,pver) ! local production of precipitation
1405 : real(r8), intent(in) :: nevapr(pcols,pver) ! local evaporation of precipitation
1406 : real(r8), intent(in) :: prodsnow(pcols,pver) ! local production of snow
1407 : real(r8), intent(in) :: evapsnow(pcols,pver) ! local evaporation of snow
1408 : real(r8), intent(in) :: ice2pr(pcols,pver) ! rate of conversion of ice to precip
1409 : real(r8), intent(in) :: liq2pr(pcols,pver) ! rate of conversion of liquid to precip
1410 : real(r8), intent(in) :: liq2snow(pcols,pver) ! rate of conversion of liquid to snow
1411 : real(r8), intent(in) :: fice (pcols,pver) ! Fractional ice content within cloud
1412 : real(r8), intent(in) :: snow_pcw(pcols)
1413 : real(r8), intent(in) :: prec_pcw(pcols)
1414 :
1415 : real(r8) hs1, qv1, ql1, qi1, qs1, qr1, fice2, pr1, w1, w2, w3, fliq, res
1416 : real(r8) w4, wl, wv, wi, wlf, wvf, wif, qif, qlf, qvf
1417 :
1418 0 : pr1 = 0
1419 0 : hs1 = 0
1420 0 : qv1 = 0
1421 0 : ql1 = 0
1422 0 : qi1 = 0
1423 0 : qs1 = 0
1424 0 : qr1 = 0
1425 0 : w1 = 0
1426 0 : wl = 0
1427 0 : wv = 0
1428 0 : wi = 0
1429 0 : wlf = 0
1430 0 : wvf = 0
1431 0 : wif = 0
1432 :
1433 :
1434 0 : write(iulog,*)
1435 0 : write(iulog,*) ' input state, t, q, l, i ', k, state1%t(i,k), state1%q(i,k,1), state1%q(i,k,ixcldliq), state1%q(i,k,ixcldice)
1436 0 : write(iulog,*) ' rain, snow, total from components before accumulation ', qr1, qs1, qr1+qs1
1437 0 : write(iulog,*) ' total precip before accumulation ', k, pr1
1438 :
1439 0 : wv = wv + state1%q(i,k,1 )*state1%pdel(i,k)/gravit
1440 0 : wl = wl + state1%q(i,k,ixcldliq)*state1%pdel(i,k)/gravit
1441 0 : wi = wi + state1%q(i,k,ixcldice)*state1%pdel(i,k)/gravit
1442 :
1443 0 : qvf = state1%q(i,k,1) + ptend%q(i,k,1)*dtime
1444 0 : qlf = state1%q(i,k,ixcldliq) + ptend%q(i,k,ixcldliq)*dtime
1445 0 : qif = state1%q(i,k,ixcldice) + ptend%q(i,k,ixcldice)*dtime
1446 :
1447 0 : if (qvf.lt.0._r8) then
1448 0 : write(iulog,*) ' qvf is negative *******', qvf
1449 : endif
1450 0 : if (qlf.lt.0._r8) then
1451 0 : write(iulog,*) ' qlf is negative *******', qlf
1452 : endif
1453 0 : if (qif.lt.0._r8) then
1454 0 : write(iulog,*) ' qif is negative *******', qif
1455 : endif
1456 0 : write(iulog,*) ' qvf, qlf, qif ', qvf, qlf, qif
1457 :
1458 0 : wvf = wvf + qvf*state1%pdel(i,k)/gravit
1459 0 : wlf = wlf + qlf*state1%pdel(i,k)/gravit
1460 0 : wif = wif + qif*state1%pdel(i,k)/gravit
1461 :
1462 0 : hs1 = hs1 + ptend%s(i,k)*state1%pdel(i,k)/gravit
1463 0 : pr1 = pr1 + state1%pdel(i,k)/gravit*(prain(i,k)-nevapr(i,k))
1464 0 : qv1 = qv1 - (qme(i,k)-nevapr(i,k))*state1%pdel(i,k)/gravit ! vdot
1465 0 : w1 = w1 + (qme(i,k)-prain(i,k))*state1%pdel(i,k)/gravit ! cdot
1466 0 : qi1 = qi1 + ((qme(i,k))*fice(i,k) -ice2pr(i,k) )*state1%pdel(i,k)/gravit ! idot
1467 0 : ql1 = ql1 + ((qme(i,k))*(1._r8-fice(i,k))-liq2pr(i,k) )*state1%pdel(i,k)/gravit ! ldot
1468 :
1469 : qr1 = qr1 &
1470 0 : + ( liq2pr(i,k)-liq2snow(i,k) & ! production of rain
1471 0 : -(nevapr(i,k)-evapsnow(i,k)) & ! rain evaporation
1472 0 : )*state1%pdel(i,k)/gravit
1473 : qs1 = qs1 &
1474 0 : + ( ice2pr(i,k) + liq2snow(i,k) & ! production of snow.Note last term has phase change
1475 0 : -evapsnow(i,k) & ! snow evaporation
1476 0 : )*state1%pdel(i,k)/gravit
1477 :
1478 0 : if (state1%t(i,k).gt.tmelt) then
1479 0 : qr1 = qr1 + qs1
1480 0 : qs1 = 0._r8
1481 : endif
1482 0 : write(iulog,*) ' rain, snow, total after accumulation ', qr1, qs1, qr1+qs1
1483 0 : write(iulog,*) ' total precip after accumulation ', k, pr1
1484 0 : write(iulog,*)
1485 0 : write(iulog,*) ' layer prain, nevapr, pdel ', prain(i,k), nevapr(i,k), state1%pdel(i,k)
1486 0 : write(iulog,*) ' layer prodsnow, ice2pr+liq2snow ', prodsnow(i,k), ice2pr(i,k)+liq2snow(i,k)
1487 0 : write(iulog,*) ' layer prain-prodsnow, liq2pr-liq2snow ', prain(i,k)-prodsnow(i,k), liq2pr(i,k)-liq2snow(i,k)
1488 0 : write(iulog,*) ' layer evapsnow, evaprain ', k, evapsnow(i,k), nevapr(i,k)-evapsnow(i,k)
1489 0 : write(iulog,*) ' layer ice2pr, liq2pr, liq2snow ', ice2pr(i,k), liq2pr(i,k), liq2snow(i,k)
1490 0 : write(iulog,*) ' layer ice2pr+liq2pr, prain ', ice2pr(i,k)+liq2pr(i,k), prain(i,k)
1491 0 : write(iulog,*)
1492 0 : write(iulog,*) ' qv1 vapor removed from col after accum (vdot) ', k, qv1
1493 0 : write(iulog,*) ' - (precip produced - vapor removed) after accum ', k, -pr1-qv1
1494 0 : write(iulog,*) ' condensate produce after accum ', k, w1
1495 0 : write(iulog,*) ' liq+ice tends accum ', k, ql1+qi1
1496 0 : write(iulog,*) ' change in total water after accum ', k, qv1+ql1+qi1
1497 0 : write(iulog,*) ' imbalance in colum after accum ', k, qs1+qr1+qv1+ql1+qi1
1498 0 : write(iulog,*) ' fice at this lev ', fice(i,k)
1499 0 : write(iulog,*)
1500 :
1501 0 : res = abs((qs1+qr1+qv1+ql1+qi1)/max(abs(qv1),abs(ql1),abs(qi1),abs(qs1),abs(qr1),1.e-36_r8))
1502 0 : write(iulog,*) ' relative residual in column method 1 ', k, res
1503 :
1504 0 : write(iulog,*) ' relative residual in column method 2 ',&
1505 0 : k, abs((qs1+qr1+qv1+ql1+qi1)/max(abs(qv1+ql1+qi1),1.e-36_r8))
1506 : ! if (abs((qs1+qr1+qv1+ql1+qi1)/(qs1+qr1+1.e-36)).gt.1.e-14) then
1507 0 : if (res.gt.1.e-14_r8) then
1508 0 : call endrun ('STRATIFORM_TEND')
1509 : endif
1510 :
1511 : ! w3 = qme(i,k) * (latvap + latice*fice(i,k)) &
1512 : ! + evapheat(i,k) + prfzheat(i,k) + meltheat(i,k)
1513 :
1514 0 : res = qs1+qr1-pr1
1515 0 : w4 = max(abs(qs1),abs(qr1),abs(pr1))
1516 0 : if (w4.gt.0._r8) then
1517 0 : if (res/w4.gt.1.e-14_r8) then
1518 0 : write(iulog,*) ' imbalance in precips calculated two ways '
1519 0 : write(iulog,*) ' res/w4, pr1, qr1, qs1, qr1+qs1 ', &
1520 0 : res/w4, pr1, qr1, qs1, qr1+qs1
1521 : ! call endrun()
1522 : endif
1523 : endif
1524 0 : if (k.eq.pver) then
1525 0 : write(iulog,*) ' pcond returned precip, rain and snow rates ', prec_pcw(i), prec_pcw(i)-snow_pcw(i), snow_pcw(i)
1526 0 : write(iulog,*) ' I calculate ', pr1, qr1, qs1
1527 : ! call endrun
1528 0 : write(iulog,*) ' byrons water check ', wv+wl+wi-pr1*dtime, wvf+wlf+wif
1529 : endif
1530 0 : write(iulog,*)
1531 :
1532 :
1533 0 : end subroutine debug_microphys_1
1534 :
1535 : !============================================================================ !
1536 : ! !
1537 : !============================================================================ !
1538 :
1539 0 : subroutine debug_microphys_2(state1,&
1540 : snow_pcw,fsaut,fsacw ,fsaci, meltheat)
1541 :
1542 0 : use ppgrid, only: pver
1543 : use physconst, only: tmelt
1544 : use physics_types, only: physics_state
1545 :
1546 : implicit none
1547 :
1548 : type(physics_state), intent(in) :: state1 ! local copy of the state variable
1549 : real(r8), intent(in) :: snow_pcw(pcols)
1550 : real(r8), intent(in) :: fsaut(pcols,pver)
1551 : real(r8), intent(in) :: fsacw(pcols,pver)
1552 : real(r8), intent(in) :: fsaci(pcols,pver)
1553 : real(r8), intent(in) :: meltheat(pcols,pver) ! heating rate due to phase change of precip
1554 :
1555 :
1556 : integer i,ncol,lchnk
1557 :
1558 :
1559 0 : ncol = state1%ncol
1560 0 : lchnk = state1%lchnk
1561 :
1562 0 : do i = 1,ncol
1563 0 : if (snow_pcw(i) .gt. 0.01_r8/8.64e4_r8 .and. state1%t(i,pver) .gt. tmelt) then
1564 0 : write(iulog,*) ' stratiform: snow, temp, ', i, lchnk, &
1565 0 : snow_pcw(i), state1%t(i,pver)
1566 0 : write(iulog,*) ' t ', state1%t(i,:)
1567 0 : write(iulog,*) ' fsaut ', fsaut(i,:)
1568 0 : write(iulog,*) ' fsacw ', fsacw(i,:)
1569 0 : write(iulog,*) ' fsaci ', fsaci(i,:)
1570 0 : write(iulog,*) ' meltheat ', meltheat(i,:)
1571 0 : call endrun ('STRATIFORM_TEND')
1572 : endif
1573 :
1574 0 : if (snow_pcw(i)*8.64e4_r8 .lt. -1.e-5_r8) then
1575 0 : write(iulog,*) ' neg snow ', snow_pcw(i)*8.64e4_r8
1576 0 : write(iulog,*) ' stratiform: snow_pcw, temp, ', i, lchnk, &
1577 0 : snow_pcw(i), state1%t(i,pver)
1578 0 : write(iulog,*) ' t ', state1%t(i,:)
1579 0 : write(iulog,*) ' fsaut ', fsaut(i,:)
1580 0 : write(iulog,*) ' fsacw ', fsacw(i,:)
1581 0 : write(iulog,*) ' fsaci ', fsaci(i,:)
1582 0 : write(iulog,*) ' meltheat ', meltheat(i,:)
1583 0 : call endrun ('STRATIFORM_TEND')
1584 : endif
1585 : end do
1586 :
1587 0 : end subroutine debug_microphys_2
1588 :
1589 : end module rk_stratiform_cam
|