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