Line data Source code
1 : module cloud_diagnostics
2 :
3 : !---------------------------------------------------------------------------------
4 : ! Purpose:
5 : !
6 : ! Put cloud physical specifications on the history tape
7 : ! Modified from code that computed cloud optics
8 : !
9 : ! Author: Byron Boville Sept 06, 2002
10 : ! Modified Oct 15, 2008
11 : !
12 : !
13 : !---------------------------------------------------------------------------------
14 :
15 : use shr_kind_mod, only: r8=>shr_kind_r8
16 : use ppgrid, only: pcols, pver,pverp
17 : use physconst, only: gravit
18 : use cam_history, only: outfld
19 : use cam_history, only: addfld, add_default, horiz_only
20 :
21 : implicit none
22 : private
23 : save
24 :
25 : public :: cloud_diagnostics_init
26 : public :: cloud_diagnostics_calc
27 : public :: cloud_diagnostics_register
28 :
29 : ! Local variables
30 : integer :: dei_idx, mu_idx, lambda_idx, iciwp_idx, iclwp_idx, cld_idx ! index into pbuf for cloud fields
31 : integer :: ixcldice, ixcldliq, rei_idx, rel_idx
32 :
33 : logical :: do_cld_diag, mg_clouds, rk_clouds, camrt_rad, spcam_m2005_clouds, spcam_sam1mom_clouds
34 : logical :: one_mom_clouds, two_mom_clouds
35 :
36 : integer :: cicewp_idx = -1
37 : integer :: cliqwp_idx = -1
38 : integer :: cldemis_idx = -1
39 : integer :: cldtau_idx = -1
40 : integer :: nmxrgn_idx = -1
41 : integer :: pmxrgn_idx = -1
42 : integer :: gb_totcldliqmr_idx = -1
43 : integer :: gb_totcldicemr_idx = -1
44 :
45 : ! Index fields for precipitation efficiency.
46 : integer :: acpr_idx, acgcme_idx, acnum_idx
47 :
48 : logical :: use_spcam
49 :
50 : contains
51 :
52 : !===============================================================================
53 1490712 : subroutine cloud_diagnostics_register
54 :
55 : use phys_control, only: phys_getopts
56 : use physics_buffer,only: pbuf_add_field, dtype_r8, dtype_i4
57 :
58 : character(len=16) :: rad_pkg, microp_pgk
59 :
60 1536 : call phys_getopts(radiation_scheme_out=rad_pkg,microp_scheme_out=microp_pgk)
61 1536 : camrt_rad = rad_pkg .eq. 'camrt'
62 1536 : rk_clouds = microp_pgk == 'RK'
63 1536 : mg_clouds = microp_pgk == 'MG'
64 1536 : spcam_m2005_clouds = microp_pgk == 'SPCAM_m2005'
65 1536 : spcam_sam1mom_clouds = microp_pgk == 'SPCAM_sam1mom'
66 1536 : one_mom_clouds = (rk_clouds .or. spcam_sam1mom_clouds)
67 1536 : two_mom_clouds = (mg_clouds .or. spcam_m2005_clouds)
68 :
69 1536 : if (one_mom_clouds) then
70 0 : call pbuf_add_field('CLDEMIS','physpkg', dtype_r8,(/pcols,pver/), cldemis_idx)
71 0 : call pbuf_add_field('CLDTAU', 'physpkg', dtype_r8,(/pcols,pver/), cldtau_idx)
72 :
73 0 : call pbuf_add_field('CICEWP', 'physpkg', dtype_r8,(/pcols,pver/), cicewp_idx)
74 0 : call pbuf_add_field('CLIQWP', 'physpkg', dtype_r8,(/pcols,pver/), cliqwp_idx)
75 :
76 0 : call pbuf_add_field('PMXRGN', 'physpkg', dtype_r8,(/pcols,pverp/), pmxrgn_idx)
77 0 : call pbuf_add_field('NMXRGN', 'physpkg', dtype_i4,(/pcols /), nmxrgn_idx)
78 1536 : else if (two_mom_clouds) then
79 : ! In cloud ice water path for radiation
80 1536 : call pbuf_add_field('ICIWP', 'global', dtype_r8,(/pcols,pver/), iciwp_idx)
81 : ! In cloud liquid water path for radiation
82 1536 : call pbuf_add_field('ICLWP', 'global', dtype_r8,(/pcols,pver/), iclwp_idx)
83 : endif
84 1536 : end subroutine cloud_diagnostics_register
85 :
86 : !===============================================================================
87 1536 : subroutine cloud_diagnostics_init(pbuf2d)
88 : !-----------------------------------------------------------------------
89 1536 : use physics_buffer,only: pbuf_get_index, pbuf_set_field, physics_buffer_desc
90 : use phys_control, only: phys_getopts
91 : use constituents, only: cnst_get_ind
92 : use cloud_cover_diags, only: cloud_cover_diags_init
93 : use time_manager, only: is_first_step
94 :
95 : implicit none
96 :
97 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
98 :
99 : !-----------------------------------------------------------------------
100 :
101 : character(len=16) :: wpunits, sampling_seq
102 : logical :: history_amwg ! output the variables used by the AMWG diag package
103 :
104 :
105 : !-----------------------------------------------------------------------
106 :
107 1536 : cld_idx = pbuf_get_index('CLD')
108 : ! grid box total cloud liquid water mixing ratio (kg/kg)
109 1536 : gb_totcldliqmr_idx = pbuf_get_index('GB_TOTCLDLIQMR')
110 : ! grid box total cloud ice water mixing ratio (kg/kg)
111 1536 : gb_totcldicemr_idx = pbuf_get_index('GB_TOTCLDICEMR')
112 :
113 1536 : call phys_getopts(use_spcam_out=use_spcam)
114 :
115 1536 : if (two_mom_clouds) then
116 :
117 : ! initialize to zero
118 1536 : if (is_first_step()) then
119 768 : call pbuf_set_field(pbuf2d, iciwp_idx, 0._r8)
120 768 : call pbuf_set_field(pbuf2d, iclwp_idx, 0._r8)
121 : end if
122 :
123 3072 : call addfld ('ICWMR', (/ 'lev' /), 'A', 'kg/kg', 'Prognostic in-cloud water mixing ratio')
124 3072 : call addfld ('ICIMR', (/ 'lev' /), 'A', 'kg/kg', 'Prognostic in-cloud ice mixing ratio' )
125 3072 : call addfld ('IWC', (/ 'lev' /), 'A', 'kg/m3', 'Grid box average ice water content' )
126 3072 : call addfld ('LWC', (/ 'lev' /), 'A', 'kg/m3', 'Grid box average liquid water content' )
127 :
128 : ! determine the add_default fields
129 1536 : call phys_getopts(history_amwg_out = history_amwg)
130 :
131 1536 : if (history_amwg) then
132 1536 : call add_default ('ICWMR', 1, ' ')
133 1536 : call add_default ('ICIMR', 1, ' ')
134 1536 : call add_default ('IWC ', 1, ' ')
135 : end if
136 :
137 1536 : dei_idx = pbuf_get_index('DEI')
138 1536 : mu_idx = pbuf_get_index('MU')
139 1536 : lambda_idx = pbuf_get_index('LAMBDAC')
140 :
141 0 : elseif (one_mom_clouds) then
142 :
143 0 : rei_idx = pbuf_get_index('REI')
144 0 : rel_idx = pbuf_get_index('REL')
145 :
146 : endif
147 :
148 1536 : call cnst_get_ind('CLDICE', ixcldice)
149 1536 : call cnst_get_ind('CLDLIQ', ixcldliq)
150 :
151 1536 : do_cld_diag = one_mom_clouds .or. two_mom_clouds
152 :
153 1536 : if (.not.do_cld_diag) return
154 :
155 1536 : if (rk_clouds) then
156 0 : wpunits = 'gram/m2'
157 0 : sampling_seq='rad_lwsw'
158 1536 : else if (two_mom_clouds .or. spcam_sam1mom_clouds) then
159 1536 : wpunits = 'kg/m2'
160 1536 : sampling_seq=''
161 : end if
162 :
163 3072 : call addfld ('ICLDIWP', (/ 'lev' /), 'A', wpunits,'In-cloud ice water path' , sampling_seq=sampling_seq)
164 : call addfld ('ICLDTWP', (/ 'lev' /), 'A',wpunits,'In-cloud cloud total water path (liquid and ice)', &
165 3072 : sampling_seq=sampling_seq)
166 :
167 : call addfld ('GCLDLWP', (/ 'lev' /), 'A',wpunits,'Grid-box cloud water path' , &
168 3072 : sampling_seq=sampling_seq)
169 : call addfld ('TGCLDCWP',horiz_only, 'A',wpunits,'Total grid-box cloud water path (liquid and ice)', &
170 1536 : sampling_seq=sampling_seq)
171 : call addfld ('TGCLDLWP',horiz_only, 'A',wpunits,'Total grid-box cloud liquid water path', &
172 1536 : sampling_seq=sampling_seq)
173 : call addfld ('TGCLDIWP',horiz_only, 'A',wpunits,'Total grid-box cloud ice water path' , &
174 1536 : sampling_seq=sampling_seq)
175 :
176 1536 : if(two_mom_clouds) then
177 3072 : call addfld ('lambda_cloud',(/ 'lev' /),'I','1/meter','lambda in cloud')
178 3072 : call addfld ('mu_cloud', (/ 'lev' /),'I','1','mu in cloud')
179 3072 : call addfld ('dei_cloud', (/ 'lev' /),'I','micrometers','ice radiative effective diameter in cloud')
180 : endif
181 :
182 1536 : if(one_mom_clouds) then
183 0 : call addfld ('rel_cloud',(/ 'lev' /),'I','1/meter','effective radius of liq in cloud', sampling_seq=sampling_seq)
184 0 : call addfld ('rei_cloud',(/ 'lev' /),'I','1','effective radius of ice in cloud', sampling_seq=sampling_seq)
185 : endif
186 :
187 3072 : call addfld ('SETLWP',(/ 'lev' /), 'A','gram/m2','Prescribed liquid water path' , sampling_seq=sampling_seq)
188 1536 : call addfld ('LWSH',horiz_only, 'A','m','Liquid water scale height' , sampling_seq=sampling_seq)
189 :
190 3072 : call addfld ('EFFCLD',(/ 'lev' /), 'A','fraction','Effective cloud fraction' , sampling_seq=sampling_seq)
191 :
192 1536 : if (camrt_rad) then
193 0 : call addfld ('EMIS', (/ 'lev' /), 'A', '1','cloud emissivity' , sampling_seq=sampling_seq)
194 : else
195 3072 : call addfld ('EMISCLD', (/ 'lev' /), 'A', '1','cloud emissivity' , sampling_seq=sampling_seq)
196 : endif
197 :
198 1536 : call cloud_cover_diags_init(sampling_seq)
199 :
200 : ! ----------------------------
201 : ! determine default variables
202 : ! ----------------------------
203 1536 : call phys_getopts( history_amwg_out = history_amwg)
204 :
205 1536 : if (history_amwg) then
206 1536 : call add_default ('TGCLDLWP', 1, ' ')
207 1536 : call add_default ('TGCLDIWP', 1, ' ')
208 1536 : call add_default ('TGCLDCWP', 1, ' ')
209 1536 : if(rk_clouds) then
210 0 : if (camrt_rad) then
211 0 : call add_default ('EMIS', 1, ' ')
212 : else
213 0 : call add_default ('EMISCLD', 1, ' ')
214 : endif
215 : endif
216 : endif
217 :
218 : return
219 1536 : end subroutine cloud_diagnostics_init
220 :
221 1489176 : subroutine cloud_diagnostics_calc(state, pbuf)
222 : !===============================================================================
223 : !
224 : ! Compute (liquid+ice) water path and cloud water/ice diagnostics
225 : ! *** soon this code will compute liquid and ice paths from input liquid and ice mixing ratios
226 : !
227 : ! **** mixes interface and physics code temporarily
228 : !-----------------------------------------------------------------------
229 1536 : use physics_types, only: physics_state
230 : use physics_buffer,only: physics_buffer_desc, pbuf_get_field, pbuf_old_tim_idx
231 : use pkg_cldoptics, only: cldovrlap, cldclw, cldems
232 : use conv_water, only: conv_water_in_rad, conv_water_4rad
233 : use radiation, only: radiation_do
234 : use cloud_cover_diags, only: cloud_cover_diags_out
235 :
236 : use ref_pres, only: top_lev=>trop_cloud_top_lev
237 :
238 : implicit none
239 :
240 : ! Arguments
241 : type(physics_state), intent(in) :: state ! state variables
242 : type(physics_buffer_desc), pointer :: pbuf(:)
243 :
244 : ! Local variables
245 :
246 1489176 : real(r8), pointer :: cld(:,:) ! cloud fraction
247 1489176 : real(r8), pointer :: iciwp(:,:) ! in-cloud cloud ice water path
248 1489176 : real(r8), pointer :: iclwp(:,:) ! in-cloud cloud liquid water path
249 1489176 : real(r8), pointer :: dei(:,:) ! effective radiative diameter of ice
250 1489176 : real(r8), pointer :: mu(:,:) ! gamma distribution for liq clouds
251 1489176 : real(r8), pointer :: lambda(:,:) ! gamma distribution for liq clouds
252 1489176 : real(r8), pointer :: rei(:,:) ! effective radiative radius of ice
253 1489176 : real(r8), pointer :: rel(:,:) ! effective radiative radius of liq
254 :
255 1489176 : real(r8), pointer :: cldemis(:,:) ! cloud emissivity
256 1489176 : real(r8), pointer :: cldtau(:,:) ! cloud optical depth
257 1489176 : real(r8), pointer :: cicewp(:,:) ! in-cloud cloud ice water path
258 1489176 : real(r8), pointer :: cliqwp(:,:) ! in-cloud cloud liquid water path
259 :
260 1489176 : integer, pointer :: nmxrgn(:) ! Number of maximally overlapped regions
261 1489176 : real(r8), pointer :: pmxrgn(:,:) ! Maximum values of pressure for each
262 :
263 1489176 : real(r8), pointer :: totg_ice(:,:) ! grid box total cloud ice mixing ratio
264 1489176 : real(r8), pointer :: totg_liq(:,:) ! grid box total cloud liquid mixing ratio
265 :
266 : integer :: itim_old
267 :
268 : real(r8) :: cwp (pcols,pver) ! in-cloud cloud (total) water path
269 : real(r8) :: gicewp(pcols,pver) ! grid-box cloud ice water path
270 : real(r8) :: gliqwp(pcols,pver) ! grid-box cloud liquid water path
271 : real(r8) :: gwp (pcols,pver) ! grid-box cloud (total) water path
272 : real(r8) :: tgicewp(pcols) ! Vertically integrated ice water path
273 : real(r8) :: tgliqwp(pcols) ! Vertically integrated liquid water path
274 : real(r8) :: tgwp (pcols) ! Vertically integrated (total) cloud water path
275 :
276 : real(r8) :: ficemr (pcols,pver) ! Ice fraction from ice and liquid mixing ratios
277 :
278 : real(r8) :: icimr(pcols,pver) ! In cloud ice mixing ratio
279 : real(r8) :: icwmr(pcols,pver) ! In cloud water mixing ratio
280 : real(r8) :: iwc(pcols,pver) ! Grid box average ice water content
281 : real(r8) :: lwc(pcols,pver) ! Grid box average liquid water content
282 :
283 : ! old data
284 : real(r8) :: tpw (pcols) ! total precipitable water
285 : real(r8) :: clwpold(pcols,pver) ! Presribed cloud liq. h2o path
286 : real(r8) :: hl (pcols) ! Liquid water scale height
287 :
288 : integer :: i,k ! loop indexes
289 : integer :: ncol, lchnk
290 : real(r8) :: rgrav
291 :
292 : real(r8) :: allcld_ice (pcols,pver) ! Convective cloud ice
293 : real(r8) :: allcld_liq (pcols,pver) ! Convective cloud liquid
294 :
295 : real(r8) :: effcld(pcols,pver) ! effective cloud=cld*emis
296 :
297 : logical :: dosw,dolw
298 :
299 : !-----------------------------------------------------------------------
300 1489176 : if (.not.do_cld_diag) return
301 :
302 1489176 : if(one_mom_clouds) then
303 0 : dosw = radiation_do('sw') ! do shortwave heating calc this timestep?
304 0 : dolw = radiation_do('lw') ! do longwave heating calc this timestep?
305 : else
306 : dosw = .true.
307 : dolw = .true.
308 : endif
309 :
310 1489176 : if (.not.(dosw .or. dolw)) return
311 :
312 1489176 : ncol = state%ncol
313 1489176 : lchnk = state%lchnk
314 :
315 1489176 : itim_old = pbuf_old_tim_idx()
316 5956704 : call pbuf_get_field(pbuf, cld_idx, cld, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
317 :
318 1489176 : call pbuf_get_field(pbuf, gb_totcldicemr_idx, totg_ice)
319 1489176 : call pbuf_get_field(pbuf, gb_totcldliqmr_idx, totg_liq)
320 :
321 1489176 : if(two_mom_clouds)then
322 :
323 1489176 : call pbuf_get_field(pbuf, iclwp_idx, iclwp )
324 1489176 : call pbuf_get_field(pbuf, iciwp_idx, iciwp )
325 1489176 : call pbuf_get_field(pbuf, dei_idx, dei )
326 1489176 : call pbuf_get_field(pbuf, mu_idx, mu )
327 1489176 : call pbuf_get_field(pbuf, lambda_idx, lambda )
328 :
329 1489176 : call outfld('dei_cloud',dei(:,:),pcols,lchnk)
330 1489176 : call outfld('mu_cloud',mu(:,:),pcols,lchnk)
331 1489176 : call outfld('lambda_cloud',lambda(:,:),pcols,lchnk)
332 :
333 0 : elseif(one_mom_clouds) then
334 :
335 0 : call pbuf_get_field(pbuf, rei_idx, rei )
336 0 : call pbuf_get_field(pbuf, rel_idx, rel )
337 :
338 0 : call outfld('rel_cloud', rel, pcols, lchnk)
339 0 : call outfld('rei_cloud', rei, pcols, lchnk)
340 :
341 0 : if (cldemis_idx>0) then
342 0 : call pbuf_get_field(pbuf, cldemis_idx, cldemis )
343 : else
344 0 : allocate(cldemis(pcols,pver))
345 : endif
346 0 : if (cldtau_idx>0) then
347 0 : call pbuf_get_field(pbuf, cldtau_idx, cldtau )
348 : else
349 0 : allocate(cldtau(pcols,pver))
350 : endif
351 :
352 : endif
353 :
354 1489176 : if (cicewp_idx>0) then
355 0 : call pbuf_get_field(pbuf, cicewp_idx, cicewp )
356 : else
357 1489176 : allocate(cicewp(pcols,pver))
358 : endif
359 1489176 : if (cliqwp_idx>0) then
360 0 : call pbuf_get_field(pbuf, cliqwp_idx, cliqwp )
361 : else
362 1489176 : allocate(cliqwp(pcols,pver))
363 : endif
364 :
365 1489176 : if (nmxrgn_idx>0) then
366 0 : call pbuf_get_field(pbuf, nmxrgn_idx, nmxrgn )
367 : else
368 1489176 : allocate(nmxrgn(pcols))
369 : endif
370 :
371 1489176 : if (pmxrgn_idx>0) then
372 0 : call pbuf_get_field(pbuf, pmxrgn_idx, pmxrgn )
373 : else
374 1489176 : allocate(pmxrgn(pcols,pverp))
375 : endif
376 :
377 : ! Compute liquid and ice water paths
378 1489176 : if(two_mom_clouds) then
379 :
380 : ! ----------------------------------------------------------- !
381 : ! Adjust in-cloud water values to take account of convective !
382 : ! in-cloud water. It is used to calculate the values of !
383 : ! iclwp and iciwp to pass to the radiation. !
384 : ! ----------------------------------------------------------- !
385 1489176 : if( conv_water_in_rad /= 0 ) then
386 1489176 : call conv_water_4rad(state, pbuf)
387 2314006344 : allcld_ice(:ncol,:) = totg_ice(:ncol,:) ! Grid-avg all cloud liquid
388 2314006344 : allcld_liq(:ncol,:) = totg_liq(:ncol,:) ! Grid-avg all cloud ice
389 : else
390 0 : allcld_liq(:ncol,top_lev:pver) = state%q(:ncol,top_lev:pver,ixcldliq) ! Grid-ave all cloud liquid
391 0 : allcld_ice(:ncol,top_lev:pver) = state%q(:ncol,top_lev:pver,ixcldice) ! " ice
392 : end if
393 :
394 : ! ------------------------------------------------------------ !
395 : ! Compute in cloud ice and liquid mixing ratios !
396 : ! Note that 'iclwp, iciwp' are used for radiation computation. !
397 : ! ------------------------------------------------------------ !
398 :
399 :
400 2355876432 : iciwp = 0._r8
401 2355876432 : iclwp = 0._r8
402 1489176 : icimr = 0._r8
403 1489176 : icwmr = 0._r8
404 1489176 : iwc = 0._r8
405 1489176 : lwc = 0._r8
406 :
407 126579960 : do k = top_lev, pver
408 2090214360 : do i = 1, ncol
409 : ! Limits for in-cloud mixing ratios consistent with MG microphysics
410 : ! in-cloud mixing ratio maximum limit of 0.005 kg/kg
411 1963634400 : icimr(i,k) = min( allcld_ice(i,k) / max(0.0001_r8,cld(i,k)),0.005_r8 )
412 1963634400 : icwmr(i,k) = min( allcld_liq(i,k) / max(0.0001_r8,cld(i,k)),0.005_r8 )
413 1963634400 : iwc(i,k) = allcld_ice(i,k) * state%pmid(i,k) / (287.15_r8*state%t(i,k))
414 1963634400 : lwc(i,k) = allcld_liq(i,k) * state%pmid(i,k) / (287.15_r8*state%t(i,k))
415 : ! Calculate total cloud water paths in each layer
416 1963634400 : iciwp(i,k) = icimr(i,k) * state%pdel(i,k) / gravit
417 2088725184 : iclwp(i,k) = icwmr(i,k) * state%pdel(i,k) / gravit
418 : end do
419 : end do
420 :
421 139982544 : do k=1,pver
422 2314006344 : do i = 1,ncol
423 2174023800 : gicewp(i,k) = iciwp(i,k)*cld(i,k)
424 2174023800 : gliqwp(i,k) = iclwp(i,k)*cld(i,k)
425 2174023800 : cicewp(i,k) = iciwp(i,k)
426 2312517168 : cliqwp(i,k) = iclwp(i,k)
427 : end do
428 : end do
429 :
430 0 : elseif(one_mom_clouds) then
431 :
432 0 : if (conv_water_in_rad /= 0) then
433 0 : call conv_water_4rad(state, pbuf)
434 0 : allcld_ice(:ncol,:) = totg_ice(:ncol,:) ! Grid-avg all cloud liquid
435 0 : allcld_liq(:ncol,:) = totg_liq(:ncol,:) ! Grid-avg all cloud ice
436 : else
437 0 : allcld_liq = state%q(:,:,ixcldliq)
438 0 : allcld_ice = state%q(:,:,ixcldice)
439 : end if
440 :
441 0 : do k=1,pver
442 0 : do i = 1,ncol
443 0 : gicewp(i,k) = allcld_ice(i,k)*state%pdel(i,k)/gravit*1000.0_r8 ! Grid box ice water path.
444 0 : gliqwp(i,k) = allcld_liq(i,k)*state%pdel(i,k)/gravit*1000.0_r8 ! Grid box liquid water path.
445 0 : cicewp(i,k) = gicewp(i,k) / max(0.01_r8,cld(i,k)) ! In-cloud ice water path.
446 0 : cliqwp(i,k) = gliqwp(i,k) / max(0.01_r8,cld(i,k)) ! In-cloud liquid water path.
447 0 : ficemr(i,k) = allcld_ice(i,k) / max(1.e-10_r8,(allcld_ice(i,k) + allcld_liq(i,k)))
448 : end do
449 : end do
450 : endif
451 :
452 : ! Determine parameters for maximum/random overlap
453 1489176 : call cldovrlap(lchnk, ncol, state%pint, cld, nmxrgn, pmxrgn)
454 :
455 1489176 : if(.not. use_spcam) then ! in spcam, these diagnostics are calcluated in crm_physics.F90
456 : ! Cloud cover diagnostics (done in radiation_tend for camrt)
457 1489176 : if (.not.camrt_rad) then
458 1489176 : call cloud_cover_diags_out(lchnk, ncol, cld, state%pmid, nmxrgn, pmxrgn )
459 : endif
460 : end if
461 :
462 24865776 : tgicewp(:ncol) = 0._r8
463 24865776 : tgliqwp(:ncol) = 0._r8
464 :
465 139982544 : do k=1,pver
466 2312517168 : tgicewp(:ncol) = tgicewp(:ncol) + gicewp(:ncol,k)
467 2314006344 : tgliqwp(:ncol) = tgliqwp(:ncol) + gliqwp(:ncol,k)
468 : end do
469 :
470 24865776 : tgwp(:ncol) = tgicewp(:ncol) + tgliqwp(:ncol)
471 2314006344 : gwp(:ncol,:pver) = gicewp(:ncol,:pver) + gliqwp(:ncol,:pver)
472 2314006344 : cwp(:ncol,:pver) = cicewp(:ncol,:pver) + cliqwp(:ncol,:pver)
473 :
474 1489176 : if(one_mom_clouds) then
475 :
476 : ! Cloud emissivity.
477 0 : call cldems(lchnk, ncol, cwp, ficemr, rei, cldemis, cldtau)
478 :
479 : ! Effective cloud cover
480 0 : do k=1,pver
481 0 : do i=1,ncol
482 0 : effcld(i,k) = cld(i,k)*cldemis(i,k)
483 : end do
484 : end do
485 :
486 0 : call outfld('EFFCLD' ,effcld , pcols,lchnk)
487 0 : if (camrt_rad) then
488 0 : call outfld('EMIS' ,cldemis, pcols,lchnk)
489 : else
490 0 : call outfld('EMISCLD' ,cldemis, pcols,lchnk)
491 : endif
492 :
493 1489176 : else if (two_mom_clouds) then
494 :
495 : ! --------------------------------------------- !
496 : ! General outfield calls for microphysics !
497 : ! --------------------------------------------- !
498 :
499 1489176 : call outfld( 'IWC' , iwc, pcols, lchnk )
500 1489176 : call outfld( 'LWC' , lwc, pcols, lchnk )
501 1489176 : call outfld( 'ICIMR' , icimr, pcols, lchnk )
502 1489176 : call outfld( 'ICWMR' , icwmr, pcols, lchnk )
503 :
504 : endif
505 :
506 1489176 : if (.not. use_spcam) then
507 : ! for spcam, these are diagnostics in crm_physics.F90
508 1489176 : call outfld('GCLDLWP' ,gwp , pcols,lchnk)
509 1489176 : call outfld('TGCLDCWP',tgwp , pcols,lchnk)
510 1489176 : call outfld('TGCLDLWP',tgliqwp, pcols,lchnk)
511 1489176 : call outfld('TGCLDIWP',tgicewp, pcols,lchnk)
512 1489176 : call outfld('ICLDTWP' ,cwp , pcols,lchnk)
513 1489176 : call outfld('ICLDIWP' ,cicewp , pcols,lchnk)
514 : endif
515 :
516 : ! Compute total preciptable water in column (in mm)
517 24865776 : tpw(:ncol) = 0.0_r8
518 1489176 : rgrav = 1.0_r8/gravit
519 139982544 : do k=1,pver
520 2314006344 : do i=1,ncol
521 2312517168 : tpw(i) = tpw(i) + state%pdel(i,k)*state%q(i,k,1)*rgrav
522 : end do
523 : end do
524 :
525 : ! Diagnostic liquid water path (old specified form)
526 :
527 1489176 : call cldclw(lchnk, ncol, state%zi, clwpold, tpw, hl)
528 1489176 : call outfld('SETLWP' ,clwpold, pcols,lchnk)
529 1489176 : call outfld('LWSH' ,hl , pcols,lchnk)
530 :
531 1489176 : if(one_mom_clouds) then
532 0 : if (cldemis_idx<0) deallocate(cldemis)
533 0 : if (cldtau_idx<0) deallocate(cldtau)
534 : endif
535 1489176 : if (cicewp_idx<0) deallocate(cicewp)
536 1489176 : if (cliqwp_idx<0) deallocate(cliqwp)
537 1489176 : if (pmxrgn_idx<0) deallocate(pmxrgn)
538 1489176 : if (nmxrgn_idx<0) deallocate(nmxrgn)
539 :
540 : return
541 2978352 : end subroutine cloud_diagnostics_calc
542 :
543 : end module cloud_diagnostics
|