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 : ! Just for the purposes of radiation. !
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, sh_cldliq1_idx, sh_cldice1_idx, rei_idx
42 :
43 : integer :: ixcldice, ixcldliq
44 :
45 : ! Namelist
46 : integer, parameter :: unset_int = huge(1)
47 :
48 : integer :: conv_water_in_rad = unset_int ! 0==> No; 1==> Yes-Arithmetic average;
49 : ! 2==> Yes-Average in emissivity.
50 : integer :: conv_water_mode
51 : real(r8) :: frac_limit
52 :
53 : !=============================================================================================
54 : contains
55 : !=============================================================================================
56 :
57 60360 : subroutine conv_water_readnl(nlfile)
58 :
59 : use namelist_utils, only: find_group_name
60 : use units, only: getunit, freeunit
61 : use mpishorthand
62 :
63 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
64 :
65 : ! Local variables
66 : integer :: unitn, ierr
67 : character(len=*), parameter :: subname = 'conv_water_readnl'
68 :
69 : real(r8) :: conv_water_frac_limit = 0._r8
70 :
71 : namelist /conv_water_nl/ conv_water_in_rad, conv_water_frac_limit
72 : !-----------------------------------------------------------------------------
73 :
74 1536 : if (masterproc) then
75 2 : unitn = getunit()
76 2 : open( unitn, file=trim(nlfile), status='old' )
77 2 : call find_group_name(unitn, 'conv_water_nl', status=ierr)
78 2 : if (ierr == 0) then
79 2 : read(unitn, conv_water_nl, iostat=ierr)
80 2 : if (ierr /= 0) then
81 0 : call endrun(subname // ':: ERROR reading namelist')
82 : end if
83 : end if
84 2 : close(unitn)
85 2 : call freeunit(unitn)
86 : end if
87 :
88 : #ifdef SPMD
89 : ! Broadcast namelist variables
90 1536 : call mpibcast(conv_water_in_rad, 1, mpiint, 0, mpicom)
91 1536 : call mpibcast(conv_water_frac_limit, 1, mpir8, 0, mpicom)
92 : #endif
93 :
94 1536 : conv_water_mode = conv_water_in_rad
95 1536 : frac_limit = conv_water_frac_limit
96 :
97 1536 : end subroutine conv_water_readnl
98 :
99 : !=============================================================================================
100 :
101 1536 : subroutine conv_water_register
102 :
103 : !---------------------------------------------------------------------- !
104 : ! !
105 : ! Register the fields in the physics buffer. !
106 : ! !
107 : !---------------------------------------------------------------------- !
108 :
109 : use constituents, only: cnst_add, pcnst
110 : use physconst, only: mwdry, cpair
111 :
112 : use physics_buffer, only : pbuf_add_field, dtype_r8
113 :
114 : !-----------------------------------------------------------------------
115 :
116 : ! these calls were already done in convect_shallow...so here I add the same fields to the physics buffer with a "1" at the end
117 : ! shallow gbm cloud liquid water (kg/kg)
118 1536 : call pbuf_add_field('SH_CLDLIQ1','physpkg',dtype_r8,(/pcols,pver/),sh_cldliq1_idx)
119 : ! shallow gbm cloud ice water (kg/kg)
120 1536 : call pbuf_add_field('SH_CLDICE1','physpkg',dtype_r8,(/pcols,pver/),sh_cldice1_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 58824 : subroutine conv_water_4rad(state, pbuf, totg_liq, totg_ice)
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 : real(r8), intent(out):: totg_ice(pcols,pver) ! Total GBA in-cloud ice
206 : real(r8), intent(out):: totg_liq(pcols,pver) ! Total GBA in-cloud liquid
207 :
208 : ! --------------- !
209 : ! Local Workspace !
210 : ! --------------- !
211 :
212 58824 : real(r8), pointer, dimension(:,:) :: pdel ! Moist pressure difference across layer
213 58824 : real(r8), pointer, dimension(:,:) :: ls_liq ! Large-scale contributions to GBA cloud liq
214 58824 : real(r8), pointer, dimension(:,:) :: ls_ice ! Large-scale contributions to GBA cloud ice
215 :
216 : ! Physics buffer fields
217 58824 : real(r8), pointer, dimension(:,:) :: ast ! Physical liquid+ice stratus cloud fraction
218 58824 : real(r8), pointer, dimension(:,:) :: sh_frac ! Shallow convective cloud fraction
219 58824 : real(r8), pointer, dimension(:,:) :: dp_frac ! Deep convective cloud fraction
220 58824 : real(r8), pointer, dimension(:,:) :: rei ! Ice effective drop size (microns)
221 :
222 58824 : real(r8), pointer, dimension(:,:) :: dp_icwmr ! Deep conv. cloud water
223 58824 : real(r8), pointer, dimension(:,:) :: sh_icwmr ! Shallow conv. cloud water
224 58824 : real(r8), pointer, dimension(:,:) :: fice ! Ice partitioning ratio
225 58824 : real(r8), pointer, dimension(:,:) :: sh_cldliq ! shallow convection gbx liq cld mixing ratio for COSP
226 58824 : real(r8), pointer, dimension(:,:) :: sh_cldice ! shallow convection gbx ice cld mixing ratio for COSP
227 :
228 : real(r8) :: conv_ice(pcols,pver) ! Convective contributions to IC cloud ice
229 : real(r8) :: conv_liq(pcols,pver) ! Convective contributions to IC cloud liquid
230 : real(r8) :: tot_ice(pcols,pver) ! Total IC ice
231 : real(r8) :: tot_liq(pcols,pver) ! Total IC liquid
232 :
233 : integer :: i,k,itim_old ! Lon, lev indices buff stuff.
234 : real(r8) :: cu_icwmr ! Convective water for this grid-box.
235 : real(r8) :: ls_icwmr ! Large-scale water for this grid-box.
236 : real(r8) :: tot_icwmr ! Large-scale water for this grid-box.
237 : real(r8) :: ls_frac ! Large-scale cloud frac for this grid-box.
238 : real(r8) :: tot0_frac, cu0_frac, dp0_frac, sh0_frac
239 : real(r8) :: kabs, kabsi, kabsl, alpha, dp0, sh0, ic_limit
240 : real(r8) :: wrk1
241 :
242 : real(r8) :: totg_ice_sh(pcols,pver) ! Grid-mean IWP from shallow convective cloud
243 : real(r8) :: totg_liq_sh(pcols,pver) ! Grid-mean LWP from shallow convective cloud
244 : real(r8) :: totg_ice_dp(pcols,pver) ! Grid-mean IWP from deep convective cloud
245 : real(r8) :: totg_liq_dp(pcols,pver) ! Grid-mean LWP from deep convective cloud
246 : real(r8) :: fresh(pcols,pver) ! Fractional occurrence of shallow cumulus
247 : real(r8) :: fredp(pcols,pver) ! Fractional occurrence of deep cumulus
248 : real(r8) :: frecu(pcols,pver) ! Fractional occurrence of cumulus
249 : real(r8) :: fretot(pcols,pver) ! Fractional occurrence of cloud
250 :
251 : integer :: lchnk
252 : integer :: ncol
253 :
254 : ! --------- !
255 : ! Parameter !
256 : ! --------- !
257 :
258 : parameter( kabsl = 0.090361_r8, ic_limit = 1.e-12_r8 )
259 : character(len=16) :: microp_scheme
260 :
261 58824 : ncol = state%ncol
262 58824 : lchnk = state%lchnk
263 58824 : pdel => state%pdel
264 58824 : ls_liq => state%q(:,:,ixcldliq)
265 58824 : ls_ice => state%q(:,:,ixcldice)
266 :
267 : ! Get microphysics option
268 58824 : call phys_getopts( microp_scheme_out = microp_scheme )
269 :
270 : ! Get convective in-cloud water and ice/water temperature partitioning.
271 :
272 58824 : call pbuf_get_field(pbuf, icwmrsh_idx, sh_icwmr )
273 58824 : call pbuf_get_field(pbuf, icwmrdp_idx, dp_icwmr )
274 58824 : call pbuf_get_field(pbuf, fice_idx, fice )
275 :
276 : ! Get convective in-cloud fraction
277 :
278 58824 : call pbuf_get_field(pbuf, sh_frac_idx, sh_frac )
279 58824 : call pbuf_get_field(pbuf, dp_frac_idx, dp_frac )
280 58824 : call pbuf_get_field(pbuf, rei_idx, rei )
281 :
282 58824 : itim_old = pbuf_old_tim_idx()
283 235296 : call pbuf_get_field(pbuf, ast_idx, ast, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
284 :
285 : ! --------------------------------------------------------------- !
286 : ! Loop through grid-boxes and determine: !
287 : ! 1. Effective mean in-cloud convective ice/liquid (deep+shallow) !
288 : ! 2. Effective mean in-cloud total ice/liquid (ls+convective) !
289 : ! --------------------------------------------------------------- !
290 :
291 58824 : fresh(:,:) = 0._r8
292 58824 : fredp(:,:) = 0._r8
293 58824 : frecu(:,:) = 0._r8
294 58824 : fretot(:,:) = 0._r8
295 :
296 5529456 : do k = 1, pver
297 91405656 : do i = 1, ncol
298 :
299 85876200 : if( sh_frac(i,k) <= frac_limit .or. sh_icwmr(i,k) <= ic_limit ) then
300 : sh0_frac = 0._r8
301 : else
302 0 : sh0_frac = sh_frac(i,k)
303 : endif
304 85876200 : if( dp_frac(i,k) <= frac_limit .or. dp_icwmr(i,k) <= ic_limit ) then
305 : dp0_frac = 0._r8
306 : else
307 1947302 : dp0_frac = dp_frac(i,k)
308 : endif
309 85876200 : cu0_frac = sh0_frac + dp0_frac
310 :
311 : ! For the moment calculate the emissivity based upon the ls clouds ice fraction
312 :
313 85876200 : wrk1 = min(1._r8,max(0._r8, ls_ice(i,k)/(ls_ice(i,k)+ls_liq(i,k)+1.e-36_r8)))
314 :
315 85876200 : if( ( cu0_frac < frac_limit ) .or. ( ( sh_icwmr(i,k) + dp_icwmr(i,k) ) < ic_limit ) ) then
316 :
317 83928898 : cu0_frac = 0._r8
318 83928898 : cu_icwmr = 0._r8
319 :
320 83928898 : ls_frac = ast(i,k)
321 83928898 : if( ls_frac < frac_limit ) then
322 : ls_frac = 0._r8
323 : ls_icwmr = 0._r8
324 : else
325 8858888 : ls_icwmr = ( ls_liq(i,k) + ls_ice(i,k) )/max(frac_limit,ls_frac) ! Convert to IC value.
326 : end if
327 :
328 : tot0_frac = ls_frac
329 : tot_icwmr = ls_icwmr
330 :
331 : else
332 :
333 : ! Select radiation constants (effective radii) for emissivity averaging.
334 :
335 1947302 : if( microp_scheme == 'RK' .or. microp_scheme == 'SPCAM_sam1mom') then
336 0 : kabsi = 0.005_r8 + 1._r8/rei(i,k)
337 : else
338 1947302 : kabsi = 0.005_r8 + 1._r8/min(max(13._r8,rei(i,k)),130._r8)
339 : endif
340 1947302 : kabs = kabsl * ( 1._r8 - wrk1 ) + kabsi * wrk1
341 1947302 : alpha = -1.66_r8*kabs*pdel(i,k)/gravit*1000.0_r8
342 :
343 : ! Selecting cumulus in-cloud water.
344 :
345 : select case (conv_water_mode) ! Type of average
346 : case (1) ! Area weighted arithmetic average
347 1947302 : cu_icwmr = ( sh0_frac * sh_icwmr(i,k) + dp0_frac*dp_icwmr(i,k))/max(frac_limit,cu0_frac)
348 : case (2)
349 0 : sh0 = exp(alpha*sh_icwmr(i,k))
350 0 : dp0 = exp(alpha*dp_icwmr(i,k))
351 0 : cu_icwmr = log((sh0_frac*sh0+dp0_frac*dp0)/max(frac_limit,cu0_frac))
352 1947302 : cu_icwmr = cu_icwmr/alpha
353 : case default ! Area weighted 'arithmetic in emissivity' average.
354 : ! call endrun ('CONV_WATER_4_RAD: Unknown option for conv_water_in_rad - exiting')
355 : end select
356 :
357 : ! Selecting total in-cloud water.
358 : ! Attribute large-scale/convective area fraction differently from default.
359 :
360 1947302 : ls_frac = ast(i,k)
361 1947302 : ls_icwmr = (ls_liq(i,k) + ls_ice(i,k))/max(frac_limit,ls_frac) ! Convert to IC value.
362 1947302 : tot0_frac = (ls_frac + cu0_frac)
363 :
364 1947302 : select case (conv_water_mode) ! Type of average
365 : case (1) ! Area weighted 'arithmetic in emissivity' average
366 1947302 : tot_icwmr = (ls_frac*ls_icwmr + cu0_frac*cu_icwmr)/max(frac_limit,tot0_frac)
367 : case (2)
368 0 : tot_icwmr = log((ls_frac*exp(alpha*ls_icwmr)+cu0_frac*exp(alpha*cu_icwmr))/max(frac_limit,tot0_frac))
369 1947302 : tot_icwmr = tot_icwmr/alpha
370 : case default ! Area weighted 'arithmetic in emissivity' average.
371 : ! call endrun ('CONV_WATER_4_RAD: Unknown option for conv_water_in_rad - exiting')
372 : end select
373 :
374 : end if
375 :
376 : ! Repartition convective cloud water into liquid and ice phase.
377 : ! Currently, this partition is made using the ice fraction of stratus condensate.
378 : ! In future, we should use ice fraction explicitly computed from the convection scheme.
379 :
380 85876200 : conv_ice(i,k) = cu_icwmr * wrk1
381 85876200 : conv_liq(i,k) = cu_icwmr * (1._r8-wrk1)
382 :
383 85876200 : tot_ice(i,k) = tot_icwmr * wrk1
384 85876200 : tot_liq(i,k) = tot_icwmr * (1._r8-wrk1)
385 :
386 85876200 : totg_ice(i,k) = tot0_frac * tot_icwmr * wrk1
387 85876200 : totg_liq(i,k) = tot0_frac * tot_icwmr * (1._r8-wrk1)
388 :
389 : ! Grid-mean convective water
390 85876200 : totg_ice_sh(i,k) = sh0_frac * sh_icwmr(i,k) * wrk1
391 85876200 : totg_ice_dp(i,k) = dp0_frac * dp_icwmr(i,k) * wrk1
392 85876200 : totg_liq_sh(i,k) = sh0_frac * sh_icwmr(i,k) * (1._r8-wrk1)
393 85876200 : totg_liq_dp(i,k) = dp0_frac * dp_icwmr(i,k) * (1._r8-wrk1)
394 85876200 : if( sh0_frac > frac_limit ) then
395 0 : fresh(i,k) = 1._r8
396 : endif
397 85876200 : if( dp0_frac > frac_limit ) then
398 1947302 : fredp(i,k) = 1._r8
399 : endif
400 85876200 : if( cu0_frac > frac_limit ) then
401 1947302 : frecu(i,k) = 1._r8
402 : endif
403 91346832 : if( tot0_frac > frac_limit ) then
404 10806190 : fretot(i,k) = 1._r8
405 : endif
406 :
407 : end do
408 : end do
409 :
410 : !add pbuff calls for COSP
411 58824 : call pbuf_get_field(pbuf, sh_cldliq1_idx, sh_cldliq )
412 58824 : call pbuf_get_field(pbuf, sh_cldice1_idx, sh_cldice )
413 :
414 182752488 : sh_cldliq(:ncol,:pver)=sh_icwmr(:ncol,:pver)*(1-fice(:ncol,:pver))*sh_frac(:ncol,:pver)
415 182752488 : sh_cldice(:ncol,:pver)=sh_icwmr(:ncol,:pver)*fice(:ncol,:pver)*sh_frac(:ncol,:pver)
416 :
417 : ! Output convective IC WMRs
418 :
419 58824 : call outfld( 'ICLMRCU ', conv_liq , pcols, lchnk )
420 58824 : call outfld( 'ICIMRCU ', conv_ice , pcols, lchnk )
421 58824 : call outfld( 'ICLMRTOT', tot_liq , pcols, lchnk )
422 58824 : call outfld( 'ICIMRTOT', tot_ice , pcols, lchnk )
423 :
424 58824 : call outfld('GCLMRDP', totg_liq_dp, pcols, lchnk)
425 58824 : call outfld('GCIMRDP', totg_ice_dp, pcols, lchnk)
426 58824 : call outfld('GCLMRSH', totg_liq_sh, pcols, lchnk)
427 58824 : call outfld('GCIMRSH', totg_ice_sh, pcols, lchnk)
428 58824 : call outfld('FRESH', fresh, pcols, lchnk)
429 58824 : call outfld('FREDP', fredp, pcols, lchnk)
430 58824 : call outfld('FRECU', frecu, pcols, lchnk)
431 58824 : call outfld('FRETOT', fretot, pcols, lchnk)
432 :
433 117648 : end subroutine conv_water_4rad
434 :
435 : end module conv_water
|