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