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