Line data Source code
1 : module conv_water
2 :
3 : ! --------------------------------------------------------------------- !
4 : ! Purpose: !
5 : ! Computes grid-box average liquid (and ice) from stratus and cumulus !
6 : ! These values used by both the radiation and the COSP diagnostics. !
7 : ! !
8 : ! Method: !
9 : ! Extract information about deep+shallow liquid and cloud fraction from !
10 : ! the physics buffer. !
11 : ! !
12 : ! Author: Rich Neale, August 2006 !
13 : ! October 2006: Allow averaging of liquid to give a linear !
14 : ! average in emissivity. !
15 : ! Andrew Gettelman October 2010 Separate module !
16 : !---------------------------------------------------------------------- !
17 :
18 : use shr_kind_mod, only: r8=>shr_kind_r8
19 : use spmd_utils, only: masterproc
20 : use ppgrid, only: pcols, pver, pverp
21 : use physconst, only: gravit, latvap, latice
22 : use cam_abortutils, only: endrun
23 :
24 : use perf_mod
25 : use cam_logfile, only: iulog
26 :
27 : implicit none
28 : private
29 : save
30 :
31 : public :: &
32 : conv_water_readnl, &
33 : conv_water_register, &
34 : conv_water_init, &
35 : conv_water_4rad, &
36 : conv_water_in_rad
37 :
38 : ! pbuf indices
39 :
40 : integer :: icwmrsh_idx, icwmrdp_idx, fice_idx, sh_frac_idx, dp_frac_idx, &
41 : ast_idx, rei_idx
42 :
43 : integer :: ixcldice, ixcldliq
44 : integer :: gb_totcldliqmr_idx, gb_totcldicemr_idx
45 :
46 : ! Namelist
47 : integer, parameter :: unset_int = huge(1)
48 :
49 : integer :: conv_water_in_rad = unset_int ! 0==> No; 1==> Yes-Arithmetic average;
50 : ! 2==> Yes-Average in emissivity.
51 : integer :: conv_water_mode
52 : real(r8) :: frac_limit
53 :
54 : !=============================================================================================
55 : contains
56 : !=============================================================================================
57 :
58 1490712 : subroutine conv_water_readnl(nlfile)
59 :
60 : use namelist_utils, only: find_group_name
61 : use units, only: getunit, freeunit
62 : use mpishorthand
63 :
64 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
65 :
66 : ! Local variables
67 : integer :: unitn, ierr
68 : character(len=*), parameter :: subname = 'conv_water_readnl'
69 :
70 : real(r8) :: conv_water_frac_limit = 0._r8
71 :
72 : namelist /conv_water_nl/ conv_water_in_rad, conv_water_frac_limit
73 : !-----------------------------------------------------------------------------
74 :
75 1536 : if (masterproc) then
76 2 : unitn = getunit()
77 2 : open( unitn, file=trim(nlfile), status='old' )
78 2 : call find_group_name(unitn, 'conv_water_nl', status=ierr)
79 2 : if (ierr == 0) then
80 2 : read(unitn, conv_water_nl, iostat=ierr)
81 2 : if (ierr /= 0) then
82 0 : call endrun(subname // ':: ERROR reading namelist')
83 : end if
84 : end if
85 2 : close(unitn)
86 2 : call freeunit(unitn)
87 : end if
88 :
89 : #ifdef SPMD
90 : ! Broadcast namelist variables
91 1536 : call mpibcast(conv_water_in_rad, 1, mpiint, 0, mpicom)
92 1536 : call mpibcast(conv_water_frac_limit, 1, mpir8, 0, mpicom)
93 : #endif
94 :
95 1536 : conv_water_mode = conv_water_in_rad
96 1536 : frac_limit = conv_water_frac_limit
97 :
98 1536 : end subroutine conv_water_readnl
99 :
100 : !=============================================================================================
101 :
102 1536 : subroutine conv_water_register
103 :
104 : !---------------------------------------------------------------------- !
105 : ! !
106 : ! Register the fields in the physics buffer. !
107 : ! !
108 : !---------------------------------------------------------------------- !
109 :
110 : use constituents, only: cnst_add, pcnst
111 : use physconst, only: mwdry, cpair
112 :
113 : use physics_buffer, only : pbuf_add_field, dtype_r8
114 :
115 : !-----------------------------------------------------------------------
116 :
117 : ! grid box total cloud liquid water mixing ratio (kg/kg)
118 1536 : call pbuf_add_field('GB_TOTCLDLIQMR', 'physpkg', dtype_r8, (/pcols,pver/), gb_totcldliqmr_idx)
119 : ! grid box total cloud ice water mixing ratio (kg/kg)
120 1536 : call pbuf_add_field('GB_TOTCLDICEMR', 'physpkg', dtype_r8, (/pcols,pver/), gb_totcldicemr_idx)
121 :
122 1536 : end subroutine conv_water_register
123 :
124 :
125 : !============================================================================ !
126 : ! !
127 : !============================================================================ !
128 :
129 1536 : subroutine conv_water_init()
130 : ! --------------------------------------------------------------------- !
131 : ! Purpose: !
132 : ! Initializes the pbuf indices required by conv_water
133 : ! --------------------------------------------------------------------- !
134 :
135 :
136 1536 : use physics_buffer, only : pbuf_get_index
137 : use cam_history, only : addfld
138 :
139 : use constituents, only: cnst_get_ind
140 :
141 : implicit none
142 :
143 1536 : call cnst_get_ind('CLDICE', ixcldice)
144 1536 : call cnst_get_ind('CLDLIQ', ixcldliq)
145 :
146 1536 : icwmrsh_idx = pbuf_get_index('ICWMRSH')
147 1536 : icwmrdp_idx = pbuf_get_index('ICWMRDP')
148 1536 : fice_idx = pbuf_get_index('FICE')
149 1536 : sh_frac_idx = pbuf_get_index('SH_FRAC')
150 1536 : dp_frac_idx = pbuf_get_index('DP_FRAC')
151 1536 : ast_idx = pbuf_get_index('AST')
152 1536 : rei_idx = pbuf_get_index('REI')
153 :
154 : ! Convective cloud water variables.
155 3072 : call addfld ('ICIMRCU', (/ 'lev' /), 'A', 'kg/kg', 'Convection in-cloud ice mixing ratio ' )
156 3072 : call addfld ('ICLMRCU', (/ 'lev' /), 'A', 'kg/kg', 'Convection in-cloud liquid mixing ratio ')
157 3072 : call addfld ('ICIMRTOT', (/ 'lev' /), 'A', 'kg/kg', 'Total in-cloud ice mixing ratio ' )
158 3072 : call addfld ('ICLMRTOT', (/ 'lev' /), 'A', 'kg/kg', 'Total in-cloud liquid mixing ratio ' )
159 :
160 3072 : call addfld ('GCLMRDP', (/ 'lev' /), 'A', 'kg/kg', 'Grid-mean deep convective LWC' )
161 3072 : call addfld ('GCIMRDP', (/ 'lev' /), 'A', 'kg/kg', 'Grid-mean deep convective IWC' )
162 3072 : call addfld ('GCLMRSH', (/ 'lev' /), 'A', 'kg/kg', 'Grid-mean shallow convective LWC' )
163 3072 : call addfld ('GCIMRSH', (/ 'lev' /), 'A', 'kg/kg', 'Grid-mean shallow convective IWC' )
164 3072 : call addfld ('FRESH', (/ 'lev' /), 'A', '1', 'Fractional occurrence of shallow cumulus with condensate')
165 3072 : call addfld ('FREDP', (/ 'lev' /), 'A', '1', 'Fractional occurrence of deep cumulus with condensate')
166 3072 : call addfld ('FRECU', (/ 'lev' /), 'A', '1', 'Fractional occurrence of cumulus with condensate')
167 3072 : call addfld ('FRETOT', (/ 'lev' /), 'A', '1', 'Fractional occurrence of cloud with condensate')
168 :
169 1536 : end subroutine conv_water_init
170 :
171 1489176 : subroutine conv_water_4rad(state, pbuf)
172 :
173 : ! --------------------------------------------------------------------- !
174 : ! Purpose: !
175 : ! Computes grid-box average liquid (and ice) from stratus and cumulus !
176 : ! Just for the purposes of radiation. !
177 : ! !
178 : ! Method: !
179 : ! Extract information about deep+shallow liquid and cloud fraction from !
180 : ! the physics buffer. !
181 : ! !
182 : ! Author: Rich Neale, August 2006 !
183 : ! October 2006: Allow averaging of liquid to give a linear !
184 : ! average in emissivity. !
185 : ! !
186 : !---------------------------------------------------------------------- !
187 :
188 :
189 1536 : use physics_buffer, only : physics_buffer_desc, pbuf_get_field, pbuf_old_tim_idx
190 :
191 : use physics_types, only: physics_state
192 : use cam_history, only: outfld
193 : use phys_control, only: phys_getopts
194 :
195 : implicit none
196 :
197 : ! ---------------------- !
198 : ! Input-Output Arguments !
199 : ! ---------------------- !
200 :
201 :
202 : type(physics_state), target, intent(in) :: state ! state variables
203 : type(physics_buffer_desc), pointer :: pbuf(:)
204 :
205 : ! --------------- !
206 : ! Local Workspace !
207 : ! --------------- !
208 :
209 1489176 : real(r8), pointer, dimension(:,:) :: pdel ! Moist pressure difference across layer
210 1489176 : real(r8), pointer, dimension(:,:) :: ls_liq ! Large-scale contributions to GBA cloud liq
211 1489176 : real(r8), pointer, dimension(:,:) :: ls_ice ! Large-scale contributions to GBA cloud ice
212 :
213 : ! Physics buffer fields
214 1489176 : real(r8), pointer, dimension(:,:) :: ast ! Physical liquid+ice stratus cloud fraction
215 1489176 : real(r8), pointer, dimension(:,:) :: sh_frac ! Shallow convective cloud fraction
216 1489176 : real(r8), pointer, dimension(:,:) :: dp_frac ! Deep convective cloud fraction
217 1489176 : real(r8), pointer, dimension(:,:) :: rei ! Ice effective drop size (microns)
218 :
219 1489176 : real(r8), pointer, dimension(:,:) :: dp_icwmr ! Deep conv. cloud water
220 1489176 : real(r8), pointer, dimension(:,:) :: sh_icwmr ! Shallow conv. cloud water
221 1489176 : real(r8), pointer, dimension(:,:) :: fice ! Ice partitioning ratio
222 :
223 1489176 : real(r8), pointer, dimension(:,:) :: totg_ice ! Grid box total cloud ice mixing ratio
224 1489176 : real(r8), pointer, dimension(:,:) :: totg_liq ! Grid box total cloud liquid mixing ratio
225 :
226 : real(r8) :: conv_ice(pcols,pver) ! Convective contributions to IC cloud ice
227 : real(r8) :: conv_liq(pcols,pver) ! Convective contributions to IC cloud liquid
228 : real(r8) :: tot_ice(pcols,pver) ! Total IC ice
229 : real(r8) :: tot_liq(pcols,pver) ! Total IC liquid
230 :
231 : integer :: i,k,itim_old ! Lon, lev indices buff stuff.
232 : real(r8) :: cu_icwmr ! Convective water for this grid-box.
233 : real(r8) :: ls_icwmr ! Large-scale water for this grid-box.
234 : real(r8) :: tot_icwmr ! Large-scale water for this grid-box.
235 : real(r8) :: ls_frac ! Large-scale cloud frac for this grid-box.
236 : real(r8) :: tot0_frac, cu0_frac, dp0_frac, sh0_frac
237 : real(r8) :: kabs, kabsi, kabsl, alpha, dp0, sh0, ic_limit
238 : real(r8) :: wrk1
239 :
240 : real(r8) :: totg_ice_sh(pcols,pver) ! Grid-mean IWP from shallow convective cloud
241 : real(r8) :: totg_liq_sh(pcols,pver) ! Grid-mean LWP from shallow convective cloud
242 : real(r8) :: totg_ice_dp(pcols,pver) ! Grid-mean IWP from deep convective cloud
243 : real(r8) :: totg_liq_dp(pcols,pver) ! Grid-mean LWP from deep convective cloud
244 : real(r8) :: fresh(pcols,pver) ! Fractional occurrence of shallow cumulus
245 : real(r8) :: fredp(pcols,pver) ! Fractional occurrence of deep cumulus
246 : real(r8) :: frecu(pcols,pver) ! Fractional occurrence of cumulus
247 : real(r8) :: fretot(pcols,pver) ! Fractional occurrence of cloud
248 :
249 : integer :: lchnk
250 : integer :: ncol
251 :
252 : ! --------- !
253 : ! Parameter !
254 : ! --------- !
255 :
256 : parameter( kabsl = 0.090361_r8, ic_limit = 1.e-12_r8 )
257 : character(len=16) :: microp_scheme
258 :
259 1489176 : ncol = state%ncol
260 1489176 : lchnk = state%lchnk
261 1489176 : pdel => state%pdel
262 1489176 : ls_liq => state%q(:,:,ixcldliq)
263 1489176 : ls_ice => state%q(:,:,ixcldice)
264 :
265 : ! Get microphysics option
266 1489176 : call phys_getopts( microp_scheme_out = microp_scheme )
267 :
268 : ! Get convective in-cloud water and ice/water temperature partitioning.
269 :
270 1489176 : call pbuf_get_field(pbuf, icwmrsh_idx, sh_icwmr )
271 1489176 : call pbuf_get_field(pbuf, icwmrdp_idx, dp_icwmr )
272 1489176 : call pbuf_get_field(pbuf, fice_idx, fice )
273 :
274 : ! Get convective in-cloud fraction
275 :
276 1489176 : call pbuf_get_field(pbuf, sh_frac_idx, sh_frac )
277 1489176 : call pbuf_get_field(pbuf, dp_frac_idx, dp_frac )
278 1489176 : call pbuf_get_field(pbuf, rei_idx, rei )
279 :
280 1489176 : itim_old = pbuf_old_tim_idx()
281 5956704 : call pbuf_get_field(pbuf, ast_idx, ast, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
282 :
283 : ! Fields computed below and stored in pbuf.
284 1489176 : call pbuf_get_field(pbuf, gb_totcldicemr_idx, totg_ice)
285 1489176 : call pbuf_get_field(pbuf, gb_totcldliqmr_idx, totg_liq)
286 :
287 : ! --------------------------------------------------------------- !
288 : ! Loop through grid-boxes and determine: !
289 : ! 1. Effective mean in-cloud convective ice/liquid (deep+shallow) !
290 : ! 2. Effective mean in-cloud total ice/liquid (ls+convective) !
291 : ! --------------------------------------------------------------- !
292 :
293 1489176 : fresh(:,:) = 0._r8
294 1489176 : fredp(:,:) = 0._r8
295 1489176 : frecu(:,:) = 0._r8
296 1489176 : fretot(:,:) = 0._r8
297 :
298 139982544 : do k = 1, pver
299 2314006344 : do i = 1, ncol
300 :
301 2174023800 : if( sh_frac(i,k) <= frac_limit .or. sh_icwmr(i,k) <= ic_limit ) then
302 : sh0_frac = 0._r8
303 : else
304 0 : sh0_frac = sh_frac(i,k)
305 : endif
306 2174023800 : if( dp_frac(i,k) <= frac_limit .or. dp_icwmr(i,k) <= ic_limit ) then
307 : dp0_frac = 0._r8
308 : else
309 55309321 : dp0_frac = dp_frac(i,k)
310 : endif
311 2174023800 : cu0_frac = sh0_frac + dp0_frac
312 :
313 : ! For the moment calculate the emissivity based upon the ls clouds ice fraction
314 :
315 2174023800 : wrk1 = min(1._r8,max(0._r8, ls_ice(i,k)/(ls_ice(i,k)+ls_liq(i,k)+1.e-36_r8)))
316 :
317 2174023800 : if( ( cu0_frac < frac_limit ) .or. ( ( sh_icwmr(i,k) + dp_icwmr(i,k) ) < ic_limit ) ) then
318 :
319 2118714479 : cu0_frac = 0._r8
320 2118714479 : cu_icwmr = 0._r8
321 :
322 2118714479 : ls_frac = ast(i,k)
323 2118714479 : if( ls_frac < frac_limit ) then
324 : ls_frac = 0._r8
325 : ls_icwmr = 0._r8
326 : else
327 231705349 : ls_icwmr = ( ls_liq(i,k) + ls_ice(i,k) )/max(frac_limit,ls_frac) ! Convert to IC value.
328 : end if
329 :
330 : tot0_frac = ls_frac
331 : tot_icwmr = ls_icwmr
332 :
333 : else
334 :
335 : ! Select radiation constants (effective radii) for emissivity averaging.
336 :
337 55309321 : if( microp_scheme == 'RK' .or. microp_scheme == 'SPCAM_sam1mom') then
338 0 : kabsi = 0.005_r8 + 1._r8/rei(i,k)
339 : else
340 55309321 : kabsi = 0.005_r8 + 1._r8/min(max(13._r8,rei(i,k)),130._r8)
341 : endif
342 55309321 : kabs = kabsl * ( 1._r8 - wrk1 ) + kabsi * wrk1
343 55309321 : alpha = -1.66_r8*kabs*pdel(i,k)/gravit*1000.0_r8
344 :
345 : ! Selecting cumulus in-cloud water.
346 :
347 : select case (conv_water_mode) ! Type of average
348 : case (1) ! Area weighted arithmetic average
349 55309321 : cu_icwmr = ( sh0_frac * sh_icwmr(i,k) + dp0_frac*dp_icwmr(i,k))/max(frac_limit,cu0_frac)
350 : case (2)
351 0 : sh0 = exp(alpha*sh_icwmr(i,k))
352 0 : dp0 = exp(alpha*dp_icwmr(i,k))
353 0 : cu_icwmr = log((sh0_frac*sh0+dp0_frac*dp0)/max(frac_limit,cu0_frac))
354 55309321 : cu_icwmr = cu_icwmr/alpha
355 : case default ! Area weighted 'arithmetic in emissivity' average.
356 : ! call endrun ('CONV_WATER_4_RAD: Unknown option for conv_water_in_rad - exiting')
357 : end select
358 :
359 : ! Selecting total in-cloud water.
360 : ! Attribute large-scale/convective area fraction differently from default.
361 :
362 55309321 : ls_frac = ast(i,k)
363 55309321 : ls_icwmr = (ls_liq(i,k) + ls_ice(i,k))/max(frac_limit,ls_frac) ! Convert to IC value.
364 55309321 : tot0_frac = (ls_frac + cu0_frac)
365 :
366 55309321 : select case (conv_water_mode) ! Type of average
367 : case (1) ! Area weighted 'arithmetic in emissivity' average
368 55309321 : tot_icwmr = (ls_frac*ls_icwmr + cu0_frac*cu_icwmr)/max(frac_limit,tot0_frac)
369 : case (2)
370 0 : tot_icwmr = log((ls_frac*exp(alpha*ls_icwmr)+cu0_frac*exp(alpha*cu_icwmr))/max(frac_limit,tot0_frac))
371 55309321 : tot_icwmr = tot_icwmr/alpha
372 : case default ! Area weighted 'arithmetic in emissivity' average.
373 : ! call endrun ('CONV_WATER_4_RAD: Unknown option for conv_water_in_rad - exiting')
374 : end select
375 :
376 : end if
377 :
378 : ! Repartition convective cloud water into liquid and ice phase.
379 : ! Currently, this partition is made using the ice fraction of stratus condensate.
380 : ! In future, we should use ice fraction explicitly computed from the convection scheme.
381 :
382 2174023800 : conv_ice(i,k) = cu_icwmr * wrk1
383 2174023800 : conv_liq(i,k) = cu_icwmr * (1._r8-wrk1)
384 :
385 2174023800 : tot_ice(i,k) = tot_icwmr * wrk1
386 2174023800 : tot_liq(i,k) = tot_icwmr * (1._r8-wrk1)
387 :
388 2174023800 : totg_ice(i,k) = tot0_frac * tot_icwmr * wrk1
389 2174023800 : totg_liq(i,k) = tot0_frac * tot_icwmr * (1._r8-wrk1)
390 :
391 : ! Grid-mean convective water
392 2174023800 : totg_ice_sh(i,k) = sh0_frac * sh_icwmr(i,k) * wrk1
393 2174023800 : totg_ice_dp(i,k) = dp0_frac * dp_icwmr(i,k) * wrk1
394 2174023800 : totg_liq_sh(i,k) = sh0_frac * sh_icwmr(i,k) * (1._r8-wrk1)
395 2174023800 : totg_liq_dp(i,k) = dp0_frac * dp_icwmr(i,k) * (1._r8-wrk1)
396 2174023800 : if( sh0_frac > frac_limit ) then
397 0 : fresh(i,k) = 1._r8
398 : endif
399 2174023800 : if( dp0_frac > frac_limit ) then
400 55309321 : fredp(i,k) = 1._r8
401 : endif
402 2174023800 : if( cu0_frac > frac_limit ) then
403 55309321 : frecu(i,k) = 1._r8
404 : endif
405 2312517168 : if( tot0_frac > frac_limit ) then
406 287014670 : fretot(i,k) = 1._r8
407 : endif
408 :
409 : end do
410 : end do
411 :
412 : ! Output convective IC WMRs
413 :
414 1489176 : call outfld( 'ICLMRCU ', conv_liq , pcols, lchnk )
415 1489176 : call outfld( 'ICIMRCU ', conv_ice , pcols, lchnk )
416 1489176 : call outfld( 'ICLMRTOT', tot_liq , pcols, lchnk )
417 1489176 : call outfld( 'ICIMRTOT', tot_ice , pcols, lchnk )
418 :
419 1489176 : call outfld('GCLMRDP', totg_liq_dp, pcols, lchnk)
420 1489176 : call outfld('GCIMRDP', totg_ice_dp, pcols, lchnk)
421 1489176 : call outfld('GCLMRSH', totg_liq_sh, pcols, lchnk)
422 1489176 : call outfld('GCIMRSH', totg_ice_sh, pcols, lchnk)
423 1489176 : call outfld('FRESH', fresh, pcols, lchnk)
424 1489176 : call outfld('FREDP', fredp, pcols, lchnk)
425 1489176 : call outfld('FRECU', frecu, pcols, lchnk)
426 1489176 : call outfld('FRETOT', fretot, pcols, lchnk)
427 :
428 2978352 : end subroutine conv_water_4rad
429 :
430 : end module conv_water
|