CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
radiation_clouds.f
1
4
5! module_radiation_clouds description !!!!!
6! ========================================================== !!!!!
7! !
8! the 'radiation_clouds.f' contains: !
9! !
10! 'module_radiation_clouds' --- compute cloud related quantities!
11! for radiation computations !
12! !
13! the following are the externally accessable subroutines: !
14! !
15! 'cld_init' --- initialization routine !
16! inputs: !
17! ( si, NLAY, imp_physics, me ) !
18! outputs: !
19! ( errflg, errmsg ) !
20! !
21! 'radiation_clouds_prop' --- radiation cloud properties !
22! obtained from various cloud schemes !
23! inputs: !
24! (plyr,plvl,tlyr,tvly,qlyr,qstl,rhly, !
25! ccnd, ncndl, cnvw, cnvc, tracer1, !
26! xlat,xlon,slmsk,dz,delp, IX, LM, NLAY, NLP1, !
27! deltaq, sup, me, icloud, kdt, !
28! ntrac, ntcw, ntiw, ntrw, ntsw, ntgl, ntclamt, !
29! imp_physics, imp_physics_nssl, imp_physics_fer_hires, !
30! imp_physics_gfdl, imp_physics_thompson, imp_physics_wsm6, !
31! imp_physics_zhao_carr, imp_physics_zhao_carr_pdf, !
32! imp_physics_mg, iovr, iovr_rand, iovr_maxrand, iovr_max, !
33! iovr_dcorr, iovr_exp, iovr_exprand, idcor, idcor_con, !
34! idcor_hogan, idcor_oreopoulos, !
35! imfdeepcnv, imfdeepcnv_gf, imfdeepcnv_c3, do_mynnedmf, lgfdlmprad, !
36! uni_cld, lmfshal, lmfdeep2, cldcov, clouds1, !
37! effrl, effri, effrr, effrs, effr_in, !
38! effrl_inout, effri_inout, effrs_inout, !
39! lwp_ex, iwp_ex, lwp_fc, iwp_fc, !
40! dzlay, latdeg, julian, yearlen, gridkm, !
41! outputs: !
42! cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, !
43! cld_rwp, cld_rerain, cld_swp, cld_resnow, !
44! clds,mtop,mbot,de_lgth,alpha) !
45! !
46! internal/external accessable subroutines: !
47! 'progcld_zhao_carr' --- zhao/moorthi prognostic cloud scheme !
48! 'progcld_zhao_carr_pdf' --- zhao/moorthi prognostic cloud + pdfcld !
49! 'progcld_gfdl_lin' --- GFDL-Lin cloud microphysics !
50! 'progcld_fer_hires' --- Ferrier-Aligo cloud microphysics !
51! 'progcld_thompson_wsm6' --- Thompson/wsm6 cloud microphysics (EMC) !
52! 'progclduni' --- MG2/3 cloud microphysics !
53! (with/without SHOC) (EMC) !
54! also used by GFDL MP (EMC) !
55! 'progcld_thompson' --- Thompson MP (added by G. Thompson) !
56! 'gethml' --- get diagnostic hi, mid, low clouds !
57! !
58! cloud property array description: !
59! cld_frac (:,:) - layer total cloud fraction !
60! cld_lwp (:,:) - layer cloud liq water path !
61! cld_reliq (:,:) - mean effective radius for liquid cloud !
62! cld_iwp (:,:) - layer cloud ice water path !
63! cld_reice (:,:) - mean effective radius for ice cloud !
64! cld_rwp (:,:) - layer rain drop water path !
65! cld_rerain(:,:) - mean effective radius for rain drop !
66! ** cld_swp (:,:) - layer snow flake water path !
67! cld_resnow(:,:) - mean effective radius for snow flake !
68! ** fu's scheme need to be normalized by snow density (g/m**3/1.0e6)!
69! !
70! external modules referenced: !
71! 'module module_microphysics' in 'module_bfmicrophysics.f' !
72! !
73! program history log: !
74! nov 1992, y.h., k.a.c, a.k. - cloud parameterization !
75! 'cldjms' patterned after slingo and slingo's work (jgr, !
76! 1992), stratiform clouds are allowed in any layer except !
77! the surface and upper stratosphere. the relative humidity !
78! criterion may cery in different model layers. !
79! mar 1993, kenneth campana - created original crhtab tables !
80! for critical rh look up references.
81! feb 1994, kenneth campana - modified to use only one table !
82! for all forecast hours. !
83! oct 1995, kenneth campana - tuned cloud rh curves !
84! rh-cld relation from tables created using mitchell-hahn !
85! tuning technique on airforce rtneph observations. !
86! nov 1995, kenneth campana - the bl relationships used !
87! below llyr, except in marine stratus regions. !
88! apr 1996, kenneth campana - save bl cld amt in cld(,5) !
89! aug 1997, kenneth campana - smooth out last bunch of bins !
90! of the rh look up tables !
91! dec 1998, s. moorthi - added prognostic cloud method !
92! apr 2003, yu-tai hou - rewritten in frotran 90 !
93! modulized form 'module_rad_clouds' from combining the original!
94! subroutines 'cldjms', 'cldprp', and 'gcljms'. and seperated !
95! prognostic and diagnostic methods into two packages. !
96! --- 2003, s. moorthi - adapted b. ferrier's prognostic !
97! cloud scheme to ncep gfs model as additional option. !
98! apr 2004, yu-tai hou - separated calculation of the !
99! averaged h,m,l,bl cloud amounts from each of the cld schemes !
100! to become an shared individule subprogram 'gethml'. !
101! apr 2005, yu-tai hou - modified cloud array and module !
102! structures. !
103! dec 2008, yu-tai hou - changed low-cld calculation, !
104! now cantains clds from sfc layer and upward to the low/mid !
105! boundary (include bl-cld). h,m,l clds domain boundaries are !
106! adjusted for better agreement with observations. !
107! jan 2011, yu-tai hou - changed virtual temperature !
108! as input variable instead of originally computed inside the !
109! two prognostic cld schemes 'progcld_zhao_carr' !
110! aug 2012, yu-tai hou - modified subroutine cld_init !
111! to pass all fixed control variables at the start. and set !
112! their correponding internal module variables to be used by !
113! module subroutines. !
114! nov 2012, yu-tai hou - modified control parameters !
115! thru module 'physparam'. !
116! apr 2013, h.lin/y.hou - corrected error in calculating !
117! llyr for different vertical indexing directions. !
118! jul 2013, r. sun/h. pan - modified to use pdf cloud and !
119! convective cloud cover and water for radiation !
120! !
121! jul 2014 s. moorthi - merging with gfs version !
122! feb 2017 a. cheng - add odepth output, effective radius input !
123! Jan 2018 S Moorthi - update to include physics from ipdv4 !
124! jun 2018 h-m lin/y-t hou - removed the legacy subroutine !
125! 'diagcld1' for diagnostic cloud scheme, added new cloud !
126! overlapping method of de-correlation length, and optimized !
127! the code structure. !
128! jul 2020, m.j. iacono - added rrtmg/mcica cloud overlap options !
129! exponential and exponential-random. each method can use !
130! either a constant or a latitude-varying and day-of-year !
131! varying decorrelation length selected with parameter "idcor". !
132! !
133! Jan 2022, Q.Liu - add subroutine radiation_clouds_prop, and !
134! move all the subroutine call "progcld*" from !
135! GFS_rrtmg_pre.F90 to this new subroutine !
136!!!!! ========================================================== !!!!!
137!!!!! end descriptions !!!!!
138!!!!! ========================================================== !!!!!
139
164
167!
168 use module_iounitdef, only : nicltun
170 & get_alpha_exper
171 use machine, only : kind_phys
172!
173 implicit none
174!
175 private
176
177! --- version tag and last revision date
178 character(40), parameter :: &
179 & VTAGCLD='NCEP-Radiation_clouds v5.1 Nov 2012 '
180! & VTAGCLD='NCEP-Radiation_clouds v5.0 Aug 2012 '
181
182! --- set constant parameters
183 real (kind=kind_phys) :: gfac,gord
184
185 integer, parameter, public :: nf_clds = 9
186 integer, parameter, public :: nk_clds = 3
187
188! pressure limits of cloud domain interfaces (low,mid,high) in mb (0.1kPa)
189 real (kind=kind_phys), save :: ptopc(nk_clds+1,2)
191
192!org data ptopc / 1050., 642., 350., 0.0, 1050., 750., 500., 0.0 /
193 data ptopc / 1050., 650., 400., 0.0, 1050., 750., 500., 0.0 /
194
195! real (kind=kind_phys), parameter :: climit = 0.01
196 real (kind=kind_phys), parameter :: climit = 0.001, climit2=0.05
197 real (kind=kind_phys), parameter :: ovcst = 1.0 - 1.0e-8
198
199 real (kind=kind_phys), parameter :: reliq_def = 10.0
200 real (kind=kind_phys), parameter :: reice_def = 50.0
201 real (kind=kind_phys), parameter :: rrain_def = 1000.0
202 real (kind=kind_phys), parameter :: rsnow_def = 250.0
203 real (kind=kind_phys), parameter :: creice_def = 25.0
204
205 real (kind=kind_phys), parameter :: cldssa_def = 0.99
206 real (kind=kind_phys), parameter :: cldasy_def = 0.84
207
208 integer :: llyr = 2
209
210 ! Default ice crystal sizes vs. temperature following Kristjansson and Mitchell
211 real (kind=kind_phys), dimension(95), parameter :: retab =(/ &
212 & 5.92779, 6.26422, 6.61973, 6.99539, 7.39234, &
213 & 7.81177, 8.25496, 8.72323, 9.21800, 9.74075, 10.2930, &
214 & 10.8765, 11.4929, 12.1440, 12.8317, 13.5581, 14.2319, &
215 & 15.0351, 15.8799, 16.7674, 17.6986, 18.6744, 19.6955, &
216 & 20.7623, 21.8757, 23.0364, 24.2452, 25.5034, 26.8125, &
217 & 27.7895, 28.6450, 29.4167, 30.1088, 30.7306, 31.2943, &
218 & 31.8151, 32.3077, 32.7870, 33.2657, 33.7540, 34.2601, &
219 & 34.7892, 35.3442, 35.9255, 36.5316, 37.1602, 37.8078, &
220 & 38.4720, 39.1508, 39.8442, 40.5552, 41.2912, 42.0635, &
221 & 42.8876, 43.7863, 44.7853, 45.9170, 47.2165, 48.7221, &
222 & 50.4710, 52.4980, 54.8315, 57.4898, 60.4785, 63.7898, &
223 & 65.5604, 71.2885, 75.4113, 79.7368, 84.2351, 88.8833, &
224 & 93.6658, 98.5739, 103.603, 108.752, 114.025, 119.424, &
225 & 124.954, 130.630, 136.457, 142.446, 148.608, 154.956, &
226 & 161.503, 168.262, 175.248, 182.473, 189.952, 197.699, &
227 & 205.728, 214.055, 222.694, 231.661, 240.971, 250.639/)
228
235
236! =================
237 contains
238! =================
239
255 subroutine cld_init &
256 & ( si, nlay, imp_physics, me, con_g, con_rd, errflg, errmsg )
257! =================================================================== !
258! !
259! abstract: cld_init is an initialization program for cloud-radiation !
260! calculations. it sets up boundary layer cloud top. !
261! !
262! !
263! inputs: !
264! si (L+1) : model vertical sigma layer interface !
265! NLAY : vertical layer number !
266! imp_physics : MP identifier !
267! me : print control flag !
268! imp_physics : cloud microphysics scheme control flag !
269! !
270! outputs: !
271! errflg : CCPP error flag !
272! errmsg : CCPP error message !
273! !
274! usage: call cld_init !
275! !
276! subroutines called: rhtable !
277! !
278! =================================================================== !
279!
280 implicit none
281
282! --- inputs:
283 integer, intent(in) :: nlay, me, imp_physics
284
285 real (kind=kind_phys), intent(in) :: si(:), con_g, con_rd
286
287! --- outputs:
288 integer, intent(out) :: errflg
289 character(len=*), intent(out) :: errmsg
290
291!
292!===> ... begin here
293!
294! Initialize CCPP error handling variables
295 errmsg = ''
296 errflg = 0
297
298 ! Initialze module parameters
299 gfac = 1.0e5/con_g
300 gord = con_g/con_rd
301
302 if (me == 0) then
303 print *, vtagcld !print out version tag
304 print *,' - Using Prognostic Cloud Method'
305 if (imp_physics == 99) then
306 print *,' --- Zhao/Carr/Sundqvist microphysics'
307 elseif (imp_physics == 98) then
308 print *,' --- zhao/carr/sundqvist + pdf cloud'
309 elseif (imp_physics == 11) then
310 print *,' --- GFDL Lin cloud microphysics'
311 elseif (imp_physics == 8) then
312 print *,' --- Thompson cloud microphysics'
313 elseif (imp_physics == 88) then
314 print *,' --- TEMPO cloud microphysics'
315 elseif (imp_physics == 6) then
316 print *,' --- WSM6 cloud microphysics'
317 elseif (imp_physics == 10) then
318 print *,' --- MG cloud microphysics'
319 elseif (imp_physics == 15) then
320 print *,' --- Ferrier-Aligo cloud microphysics'
321 elseif (imp_physics == 17) then
322 print *,' --- NSSL cloud microphysics'
323 else
324 print *,' !!! ERROR in cloud microphysc specification!!!', &
325 & ' imp_physics (NP3D) =',imp_physics
326 errflg = 1
327 errmsg = 'ERROR(cld_init): cloud mp specification is not'// &
328 & ' valid'
329 return
330 endif
331 endif
332!
333!...................................
334 end subroutine cld_init
335!-----------------------------------
336
341 & ( plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, & ! --- inputs:
342 & ccnd, ncndl, cnvw, cnvc, tracer1, &
343 & xlat, xlon, slmsk, dz, delp, ix, lm, nlay, nlp1, &
344 & deltaq, sup, dcorr_con, me, icloud, kdt, &
345 & ntrac, ntcw, ntiw, ntrw, ntsw, ntgl, ntclamt, &
346 & imp_physics, imp_physics_nssl, imp_physics_fer_hires, &
347 & imp_physics_gfdl, imp_physics_thompson, imp_physics_wsm6, &
348 & imp_physics_tempo, &
349 & imp_physics_zhao_carr, imp_physics_zhao_carr_pdf, &
350 & imp_physics_mg, iovr, iovr_rand, iovr_maxrand, iovr_max, &
351 & iovr_dcorr, iovr_exp, iovr_exprand, idcor, idcor_con, &
352 & idcor_hogan, idcor_oreopoulos, lcrick, lcnorm, &
353 & imfdeepcnv, imfdeepcnv_gf, imfdeepcnv_c3, &
354 & do_mynnedmf, lgfdlmprad, xr_cnvcld, &
355 & uni_cld, lmfshal, lmfdeep2, cldcov, clouds1, &
356 & effrl, effri, effrr, effrs, effr_in, &
357 & effrl_inout, effri_inout, effrs_inout, &
358 & lwp_ex, iwp_ex, lwp_fc, iwp_fc, &
359 & dzlay, latdeg, julian, yearlen, gridkm, top_at_1, si, &
360 & xr_con, xr_exp, con_ttp, con_pi, con_g, con_rd, con_thgni, &
361 & cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, & ! --- outputs:
362 & cld_rwp, cld_rerain, cld_swp, cld_resnow, &
363 & clds, mtop, mbot, de_lgth, alpha &
364 & )
365
366! ================= subprogram documentation block ================ !
367! !
368! subprogram: radiation_clouds_prop computes cloud related quantities using !
369! zhao/moorthi's prognostic cloud microphysics scheme. !
370! !
371! abstract: this program computes cloud fractions from cloud !
372! condensates, calculates liquid/ice cloud droplet effective radius, !
373! and computes the low, mid, high, total and boundary layer cloud !
374! fractions and the vertical indices of low, mid, and high cloud !
375! top and base. the three vertical cloud domains are set up in the !
376! initial subroutine "cld_init". !
377! !
378! usage: call radiation_clouds_prop !
379! !
380! subprograms called: !
381! !
382! 'progcld_zhao_carr' --- zhao/moorthi prognostic cloud scheme !
383! 'progcld_zhao_carr_pdf' --- zhao/moorthi prognostic cloud + pdfcld !
384! 'progcld_gfdl_lin' --- GFDL-Lin cloud microphysics !
385! 'progcld_fer_hires' --- Ferrier-Aligo cloud microphysics !
386! 'progcld_thompson_wsm6' --- Thompson/wsm6 cloud microphysics (EMC) !
387! 'progclduni' --- MG cloud microphysics !
388! --- GFDL cloud microphysics (EMC) !
389! --- Thompson + MYNN PBL (or GF convection) !
390! 'progcld_thompson' --- Thompson MP (added by G. Thompson) !
391! 'gethml' --- get diagnostic hi, mid, low clouds !
392! !
393! attributes: !
394! language: fortran 90 !
395! machine: ibm-sp, sgi !
396! !
397! !
398! ==================== definition of variables ==================== !
399! !
400! input variables: !
401! plyr (IX,NLAY) : model layer mean pressure in mb (100Pa) !
402! plvl (IX,NLP1) : model level pressure in mb (100Pa) !
403! tlyr (IX,NLAY) : model layer mean temperature in k !
404! tvly (IX,NLAY) : model layer virtual temperature in k !
405! qlyr (IX,NLAY) : layer specific humidity in gm/gm !
406! qstl (IX,NLAY) : layer saturate humidity in gm/gm !
407! rhly (IX,NLAY) : layer relative humidity (=qlyr/qstl) !
408! ccnd (IX,NLAY,ncndl) : layer cloud condensate amount !
409! water, ice, rain, snow (+ graupel) !
410! ncndl : number of layer cloud condensate types (max of 4) !
411! cnvw (IX,NLAY) : layer convective cloud condensate !
412! cnvc (IX,NLAY) : layer convective cloud cover !
413! tracer1 (IX,NLAY,1:ntrac-1) : all tracers (except sphum) !
414! xlat (IX) : grid latitude in radians, default to pi/2 -> -pi/2!
415! range, otherwise see in-line comment !
416! xlon (IX) : grid longitude in radians (not used) !
417! slmsk (IX) : sea/land mask array (sea:0,land:1,sea-ice:2) !
418! dz (ix,nlay) : layer thickness (km) !
419! delp (ix,nlay) : model layer pressure thickness in mb (100Pa) !
420! IX : horizontal dimention !
421! LM,NLAY,NLP1 : vertical layer/level dimensions !
422! deltaq (ix,nlay), half total water distribution width !
423! sup supersaturation !
424! me print control flag !
425! icloud : cloud effect to the optical depth in radiation !
426! kdt : current time step index !
427! ntrac number of tracers (Model%ntrac) !
428! ntcw tracer index for cloud liquid water (Model%ntcw) !
429! ntiw tracer index for cloud ice water (Model%ntiw) !
430! ntrw tracer index for rain water (Model%ntrw) !
431! ntsw tracer index for snow water (Model%ntsw) !
432! ntgl tracer index for graupel (Model%ntgl) !
433! ntclamt tracer index for cloud amount (Model%ntclamt) !
434! imp_physics : cloud microphysics scheme control flag !
435! imp_physics_nssl : NSSL microphysics !
436! imp_physics_fer_hires : Ferrier-Aligo microphysics scheme !
437! imp_physics_gfdl : GFDL microphysics scheme !
438! imp_physics_thompson : Thompson microphysics scheme !
439! imp_physics_wsm6 : WSMG microphysics scheme !
440! imp_physics_zhao_carr : Zhao-Carr microphysics scheme !
441! imp_physics_zhao_carr_pdf : Zhao-Carr microphysics scheme with PDF clouds
442! imp_physics_mg : Morrison-Gettelman microphysics scheme !
443! iovr : choice of cloud-overlap !
444! iovr_rand : flag of cloud-overlap: random (=0) !
445! iovr_maxrand : flag of cloud-overlap: maximum random (=1) !
446! iovr_max : flag of cloud-overlap: maximum (=2) !
447! iovr_dcorr : flag of cloud-overlap: decorrelation length(=3) !
448! iovr_exp : flag of cloud-overlap: exponential (=4) !
449! iovr_exprand : flag of cloud-overlap: exponential random (=5) !
450! idcor : choice for decorrelation-length !
451! idcor_con : flag for decorrelation-length: Use constant value (=0)
452! idcor_hogan : flag for decorrelation-length: (=1) !
453! idcor_oreopoulos: flag for decorrelation-length: (=2) !
454! imfdeepcnv : flag for mass-flux deep convection scheme !
455! imfdeepcnv_gf : flag for scale- & aerosol-aware Grell-Freitas scheme (GSD)
456! imfdeepcnv_c3 : flag for unified convection scheme
457! do_mynnedmf : flag for MYNN-EDMF !
458! lgfdlmprad : flag for GFDLMP radiation interaction !
459! uni_cld : logical - true for cloud fraction from shoc !
460! lmfshal : logical - true for mass flux shallow convection !
461! lmfdeep2 : logical - true for mass flux deep convection !
462! top_at_1 : logical - true if ordered from toa-2-sfc !
463! cldcov : layer cloud fraction (used when uni_cld=.true. !
464! clouds1 : layer total cloud fraction
465! effrl, : effective radius for liquid water
466! effri, : effective radius for ice water
467! effrr, : effective radius for rain water
468! effrs, : effective radius for snow water
469! effr_in, : flag to use effective radii of cloud species in radiation
470! effrl_inout, : eff. radius of cloud liquid water particle
471! effri_inout, : eff. radius of cloud ice water particle
472! effrs_inout : effective radius of cloud snow particle
473! lwp_ex : total liquid water path from explicit microphysics
474! iwp_ex : total ice water path from explicit microphysics
475! lwp_fc : total liquid water path from cloud fraction scheme
476! iwp_fc : total ice water path from cloud fraction scheme
477! dzlay(ix,nlay) : thickness between model layer centers (km) !
478! latdeg(ix) : latitude (in degrees 90 -> -90) !
479! julian : day of the year (fractional julian day) !
480! yearlen : current length of the year (365/366 days) !
481! gridkm : grid length in km !
482! lmfshal : mass-flux shallow conv scheme flag !
483! lmfdeep2 : scale-aware mass-flux deep conv scheme flag !
484! lcrick : control flag for eliminating CRICK !
485! =t: apply layer smoothing to eliminate CRICK !
486! =f: do not apply layer smoothing !
487! lcnorm : control flag for in-cld condensate !
488! =t: normalize cloud condensate !
489! =f: not normalize cloud condensate !
490! !
491! output variables: !
492! cloud profiles: !
493! cld_frac (:,:) - layer total cloud fraction !
494! cld_lwp (:,:) - layer cloud liq water path (g/m**2) !
495! cld_reliq (:,:) - mean eff radius for liq cloud (micron) !
496! cld_iwp (:,:) - layer cloud ice water path (g/m**2) !
497! cld_reice (:,:) - mean eff radius for ice cloud (micron) !
498! cld_rwp (:,:) - layer rain drop water path not assigned !
499! cld_rerain(:,:) - mean eff radius for rain drop (micron) !
500! *** cld_swp (:,:) - layer snow flake water path not assigned !
501! cld_resnow(:,:) - mean eff radius for snow flake (micron) !
502! *** fu's scheme need to be normalized by snow density (g/m**3/1.0e6) !
503! clds (IX,5) : fraction of clouds for low, mid, hi, tot, bl !
504! mtop (IX,3) : vertical indices for low, mid, hi cloud tops !
505! mbot (IX,3) : vertical indices for low, mid, hi cloud bases !
506! de_lgth(ix) : clouds decorrelation length (km) !
507! alpha(ix,nlay) : alpha decorrelation parameter !
508! !
509! ==================== end of description ===================== !
510 implicit none
511
512! --- inputs
513 integer, intent(in) :: ix, lm, nlay, nlp1, me, ncndl, icloud
514 integer, intent(in) :: ntrac, ntcw, ntiw, ntrw, ntsw, ntgl, &
515 & ntclamt
516 integer, intent(in) :: kdt, imfdeepcnv, imfdeepcnv_gf, &
517 & imfdeepcnv_c3
518 integer, intent(in) :: &
519 & imp_physics, ! Flag for MP scheme
520 & imp_physics_nssl, ! Flag for NSSL scheme
521 & imp_physics_fer_hires, ! Flag for fer-hires scheme
522 & imp_physics_gfdl, ! Flag for gfdl scheme
523 & imp_physics_thompson, ! Flag for thompsonscheme
524 & imp_physics_tempo, ! Flag for TEMPO scheme
525 & imp_physics_wsm6, ! Flag for wsm6 scheme
526 & imp_physics_zhao_carr, ! Flag for zhao-carr scheme
527 & imp_physics_zhao_carr_pdf, ! Flag for zhao-carr+PDF scheme
528 & imp_physics_mg ! Flag for MG scheme
529
530 integer, intent(in) :: &
531 & iovr, !
532 & iovr_rand, ! Flag for random cloud overlap method
533 & iovr_maxrand, ! Flag for maximum-random cloud overlap method
534 & iovr_max, ! Flag for maximum cloud overlap method
535 & iovr_dcorr, ! Flag for decorrelation-length cloud overlap method
536 & iovr_exp, ! Flag for exponential cloud overlap method
537 & iovr_exprand, ! Flag for exponential-random cloud overlap method
538 & idcor,
539 & idcor_con,
540 & idcor_hogan,
541 & idcor_oreopoulos
542
543
544 logical, intent(in) :: uni_cld, lmfshal, lmfdeep2, effr_in, &
545 & do_mynnedmf, lgfdlmprad, top_at_1, lcrick, lcnorm, &
546 & xr_cnvcld
547
548 real (kind=kind_phys), dimension(:,:,:), intent(in) :: ccnd, &
549 & tracer1
550 real (kind=kind_phys), dimension(:,:), intent(in) :: plvl, plyr, &
551 & tlyr, tvly, qlyr, qstl, rhly, cnvw, cnvc, cldcov, &
552 & delp, dz, effrl, effri, effrr, effrs, dzlay, clouds1
553
554 real (kind=kind_phys), intent(in) :: sup, dcorr_con, con_ttp, &
555 & con_pi, con_g, con_rd, con_thgni, xr_con, xr_exp
556 real (kind=kind_phys), dimension(:), intent(in) :: xlat, xlon, &
557 & slmsk, si
558
559 real(kind=kind_phys), dimension(:), intent(in) :: latdeg, gridkm
560 real(kind=kind_phys), intent(in) :: julian
561 integer, intent(in) :: yearlen
562
563! --- inout
564 real(kind=kind_phys),dimension(:,:) :: deltaq
565 real(kind=kind_phys),dimension(:,:),intent(inout) :: &
566 & effrl_inout, effri_inout, effrs_inout
567 real(kind=kind_phys), dimension(:), intent(inout) :: &
568 & lwp_ex, iwp_ex, lwp_fc, iwp_fc
569
570! --- outputs
571
572 real (kind=kind_phys), dimension(:,:), intent(out) :: &
573 & cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, &
574 & cld_rwp, cld_rerain, cld_swp, cld_resnow
575 real (kind=kind_phys), dimension(:,:), intent(out) :: clds
576 real (kind=kind_phys), dimension(:), intent(out) :: de_lgth
577 real (kind=kind_phys), dimension(:,:), intent(out) :: alpha
578
579 integer, dimension(:,:), intent(out) :: mtop,mbot
580
581! --- local variables:
582 real (kind=kind_phys), dimension(IX,NLAY) :: cldtot, cldcnv, &
583 & cwp, cip, crp, csp, rew, rei, res, rer, tem2d, clwf
584
585 real (kind=kind_phys) :: ptop1(ix,nk_clds+1), rxlat(ix)
586
587 real (kind=kind_phys) :: clwmin, clwm, clwt, onemrh, value, &
588 & tem1, tem2, tem3
589
590 integer :: i, k, id, nf
591
592!
593!===> ... begin here
594!
595 if (me == 0 .and. kdt == 1) then &
596 print*, 'in radiation_clouds_prop=', imp_physics, uni_cld, &
597 & ncndl, lgfdlmprad, do_mynnedmf, imfdeepcnv, kdt
598 end if
599
600 do k = 1, nlay
601 do i = 1, ix
602 cld_frac(i,k) = 0.0
603 cld_lwp(i,k) = 0.0
604 cld_reliq(i,k) = 0.0
605 cld_iwp(i,k) = 0.0
606 cld_reice(i,k) = 0.0
607 cld_rwp(i,k) = 0.0
608 cld_rerain(i,k) = 0.0
609 cld_swp(i,k) = 0.0
610 cld_resnow(i,k) = 0.0
611 enddo
612 enddo
613
614 do k = 1, nlay
615 do i = 1, ix
616 cldtot(i,k) = 0.0
617 cldcnv(i,k) = 0.0
618 end do
619 end do
620
621 if (imp_physics == imp_physics_zhao_carr .or. &
622 & imp_physics == imp_physics_mg) then ! zhao/moorthi's p
623 ! or unified cloud and/or with MG microphysics
624
625 if (uni_cld .and. ncndl >= 2) then
626 call progclduni (plyr, plvl, tlyr, tvly, ccnd, ncndl, & ! --- inputs
627 & xlat, xlon, slmsk, dz, delp, &
628 & ix, nlay, nlp1, cldcov, &
629 & effrl, effri, effrr, effrs, effr_in, &
630 & dzlay, &
631 & cldtot, cldcnv, lcrick, lcnorm, con_ttp, & ! inout
632 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
633 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
634 & cld_resnow)
635 else
636 call progcld_zhao_carr (plyr ,plvl, tlyr, tvly, qlyr, & ! --- inputs
637 & qstl, rhly, ccnd(1:ix,1:nlay,1), xlat, xlon, &
638 & slmsk, dz, delp, ix, nlay, nlp1, uni_cld, &
639 & lmfshal, lmfdeep2, xr_con, xr_exp, &
640 & cldcov, effrl, effri, effrr, effrs, effr_in, &
641 & dzlay, &
642 & cldtot, cldcnv, lcrick, lcnorm, con_ttp, & ! inout
643 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
644 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
645 & cld_resnow)
646 endif
647
648 elseif(imp_physics == imp_physics_zhao_carr_pdf) then ! zhao/moorthi's prognostic cloud+pdfcld
649
650 call progcld_zhao_carr_pdf (plyr, plvl, tlyr, tvly, qlyr, & ! --- inputs
651 & qstl, rhly, ccnd(1:ix,1:nlay,1), cnvw, cnvc, &
652 & xlat, xlon, slmsk, dz, delp, ix, nlay, nlp1, &
653 & deltaq, sup, kdt, me, dzlay, &
654 & cldtot, cldcnv, lcrick, lcnorm, con_thgni, & ! inout
655 & con_ttp, cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
656 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
657 & cld_resnow)
658
659 elseif (imp_physics == imp_physics_gfdl) then ! GFDL cloud scheme
660
661 if (.not. lgfdlmprad) then
662 call progcld_gfdl_lin (plyr, plvl, tlyr, tvly, qlyr, & ! --- inputs
663 & qstl, rhly, ccnd(1:ix,1:nlay,1), cnvw, cnvc, &
664 & xlat, xlon, slmsk, cldcov, dz, delp, &
665 & ix, nlay, nlp1, dzlay, &
666 & cldtot, cldcnv, lcrick, lcnorm, con_ttp, & ! inout
667 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
668 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
669 & cld_resnow)
670 else
671
672 call progclduni (plyr, plvl, tlyr, tvly, ccnd, ncndl, xlat, & ! --- inputs
673 & xlon, slmsk, dz,delp, ix, nlay, nlp1, cldcov, &
674 & effrl, effri, effrr, effrs, effr_in, &
675 & dzlay, &
676 & cldtot, cldcnv, lcrick, lcnorm, con_ttp, & ! inout
677 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
678 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
679 & cld_resnow)
680 endif
681
682
683 elseif(imp_physics == imp_physics_fer_hires) then
684 if (kdt == 1) then
685 effrl_inout(:,:) = 10.
686 effri_inout(:,:) = 50.
687 effrs_inout(:,:) = 250.
688 endif
689
690 call progcld_fer_hires (plyr,plvl,tlyr,tvly,qlyr,qstl,rhly, & ! --- inputs
691 & tracer1,xlat,xlon,slmsk,dz,delp, &
692 & ntrac-1, ntcw-1,ntiw-1,ntrw-1, &
693 & ix,nlay,nlp1, icloud, xr_con, xr_exp, uni_cld,&
694 & lmfshal, lmfdeep2, &
695 & cldcov(:,1:nlay),effrl_inout(:,:), &
696 & effri_inout(:,:), effrs_inout(:,:), &
697 & dzlay, &
698 & cldtot, cldcnv, lcnorm, & ! inout
699 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
700 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
701 & cld_resnow)
702
703 elseif ( imp_physics == imp_physics_nssl ) then ! NSSL MP
704
705 if(do_mynnedmf .or. imfdeepcnv == imfdeepcnv_gf .or. &
706 & imfdeepcnv == imfdeepcnv_c3) then ! MYNN PBL or GF or unified conv
707 !-- MYNN PBL or convective GF
708 !-- use cloud fractions with SGS clouds
709 do k=1,nlay
710 do i=1,ix
711 cld_frac(i,k) = clouds1(i,k)
712 enddo
713 enddo
714
715 ! --- use clduni with the NSSL microphysics.
716 ! --- make sure that effr_in=.true. in the input.nml!
717 call progclduni (plyr, plvl, tlyr, tvly, ccnd, ncndl, & ! --- inputs
718 & xlat, xlon, slmsk, dz, delp, ix, nlay, nlp1, &
719 & cld_frac, &
720 & effrl, effri, effrr, effrs, effr_in , &
721 & dzlay, &
722 & cldtot, cldcnv, lcrick, lcnorm, con_ttp, & ! inout
723 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
724 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
725 & cld_resnow)
726 else
727 ! MYNN PBL or GF convective are not used
728 call progcld_thompson_wsm6 (plyr,plvl,tlyr,qlyr,qstl, & ! --- inputs
729 & rhly,tracer1,xlat,xlon,slmsk,dz,delp, &
730 & ntrac-1, ntcw-1,ntiw-1,ntrw-1, &
731 & ntsw-1,ntgl-1,con_ttp,xr_cnvcld, &
732 & ix, nlay, nlp1, xr_con, xr_exp, uni_cld, &
733 & lmfshal, lmfdeep2, &
734 & cldcov(:,1:nlay), cnvw, effrl_inout, &
735 & effri_inout, effrs_inout, &
736 & lwp_ex, iwp_ex, lwp_fc, iwp_fc, &
737 & dzlay, &
738 & cldtot, cldcnv, lcnorm, & ! inout
739 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
740 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
741 & cld_resnow)
742 endif ! MYNN PBL or GF
743
744 elseif(imp_physics == imp_physics_thompson &
745 & .or. imp_physics == imp_physics_tempo) then ! Thompson/TEMPO MP
746
747 if(do_mynnedmf .or. imfdeepcnv == imfdeepcnv_gf &
748 & .or. imfdeepcnv == imfdeepcnv_c3) then ! MYNN PBL or GF conv
749
750 if (icloud == 3) then
751 call progcld_thompson (plyr,plvl,tlyr,qlyr,qstl,rhly, & ! --- inputs
752 & tracer1,xlat,xlon,slmsk,dz,delp, &
753 & ntrac-1, ntcw-1,ntiw-1,ntrw-1, &
754 & ntsw-1,ntgl-1, &
755 & ix, lm, nlp1, uni_cld, lmfshal, lmfdeep2, &
756 & cldcov(:,1:lm), effrl, effri, effrs, &
757 & lwp_ex, iwp_ex, lwp_fc, iwp_fc, &
758 & dzlay, gridkm, top_at_1, &
759 & cldtot, cldcnv, & ! inout
760 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
761 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
762 & cld_resnow)
763 else
764
765 !-- MYNN PBL or convective GF
766 !-- use cloud fractions with SGS clouds
767 do k=1,nlay
768 do i=1,ix
769 cld_frac(i,k) = clouds1(i,k)
770 enddo
771 enddo
772
773 ! --- use clduni as with the GFDL microphysics.
774 ! --- make sure that effr_in=.true. in the input.nml!
775 call progclduni (plyr, plvl, tlyr, tvly, ccnd, ncndl, & ! --- inputs
776 & xlat, xlon, slmsk, dz, delp, ix, nlay, nlp1, &
777 & cld_frac, &
778 & effrl, effri, effrr, effrs, effr_in , &
779 & dzlay, &
780 & cldtot, cldcnv, lcrick, lcnorm, con_ttp, & ! inout
781 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
782 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
783 & cld_resnow)
784 endif
785
786 else
787 ! MYNN PBL or GF convective are not used
788
789 if (icloud == 3) then
790 call progcld_thompson (plyr,plvl,tlyr,qlyr,qstl,rhly, & ! --- inputs
791 & tracer1,xlat,xlon,slmsk,dz,delp, &
792 & ntrac-1, ntcw-1,ntiw-1,ntrw-1, &
793 & ntsw-1,ntgl-1, &
794 & ix, lm, nlp1, uni_cld, lmfshal, lmfdeep2, &
795 & cldcov(:,1:lm), effrl, effri, effrs, &
796 & lwp_ex, iwp_ex, lwp_fc, iwp_fc, &
797 & dzlay, gridkm, top_at_1, &
798 & cldtot, cldcnv, & ! inout
799 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
800 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
801 & cld_resnow)
802
803 else
804 call progcld_thompson_wsm6 (plyr,plvl,tlyr,qlyr,qstl, & ! --- inputs
805 & rhly,tracer1,xlat,xlon,slmsk,dz,delp, &
806 & ntrac-1, ntcw-1,ntiw-1,ntrw-1, &
807 & ntsw-1,ntgl-1,con_ttp,xr_cnvcld, &
808 & ix, nlay, nlp1, xr_con, xr_exp, uni_cld, &
809 & lmfshal, lmfdeep2, &
810 & cldcov(:,1:nlay), cnvw, effrl, effri, effrs, &
811 & lwp_ex, iwp_ex, lwp_fc, iwp_fc, &
812 & dzlay, &
813 & cldtot, cldcnv, lcnorm, & ! inout
814 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
815 & cld_reice,cld_rwp, cld_rerain,cld_swp, &
816 & cld_resnow)
817 endif
818 endif ! MYNN PBL or GF
819
820 endif ! end if_imp_physics
821
824! ptopc(k,i): top presure of each cld domain (k=1-4 are sfc,L,m,h;
825! --- i=1,2 are low-lat (<45 degree) and pole regions)
826
827 do i =1, ix
828 rxlat(i) = abs( xlat(i) / con_pi ) ! if xlat in pi/2 -> -pi/2 range
829! rxlat(i) = abs(0.5 - xlat(i)/con_pi) ! if xlat in 0 -> pi range
830 enddo
831
832 do id = 1, 4
833 tem1 = ptopc(id,2) - ptopc(id,1)
834
835 do i =1, ix
836 ptop1(i,id) = ptopc(id,1) + tem1*max( 0.0, 4.0*rxlat(i)-1.0 )
837 enddo
838 enddo
839
840 ! Compute cloud decorrelation length
841 if (idcor == idcor_hogan) then
842 call cmp_dcorr_lgth(ix, xlat, con_pi, de_lgth)
843 endif
844 if (idcor == idcor_oreopoulos) then
845 call cmp_dcorr_lgth(ix, latdeg, julian, yearlen, de_lgth)
846 endif
847 if (idcor == idcor_con) then
848 de_lgth(:) = dcorr_con
849 endif
850
851 ! Call subroutine get_alpha_exper to define alpha parameter for exponential cloud overlap options
852 if ( iovr == iovr_dcorr .or. iovr == iovr_exp &
853 & .or. iovr == iovr_exprand) then
854 call get_alpha_exper(ix, nlay, iovr, iovr_exprand, dzlay, &
855 & de_lgth, cld_frac, alpha)
856 else
857 de_lgth(:) = 0.
858 alpha(:,:) = 0.
859 endif
860
864! --- compute low, mid, high, total, and boundary layer cloud fractions
865! and clouds top/bottom layer indices for low, mid, and high clouds.
866! The three cloud domain boundaries are defined by ptopc. The cloud
867! overlapping method is defined by control flag 'iovr', which may
868! be different for lw and sw radiation programs.
869
870 call gethml &
871! --- inputs:
872 & ( plyr, ptop1, cldtot, cldcnv, dz, de_lgth, alpha, &
873 & ix, nlay, iovr, iovr_rand, iovr_maxrand, iovr_max, &
874 & iovr_dcorr, iovr_exp, iovr_exprand, top_at_1, si, &
875! --- outputs:
876 & clds, mtop, mbot &
877 & )
878
879!...................................
880 end subroutine radiation_clouds_prop
881
885 subroutine progcld_zhao_carr &
886 & ( plyr,plvl,tlyr,tvly,qlyr,qstl,rhly,clw, & ! --- inputs:
887 & xlat,xlon,slmsk,dz,delp, ix, nlay, nlp1, &
888 & uni_cld, lmfshal, lmfdeep2, xr_con, xr_exp, cldcov, &
889 & effrl,effri,effrr,effrs,effr_in, &
890 & dzlay, cldtot, cldcnv, lcrick, lcnorm, con_ttp, &
891 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
892 & cld_reice,cld_rwp, cld_rerain,cld_swp, cld_resnow &
893 & )
894
895! ================= subprogram documentation block ================ !
896! !
897! subprogram: progcld_zhao_carr computes cloud related quantities using !
898! zhao/moorthi's prognostic cloud microphysics scheme. !
899! !
900! abstract: this program computes cloud fractions from cloud !
901! condensates, calculates liquid/ice cloud droplet effective radius, !
902! and computes the low, mid, high, total and boundary layer cloud !
903! fractions and the vertical indices of low, mid, and high cloud !
904! top and base. the three vertical cloud domains are set up in the !
905! initial subroutine "cld_init". !
906! !
907! usage: call progcld_zhao_carr !
908! !
909! attributes: !
910! language: fortran 90 !
911! machine: ibm-sp, sgi !
912! !
913! !
914! ==================== definition of variables ==================== !
915! !
916! input variables: !
917! plyr (IX,NLAY) : model layer mean pressure in mb (100Pa) !
918! plvl (IX,NLP1) : model level pressure in mb (100Pa) !
919! tlyr (IX,NLAY) : model layer mean temperature in k !
920! tvly (IX,NLAY) : model layer virtual temperature in k !
921! qlyr (IX,NLAY) : layer specific humidity in gm/gm !
922! qstl (IX,NLAY) : layer saturate humidity in gm/gm !
923! rhly (IX,NLAY) : layer relative humidity (=qlyr/qstl) !
924! clw (IX,NLAY) : layer cloud condensate amount !
925! xlat (IX) : grid latitude in radians, default to pi/2 -> -pi/2!
926! range, otherwise see in-line comment !
927! xlon (IX) : grid longitude in radians (not used) !
928! slmsk (IX) : sea/land mask array (sea:0,land:1,sea-ice:2) !
929! dz (ix,nlay) : layer thickness (km) !
930! delp (ix,nlay) : model layer pressure thickness in mb (100Pa) !
931! IX : horizontal dimention !
932! NLAY,NLP1 : vertical layer/level dimensions !
933! uni_cld : logical - true for cloud fraction from shoc !
934! lmfshal : logical - true for mass flux shallow convection !
935! lmfdeep2 : logical - true for mass flux deep convection !
936! cldcov : layer cloud fraction (used when uni_cld=.true. !
937! effrl : effective radius for liquid water
938! effri : effective radius for ice water
939! effrr : effective radius for rain water
940! effrs : effective radius for snow water
941! effr_in : logical, if .true. use input effective radii
942! dzlay(ix,nlay) : thickness between model layer centers (km) !
943! lmfshal : mass-flux shallow conv scheme flag !
944! lmfdeep2 : scale-aware mass-flux deep conv scheme flag !
945! lcrick : control flag for eliminating CRICK !
946! =t: apply layer smoothing to eliminate CRICK !
947! =f: do not apply layer smoothing !
948! lcnorm : control flag for in-cld condensate !
949! =t: normalize cloud condensate !
950! =f: not normalize cloud condensate !
951! !
952! output variables: !
953! cloud profiles: !
954! cld_frac (:,:) - layer total cloud fraction !
955! cld_lwp (:,:) - layer cloud liq water path (g/m**2) !
956! cld_reliq (:,:) - mean eff radius for liq cloud (micron) !
957! cld_iwp (:,:) - layer cloud ice water path (g/m**2) !
958! cld_reice (:,:) - mean eff radius for ice cloud (micron) !
959! cld_rwp (:,:) - layer rain drop water path not assigned !
960! cld_rerain(:,:) - mean eff radius for rain drop (micron) !
961! *** cld_swp (:,:) - layer snow flake water path not assigned !
962! cld_resnow(:,:) - mean eff radius for snow flake (micron) !
963! !
964! ==================== end of description ===================== !
965!
966 implicit none
967
968! --- inputs
969 integer, intent(in) :: ix, nlay, nlp1
970
971 logical, intent(in) :: uni_cld, lmfshal, lmfdeep2, effr_in, &
972 & lcrick, lcnorm
973
974 real (kind=kind_phys), dimension(:,:), intent(in) :: plvl, plyr, &
975 & tlyr, tvly, qlyr, qstl, rhly, clw, cldcov, delp, dz, &
976 & effrl, effri, effrr, effrs, dzlay
977
978 real (kind=kind_phys), dimension(:), intent(in) :: xlat, xlon, &
979 & slmsk
980 real (kind=kind_phys), intent(in) :: con_ttp, xr_con, xr_exp
981
982! --- inputs/outputs
983
984 real (kind=kind_phys), dimension(:,:), intent(inout) :: &
985 & cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, &
986 & cld_rwp, cld_rerain, cld_swp, cld_resnow
987
988! --- local variables:
989 real (kind=kind_phys), dimension(IX,NLAY) :: cldtot, cldcnv, &
990 & cwp, cip, crp, csp, rew, rei, res, rer, tem2d, clwf
991
992 real (kind=kind_phys) :: clwmin, clwm, clwt, onemrh, value, &
993 & tem1, tem2, tem3
994
995 integer :: i, k, id, nf
996
997!
998!===> ... begin here
999!
1001 if(effr_in) then
1002 do k = 1, nlay
1003 do i = 1, ix
1004 cldtot(i,k) = 0.0
1005 cldcnv(i,k) = 0.0
1006 cwp(i,k) = 0.0
1007 cip(i,k) = 0.0
1008 crp(i,k) = 0.0
1009 csp(i,k) = 0.0
1010 rew(i,k) = effrl(i,k)
1011 rei(i,k) = effri(i,k)
1012 rer(i,k) = effrr(i,k)
1013 res(i,k) = effrs(i,k)
1014 tem2d(i,k) = min(1.0, max(0.0,(con_ttp-tlyr(i,k))*0.05))
1015 clwf(i,k) = 0.0
1016 enddo
1017 enddo
1018 else
1019 do k = 1, nlay
1020 do i = 1, ix
1021 cldtot(i,k) = 0.0
1022 cldcnv(i,k) = 0.0
1023 cwp(i,k) = 0.0
1024 cip(i,k) = 0.0
1025 crp(i,k) = 0.0
1026 csp(i,k) = 0.0
1027 rew(i,k) = reliq_def ! default liq radius to 10 micron
1028 rei(i,k) = reice_def ! default ice radius to 50 micron
1029 rer(i,k) = rrain_def ! default rain radius to 1000 micron
1030 res(i,k) = rsnow_def ! default snow radius to 250 micron
1031 tem2d(i,k) = min(1.0, max(0.0, (con_ttp-tlyr(i,k))*0.05))
1032 clwf(i,k) = 0.0
1033 enddo
1034 enddo
1035 endif
1036!
1037 if ( lcrick ) then
1038 do i = 1, ix
1039 clwf(i,1) = 0.75*clw(i,1) + 0.25*clw(i,2)
1040 clwf(i,nlay) = 0.75*clw(i,nlay) + 0.25*clw(i,nlay-1)
1041 enddo
1042 do k = 2, nlay-1
1043 do i = 1, ix
1044 clwf(i,k) = 0.25*clw(i,k-1) + 0.5*clw(i,k) + 0.25*clw(i,k+1)
1045 enddo
1046 enddo
1047 else
1048 do k = 1, nlay
1049 do i = 1, ix
1050 clwf(i,k) = clw(i,k)
1051 enddo
1052 enddo
1053 endif
1054
1056
1057 do k = 1, nlay
1058 do i = 1, ix
1059 clwt = max(0.0, clwf(i,k)) * gfac * delp(i,k)
1060 cip(i,k) = clwt * tem2d(i,k)
1061 cwp(i,k) = clwt - cip(i,k)
1062 enddo
1063 enddo
1064
1066
1067 if(.not. effr_in) then
1068 do i = 1, ix
1069 if (nint(slmsk(i)) == 1) then
1070 do k = 1, nlay
1071 rew(i,k) = 5.0 + 5.0 * tem2d(i,k)
1072 enddo
1073 endif
1074 enddo
1075 endif
1076
1077 if (uni_cld) then ! use unified sgs clouds generated outside
1078 do k = 1, nlay
1079 do i = 1, ix
1080 cldtot(i,k) = cldcov(i,k)
1081 enddo
1082 enddo
1083
1084 else
1085
1087
1088
1089 if (.not. lmfshal) then
1091 & ( ix, nlay, xr_con, xr_exp, plyr, clwf, rhly, qstl, & ! --- inputs
1092 & cldtot ) & ! --- outputs
1093 else
1095 & ( ix, nlay, lmfdeep2, xr_con, xr_exp, plyr, clwf, rhly, &
1096 & qstl, & ! --- inputs
1097 & cldtot )
1098 endif
1099
1100 endif ! if (uni_cld) then
1101
1102 do k = 1, nlay
1103 do i = 1, ix
1104 if (cldtot(i,k) < climit) then
1105 cldtot(i,k) = 0.0
1106 cwp(i,k) = 0.0
1107 cip(i,k) = 0.0
1108 crp(i,k) = 0.0
1109 csp(i,k) = 0.0
1110 endif
1111 enddo
1112 enddo
1113
1114 if ( lcnorm ) then
1115 do k = 1, nlay
1116 do i = 1, ix
1117 if (cldtot(i,k) >= climit) then
1118 tem1 = 1.0 / max(climit2, cldtot(i,k))
1119 cwp(i,k) = cwp(i,k) * tem1
1120 cip(i,k) = cip(i,k) * tem1
1121 crp(i,k) = crp(i,k) * tem1
1122 csp(i,k) = csp(i,k) * tem1
1123 endif
1124 enddo
1125 enddo
1126 endif
1127
1130
1131 if(.not.effr_in) then
1132 do k = 1, nlay
1133 do i = 1, ix
1134 tem2 = tlyr(i,k) - con_ttp
1135
1136 if (cip(i,k) > 0.0) then
1137 tem3 = gord * cip(i,k) * plyr(i,k) / (delp(i,k)*tvly(i,k))
1138
1139 if (tem2 < -50.0) then
1140 rei(i,k) = (1250.0/9.917) * tem3 ** 0.109
1141 elseif (tem2 < -40.0) then
1142 rei(i,k) = (1250.0/9.337) * tem3 ** 0.08
1143 elseif (tem2 < -30.0) then
1144 rei(i,k) = (1250.0/9.208) * tem3 ** 0.055
1145 else
1146 rei(i,k) = (1250.0/9.387) * tem3 ** 0.031
1147 endif
1148! rei(i,k) = max(20.0, min(rei(i,k), 300.0))
1149! rei(i,k) = max(10.0, min(rei(i,k), 100.0))
1150 rei(i,k) = max(10.0, min(rei(i,k), 150.0))
1151! rei(i,k) = max(5.0, min(rei(i,k), 130.0))
1152 endif
1153 enddo
1154 enddo
1155 endif
1156
1157!
1158 do k = 1, nlay
1159 do i = 1, ix
1160 cld_frac(i,k) = cldtot(i,k)
1161 cld_lwp(i,k) = cwp(i,k)
1162 cld_reliq(i,k) = rew(i,k)
1163 cld_iwp(i,k) = cip(i,k)
1164 cld_reice(i,k) = rei(i,k)
1165! cld_rwp(i,k) = 0.0
1166 cld_rerain(i,k) = rer(i,k)
1167! cld_swp(i,k) = 0.0
1168 cld_resnow(i,k) = res(i,k)
1169 enddo
1170 enddo
1171!
1172!...................................
1173 end subroutine progcld_zhao_carr
1174!-----------------------------------
1175!-----------------------------------
1176
1181 & ( plyr,plvl,tlyr,tvly,qlyr,qstl,rhly,clw,cnvw,cnvc, & ! --- inputs:
1182 & xlat,xlon,slmsk, dz, delp, &
1183 & ix, nlay, nlp1, &
1184 & deltaq,sup,kdt,me, &
1185 & dzlay, cldtot, cldcnv, lcrick, lcnorm, con_thgni, con_ttp, &
1186 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
1187 & cld_reice,cld_rwp, cld_rerain,cld_swp, cld_resnow &
1188 & )
1189
1190! ================= subprogram documentation block ================ !
1191! !
1192! subprogram: progcld_zhao_carr_pdf computes cloud related quantities using !
1193! zhao/moorthi's prognostic cloud microphysics scheme. !
1194! !
1195! abstract: this program computes cloud fractions from cloud !
1196! condensates, calculates liquid/ice cloud droplet effective radius, !
1197! and computes the low, mid, high, total and boundary layer cloud !
1198! fractions and the vertical indices of low, mid, and high cloud !
1199! top and base. the three vertical cloud domains are set up in the !
1200! initial subroutine "cld_init". !
1201! !
1202! usage: call progcld_zhao_carr_pdf !
1203! !
1204! attributes: !
1205! language: fortran 90 !
1206! machine: ibm-sp, sgi !
1207! !
1208! !
1209! ==================== defination of variables ==================== !
1210! !
1211! input variables: !
1212! plyr (ix,nlay) : model layer mean pressure in mb (100pa) !
1213! plvl (ix,nlp1) : model level pressure in mb (100pa) !
1214! tlyr (ix,nlay) : model layer mean temperature in k !
1215! tvly (ix,nlay) : model layer virtual temperature in k !
1216! qlyr (ix,nlay) : layer specific humidity in gm/gm !
1217! qstl (ix,nlay) : layer saturate humidity in gm/gm !
1218! rhly (ix,nlay) : layer relative humidity (=qlyr/qstl) !
1219! clw (ix,nlay) : layer cloud condensate amount !
1220! xlat (ix) : grid latitude in radians, default to pi/2 -> -pi/2!
1221! range, otherwise see in-line comment !
1222! xlon (ix) : grid longitude in radians (not used) !
1223! slmsk (ix) : sea/land mask array (sea:0,land:1,sea-ice:2) !
1224! dz (ix,nlay) : layer thickness (km) !
1225! delp (ix,nlay) : model layer pressure thickness in mb (100Pa) !
1226! ix : horizontal dimention !
1227! nlay,nlp1 : vertical layer/level dimensions !
1228! cnvw (ix,nlay) : layer convective cloud condensate !
1229! cnvc (ix,nlay) : layer convective cloud cover !
1230! deltaq(ix,nlay) : half total water distribution width !
1231! sup : supersaturation !
1232! dzlay(ix,nlay) : thickness between model layer centers (km) !
1233! lcrick : control flag for eliminating crick !
1234! =t: apply layer smoothing to eliminate crick !
1235! =f: do not apply layer smoothing !
1236! lcnorm : control flag for in-cld condensate !
1237! =t: normalize cloud condensate !
1238! =f: not normalize cloud condensate !
1239! !
1240! output variables: !
1241! cloud profiles: !
1242! cld_frac (:,:) - layer total cloud fraction !
1243! cld_lwp (:,:) - layer cloud liq water path (g/m**2) !
1244! cld_reliq (:,:) - mean eff radius for liq cloud (micron) !
1245! cld_iwp (:,:) - layer cloud ice water path (g/m**2) !
1246! cld_reice (:,:) - mean eff radius for ice cloud (micron) !
1247! cld_rwp (:,:) - layer rain drop water path not assigned !
1248! cld_rerain(:,:) - mean eff radius for rain drop (micron) !
1249! *** cld_swp (:,:) - layer snow flake water path not assigned !
1250! cld_resnow(:,:) - mean eff radius for snow flake (micron) !
1251! !
1252! ==================== end of description ===================== !
1253!
1254 implicit none
1255
1256! --- inputs
1257 integer, intent(in) :: ix, nlay, nlp1,kdt
1258 logical, intent(in) :: lcrick, lcnorm
1259 real (kind=kind_phys), dimension(:,:), intent(in) :: plvl, plyr, &
1260 & tlyr, tvly, qlyr, qstl, rhly, clw, dz, delp, dzlay
1261! & tlyr, tvly, qlyr, qstl, rhly, clw, cnvw, cnvc
1262! real (kind=kind_phys), dimension(:,:), intent(in) :: deltaq
1263 real (kind=kind_phys), intent(in) :: con_thgni, con_ttp
1264 real (kind=kind_phys), dimension(:,:) :: deltaq, cnvw, cnvc
1265 real (kind=kind_phys) qtmp,qsc,rhs
1266 real (kind=kind_phys), intent(in) :: sup
1267 real (kind=kind_phys), parameter :: epsq = 1.0e-12
1268
1269 real (kind=kind_phys), dimension(:), intent(in) :: xlat, xlon, &
1270 & slmsk
1271 integer :: me
1272
1273! --- inputs/outputs
1274
1275 real (kind=kind_phys), dimension(:,:), intent(inout) :: &
1276 & cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, &
1277 & cld_rwp, cld_rerain, cld_swp, cld_resnow
1278
1279! --- local variables:
1280 real (kind=kind_phys), dimension(ix,nlay) :: cldtot, cldcnv, &
1281 & cwp, cip, crp, csp, rew, rei, res, rer, tem2d, clwf
1282
1283 real (kind=kind_phys) :: clwmin, clwm, clwt, onemrh, value, &
1284 & tem1, tem2, tem3
1285
1286 integer :: i, k, id, nf
1287
1288!
1289!===> ... begin here
1290!
1291 do k = 1, nlay
1292 do i = 1, ix
1293 cldtot(i,k) = 0.0
1294 cldcnv(i,k) = 0.0
1295 cwp(i,k) = 0.0
1296 cip(i,k) = 0.0
1297 crp(i,k) = 0.0
1298 csp(i,k) = 0.0
1299 rew(i,k) = reliq_def ! default liq radius to 10 micron
1300 rei(i,k) = reice_def ! default ice radius to 50 micron
1301 rer(i,k) = rrain_def ! default rain radius to 1000 micron
1302 res(i,k) = rsnow_def ! default snow radius to 250 micron
1303 tem2d(i,k) = min( 1.0, max( 0.0, (con_ttp-tlyr(i,k))*0.05 ) )
1304 clwf(i,k) = 0.0
1305 enddo
1306 enddo
1307!
1308 if ( lcrick ) then
1309 do i = 1, ix
1310 clwf(i,1) = 0.75*clw(i,1) + 0.25*clw(i,2)
1311 clwf(i,nlay) = 0.75*clw(i,nlay) + 0.25*clw(i,nlay-1)
1312 enddo
1313 do k = 2, nlay-1
1314 do i = 1, ix
1315 clwf(i,k) = 0.25*clw(i,k-1) + 0.5*clw(i,k) + 0.25*clw(i,k+1)
1316 enddo
1317 enddo
1318 else
1319 do k = 1, nlay
1320 do i = 1, ix
1321 clwf(i,k) = clw(i,k)
1322 enddo
1323 enddo
1324 endif
1325
1326 if(kdt==1) then
1327 do k = 1, nlay
1328 do i = 1, ix
1329 deltaq(i,k) = (1.-0.95)*qstl(i,k)
1330 enddo
1331 enddo
1332 endif
1333
1335
1336 do k = 1, nlay
1337 do i = 1, ix
1338 clwt = max(0.0,(clwf(i,k)+cnvw(i,k))) * gfac * delp(i,k)
1339 cip(i,k) = clwt * tem2d(i,k)
1340 cwp(i,k) = clwt - cip(i,k)
1341 enddo
1342 enddo
1343
1345
1346 do i = 1, ix
1347 if (nint(slmsk(i)) == 1) then
1348 do k = 1, nlay
1349 rew(i,k) = 5.0 + 5.0 * tem2d(i,k)
1350 enddo
1351 endif
1352 enddo
1353
1355
1356 do k = 1, nlay
1357 do i = 1, ix
1358 tem1 = tlyr(i,k) - 273.16
1359 if(tem1 < (con_thgni - 273.16)) then ! for pure ice, has to be consistent with gscond
1360 qsc = sup * qstl(i,k)
1361 rhs = sup
1362 else
1363 qsc = qstl(i,k)
1364 rhs = 1.0
1365 endif
1366 if(rhly(i,k) >= rhs) then
1367 cldtot(i,k) = 1.0
1368 else
1369 qtmp = qlyr(i,k) + clwf(i,k) - qsc
1370 if(deltaq(i,k) > epsq) then
1371! if(qtmp <= -deltaq(i,k) .or. cwmik < epsq) then
1372 if(qtmp <= -deltaq(i,k)) then
1373 cldtot(i,k) = 0.0
1374 elseif(qtmp >= deltaq(i,k)) then
1375 cldtot(i,k) = 1.0
1376 else
1377 cldtot(i,k) = 0.5*qtmp/deltaq(i,k) + 0.5
1378 cldtot(i,k) = max(cldtot(i,k),0.0)
1379 cldtot(i,k) = min(cldtot(i,k),1.0)
1380 endif
1381 else
1382 if(qtmp > 0.) then
1383 cldtot(i,k) = 1.0
1384 else
1385 cldtot(i,k) = 0.0
1386 endif
1387 endif
1388 endif
1389 cldtot(i,k) = cnvc(i,k) + (1-cnvc(i,k))*cldtot(i,k)
1390 cldtot(i,k) = max(cldtot(i,k),0.)
1391 cldtot(i,k) = min(cldtot(i,k),1.)
1392
1393 enddo
1394 enddo
1395
1396 do k = 1, nlay
1397 do i = 1, ix
1398 if (cldtot(i,k) < climit) then
1399 cldtot(i,k) = 0.0
1400 cwp(i,k) = 0.0
1401 cip(i,k) = 0.0
1402 crp(i,k) = 0.0
1403 csp(i,k) = 0.0
1404 endif
1405 enddo
1406 enddo
1407
1408 if ( lcnorm ) then
1409 do k = 1, nlay
1410 do i = 1, ix
1411 if (cldtot(i,k) >= climit) then
1412 tem1 = 1.0 / max(climit2, cldtot(i,k))
1413 cwp(i,k) = cwp(i,k) * tem1
1414 cip(i,k) = cip(i,k) * tem1
1415 crp(i,k) = crp(i,k) * tem1
1416 csp(i,k) = csp(i,k) * tem1
1417 endif
1418 enddo
1419 enddo
1420 endif
1421
1424
1425 do k = 1, nlay
1426 do i = 1, ix
1427 tem2 = tlyr(i,k) - con_ttp
1428
1429 if (cip(i,k) > 0.0) then
1430! tem3 = gord * cip(i,k) * (plyr(i,k)/delp(i,k)) / tvly(i,k)
1431 tem3 = gord * cip(i,k) * plyr(i,k) / (delp(i,k)*tvly(i,k))
1432
1433 if (tem2 < -50.0) then
1434 rei(i,k) = (1250.0/9.917) * tem3 ** 0.109
1435 elseif (tem2 < -40.0) then
1436 rei(i,k) = (1250.0/9.337) * tem3 ** 0.08
1437 elseif (tem2 < -30.0) then
1438 rei(i,k) = (1250.0/9.208) * tem3 ** 0.055
1439 else
1440 rei(i,k) = (1250.0/9.387) * tem3 ** 0.031
1441 endif
1442! rei(i,k) = max(20.0, min(rei(i,k), 300.0))
1443! rei(i,k) = max(10.0, min(rei(i,k), 100.0))
1444 rei(i,k) = max(10.0, min(rei(i,k), 150.0))
1445! rei(i,k) = max(5.0, min(rei(i,k), 130.0))
1446 endif
1447 enddo
1448 enddo
1449
1450!
1451 do k = 1, nlay
1452 do i = 1, ix
1453 cld_frac(i,k) = cldtot(i,k)
1454 cld_lwp(i,k) = cwp(i,k)
1455 cld_reliq(i,k) = rew(i,k)
1456 cld_iwp(i,k) = cip(i,k)
1457 cld_reice(i,k) = rei(i,k)
1458! cld_rwp(i,k) = 0.0
1459 cld_rerain(i,k) = rer(i,k)
1460! cld_swp(i,k) = 0.0
1461 cld_resnow(i,k) = res(i,k)
1462 enddo
1463 enddo
1464!
1465!...................................
1466 end subroutine progcld_zhao_carr_pdf
1467!-----------------------------------
1468
1469
1470!-----------------------------------
1474 subroutine progcld_gfdl_lin &
1475 & ( plyr,plvl,tlyr,tvly,qlyr,qstl,rhly,clw,cnvw,cnvc, & ! --- inputs:
1476 & xlat,xlon,slmsk,cldtot, dz, delp, &
1477 & ix, nlay, nlp1, &
1478 & dzlay, cldtot1, cldcnv, lcrick, lcnorm, con_ttp, &
1479 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
1480 & cld_reice,cld_rwp, cld_rerain,cld_swp, cld_resnow &
1481 & )
1482
1483! ================= subprogram documentation block ================ !
1484! !
1485! subprogram: progcld_gfdl_lin computes cloud related quantities using !
1486! GFDL Lin MP prognostic cloud microphysics scheme. !
1487! !
1488! abstract: this program computes cloud fractions from cloud !
1489! condensates, calculates liquid/ice cloud droplet effective radius, !
1490! and computes the low, mid, high, total and boundary layer cloud !
1491! fractions and the vertical indices of low, mid, and high cloud !
1492! top and base. the three vertical cloud domains are set up in the !
1493! initial subroutine "cld_init". !
1494! !
1495! usage: call progcld_gfdl_lin !
1496! !
1497! attributes: !
1498! language: fortran 90 !
1499! machine: ibm-sp, sgi !
1500! !
1501! !
1502! ==================== definition of variables ==================== !
1503! !
1504! input variables: !
1505! plyr (IX,NLAY) : model layer mean pressure in mb (100Pa) !
1506! plvl (IX,NLP1) : model level pressure in mb (100Pa) !
1507! tlyr (IX,NLAY) : model layer mean temperature in k !
1508! tvly (IX,NLAY) : model layer virtual temperature in k !
1509! qlyr (IX,NLAY) : layer specific humidity in gm/gm !
1510! qstl (IX,NLAY) : layer saturate humidity in gm/gm !
1511! rhly (IX,NLAY) : layer relative humidity (=qlyr/qstl) !
1512! clw (IX,NLAY) : layer cloud condensate amount !
1513! cnvw (IX,NLAY) : layer convective cloud condensate !
1514! cnvc (IX,NLAY) : layer convective cloud cover !
1515! xlat (IX) : grid latitude in radians, default to pi/2 -> -pi/2!
1516! range, otherwise see in-line comment !
1517! xlon (IX) : grid longitude in radians (not used) !
1518! slmsk (IX) : sea/land mask array (sea:0,land:1,sea-ice:2) !
1519! dz (ix,nlay) : layer thickness (km) !
1520! delp (ix,nlay) : model layer pressure thickness in mb (100Pa) !
1521! IX : horizontal dimention !
1522! NLAY,NLP1 : vertical layer/level dimensions !
1523! dzlay(ix,nlay) : thickness between model layer centers (km) !
1524! lcrick : control flag for eliminating CRICK !
1525! =t: apply layer smoothing to eliminate CRICK !
1526! =f: do not apply layer smoothing !
1527! lcnorm : control flag for in-cld condensate !
1528! =t: normalize cloud condensate !
1529! =f: not normalize cloud condensate !
1530! !
1531! output variables: !
1532! cloud profiles: !
1533! cld_frac (:,:) - layer total cloud fraction !
1534! cld_lwp (:,:) - layer cloud liq water path (g/m**2) !
1535! cld_reliq (:,:) - mean eff radius for liq cloud (micron) !
1536! cld_iwp (:,:) - layer cloud ice water path (g/m**2) !
1537! cld_reice (:,:) - mean eff radius for ice cloud (micron) !
1538! cld_rwp (:,:) - layer rain drop water path not assigned !
1539! cld_rerain(:,:) - mean eff radius for rain drop (micron) !
1540! *** cld_swp (:,:) - layer snow flake water path not assigned !
1541! cld_resnow(:,:) - mean eff radius for snow flake (micron) !
1542! !
1543! ==================== end of description ===================== !
1544!
1545 implicit none
1546
1547! --- inputs
1548 integer, intent(in) :: ix, nlay, nlp1
1549 logical, intent(in) :: lcrick, lcnorm
1550 real (kind=kind_phys), dimension(:,:), intent(in) :: plvl, plyr, &
1551 & tlyr, tvly, qlyr, qstl, rhly, clw, cldtot, cnvw, cnvc, &
1552 & delp, dz, dzlay
1553 real (kind=kind_phys) :: con_ttp
1554
1555 real (kind=kind_phys), dimension(:), intent(in) :: xlat, xlon, &
1556 & slmsk
1557
1558 real (kind=kind_phys), dimension(:,:), intent(inout) :: cldtot1
1559
1560! --- inputs/outputs
1561
1562 real (kind=kind_phys), dimension(:,:), intent(inout) :: &
1563 & cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, &
1564 & cld_rwp, cld_rerain, cld_swp, cld_resnow
1565
1566! --- local variables:
1567 real (kind=kind_phys), dimension(IX,NLAY) :: cldcnv, &
1568 & cwp, cip, crp, csp, rew, rei, res, rer, tem2d, clwf
1569
1570 real (kind=kind_phys) :: clwmin, clwm, clwt, onemrh, value, &
1571 & tem1, tem2, tem3
1572
1573 integer :: i, k, id, nf
1574
1575!
1576!===> ... begin here
1577!
1579 do k = 1, nlay
1580 do i = 1, ix
1581 cldcnv(i,k) = 0.0
1582 cwp(i,k) = 0.0
1583 cip(i,k) = 0.0
1584 crp(i,k) = 0.0
1585 csp(i,k) = 0.0
1586 rew(i,k) = reliq_def
1587 rei(i,k) = reice_def
1588 rer(i,k) = rrain_def
1589 res(i,k) = rsnow_def
1590 tem2d(i,k) = min( 1.0, max( 0.0, (con_ttp-tlyr(i,k))*0.05 ) )
1591 clwf(i,k) = 0.0
1592 enddo
1593 enddo
1594!
1595 if ( lcrick ) then
1596 do i = 1, ix
1597 clwf(i,1) = 0.75*clw(i,1) + 0.25*clw(i,2)
1598 clwf(i,nlay) = 0.75*clw(i,nlay) + 0.25*clw(i,nlay-1)
1599 enddo
1600 do k = 2, nlay-1
1601 do i = 1, ix
1602 clwf(i,k) = 0.25*clw(i,k-1) + 0.5*clw(i,k) + 0.25*clw(i,k+1)
1603 enddo
1604 enddo
1605 else
1606 do k = 1, nlay
1607 do i = 1, ix
1608 clwf(i,k) = clw(i,k)
1609 enddo
1610 enddo
1611 endif
1612
1614
1615 do k = 1, nlay
1616 do i = 1, ix
1617 clwt = max(0.0,(clwf(i,k)+cnvw(i,k))) * gfac * delp(i,k)
1618 cip(i,k) = clwt * tem2d(i,k)
1619 cwp(i,k) = clwt - cip(i,k)
1620 enddo
1621 enddo
1622
1624
1625 do i = 1, ix
1626 if (nint(slmsk(i)) == 1) then
1627 do k = 1, nlay
1628 rew(i,k) = 5.0 + 5.0 * tem2d(i,k)
1629 enddo
1630 endif
1631 enddo
1632
1633 do k = 1, nlay
1634 do i = 1, ix
1635 if (cldtot(i,k) < climit) then
1636 cwp(i,k) = 0.0
1637 cip(i,k) = 0.0
1638 crp(i,k) = 0.0
1639 csp(i,k) = 0.0
1640 endif
1641 enddo
1642 enddo
1643
1644 if ( lcnorm ) then
1645 do k = 1, nlay
1646 do i = 1, ix
1647 if (cldtot(i,k) >= climit) then
1648 tem1 = 1.0 / max(climit2, cldtot(i,k))
1649 cwp(i,k) = cwp(i,k) * tem1
1650 cip(i,k) = cip(i,k) * tem1
1651 crp(i,k) = crp(i,k) * tem1
1652 csp(i,k) = csp(i,k) * tem1
1653 endif
1654 enddo
1655 enddo
1656 endif
1657
1660
1661 do k = 1, nlay
1662 do i = 1, ix
1663 tem2 = tlyr(i,k) - con_ttp
1664
1665 if (cip(i,k) > 0.0) then
1666 tem3 = gord * cip(i,k) * plyr(i,k) / (delp(i,k)*tvly(i,k))
1667
1668 if (tem2 < -50.0) then
1669 rei(i,k) = (1250.0/9.917) * tem3 ** 0.109
1670 elseif (tem2 < -40.0) then
1671 rei(i,k) = (1250.0/9.337) * tem3 ** 0.08
1672 elseif (tem2 < -30.0) then
1673 rei(i,k) = (1250.0/9.208) * tem3 ** 0.055
1674 else
1675 rei(i,k) = (1250.0/9.387) * tem3 ** 0.031
1676 endif
1677! rei(i,k) = max(20.0, min(rei(i,k), 300.0))
1678! rei(i,k) = max(10.0, min(rei(i,k), 100.0))
1679 rei(i,k) = max(10.0, min(rei(i,k), 150.0))
1680! rei(i,k) = max(5.0, min(rei(i,k), 130.0))
1681 endif
1682 enddo
1683 enddo
1684
1685 do k = 1, nlay
1686 do i = 1, ix
1687 cldtot1(i,k) = cldtot(i,k)
1688 enddo
1689 enddo
1690
1691!
1692 do k = 1, nlay
1693 do i = 1, ix
1694 cld_frac(i,k) = cldtot(i,k)
1695 cld_lwp(i,k) = cwp(i,k)
1696 cld_reliq(i,k) = rew(i,k)
1697 cld_iwp(i,k) = cip(i,k)
1698 cld_reice(i,k) = rei(i,k)
1699! cld_rwp(i,k) = 0.0
1700 cld_rerain(i,k) = rer(i,k)
1701! cld_swp(i,k) = 0.0
1702 cld_resnow(i,k) = res(i,k)
1703 enddo
1704 enddo
1705!
1706!...................................
1707 end subroutine progcld_gfdl_lin
1708!-----------------------------------
1709
1710!-----------------------------------
1714 & ( plyr,plvl,tlyr,tvly,qlyr,qstl,rhly,clw, & ! --- inputs:
1715 & xlat,xlon,slmsk,dz,delp, &
1716 & ntrac,ntcw,ntiw,ntrw, &
1717 & ix, nlay, nlp1, icloud, xr_con, xr_exp, &
1718 & uni_cld, lmfshal, lmfdeep2, cldcov, &
1719 & re_cloud,re_ice,re_snow, &
1720 & dzlay, cldtot, cldcnv, lcnorm, &
1721 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
1722 & cld_reice,cld_rwp, cld_rerain,cld_swp, cld_resnow &
1723 & )
1724
1725! ================= subprogram documentation block ================ !
1726! !
1727! subprogram: progcld_fer_hires computes cloud related quantities using !
1728! Ferrier-Aligo cloud microphysics scheme. !
1729! !
1730! abstract: this program computes cloud fractions from cloud !
1731! condensates, !
1732! and computes the low, mid, high, total and boundary layer cloud !
1733! fractions and the vertical indices of low, mid, and high cloud !
1734! top and base. the three vertical cloud domains are set up in the !
1735! initial subroutine "cld_init". !
1736! !
1737! usage: call progcld_fer_hires !
1738! !
1739! attributes: !
1740! language: fortran 90 !
1741! machine: ibm-sp, sgi !
1742! !
1743! !
1744! ==================== definition of variables ==================== !
1745! !
1746! input variables: !
1747! plyr (IX,NLAY) : model layer mean pressure in mb (100Pa) !
1748! plvl (IX,NLP1) : model level pressure in mb (100Pa) !
1749! tlyr (IX,NLAY) : model layer mean temperature in k !
1750! tvly (IX,NLAY) : model layer virtual temperature in k !
1751! qlyr (IX,NLAY) : layer specific humidity in gm/gm !
1752! qstl (IX,NLAY) : layer saturate humidity in gm/gm !
1753! rhly (IX,NLAY) : layer relative humidity (=qlyr/qstl) !
1754! clw (IX,NLAY,ntrac) : layer cloud condensate amount !
1755! xlat (IX) : grid latitude in radians, default to pi/2 -> -pi/2!
1756! range, otherwise see in-line comment !
1757! xlon (IX) : grid longitude in radians (not used) !
1758! slmsk (IX) : sea/land mask array (sea:0,land:1,sea-ice:2) !
1759! dz (ix,nlay) : layer thickness (km) !
1760! delp (ix,nlay) : model layer pressure thickness in mb (100Pa) !
1761! IX : horizontal dimention !
1762! NLAY,NLP1 : vertical layer/level dimensions !
1763! icloud : cloud effect to the optical depth in radiation !
1764! uni_cld : logical - true for cloud fraction from shoc !
1765! lmfshal : logical - true for mass flux shallow convection !
1766! lmfdeep2 : logical - true for mass flux deep convection !
1767! cldcov : layer cloud fraction (used when uni_cld=.true. !
1768! dzlay(ix,nlay) : thickness between model layer centers (km) !
1769! lmfshal : mass-flux shallow conv scheme flag !
1770! lmfdeep2 : scale-aware mass-flux deep conv scheme flag !
1771! lcrick : control flag for eliminating CRICK !
1772! =t: apply layer smoothing to eliminate CRICK !
1773! =f: do not apply layer smoothing !
1774! lcnorm : control flag for in-cld condensate !
1775! =t: normalize cloud condensate !
1776! =f: not normalize cloud condensate !
1777! !
1778! output variables: !
1779! cloud profiles: !
1780! cld_frac (:,:) - layer total cloud fraction !
1781! cld_lwp (:,:) - layer cloud liq water path (g/m**2) !
1782! cld_reliq (:,:) - mean eff radius for liq cloud (micron) !
1783! cld_iwp (:,:) - layer cloud ice water path (g/m**2) !
1784! cld_reice (:,:) - mean eff radius for ice cloud (micron) !
1785! cld_rwp (:,:) - layer rain drop water path not assigned !
1786! cld_rerain(:,:) - mean eff radius for rain drop (micron) !
1787! *** cld_swp (:,:) - layer snow flake water path not assigned !
1788! cld_resnow(:,:) - mean eff radius for snow flake (micron) !
1789! !
1790! ==================== end of description ===================== !
1791!
1792 implicit none
1793
1794! --- inputs
1795 integer, intent(in) :: ix, nlay, nlp1, icloud
1796 integer, intent(in) :: ntrac, ntcw, ntiw, ntrw
1797
1798 logical, intent(in) :: uni_cld, lmfshal, lmfdeep2, lcnorm
1799
1800 real (kind=kind_phys), intent(in) :: xr_con, xr_exp
1801 real (kind=kind_phys), dimension(:,:), intent(in) :: plvl, plyr, &
1802 & tlyr, tvly, qlyr, qstl, rhly, cldcov, delp, dz, dzlay
1803
1804 real (kind=kind_phys), dimension(:,:), intent(inout) :: &
1805 & re_cloud, re_ice, re_snow
1806
1807 real (kind=kind_phys), dimension(:,:,:), intent(in) :: clw
1808
1809 real (kind=kind_phys), dimension(:), intent(in) :: xlat, xlon, &
1810 & slmsk
1811
1812! --- inputs/outputs
1813
1814 real (kind=kind_phys), dimension(:,:), intent(inout) :: &
1815 & cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, &
1816 & cld_rwp, cld_rerain, cld_swp, cld_resnow
1817
1818! --- local variables:
1819 real (kind=kind_phys), dimension(IX,NLAY) :: cldtot, cldcnv, &
1820 & cwp, cip, crp, csp, rew, rei, res, rer, tem2d, clwf
1821
1822 real (kind=kind_phys) :: clwmin, clwm, clwt, onemrh, value, &
1823 & tem1, tem2, tem3
1824
1825 integer :: i, k, id, nf
1826
1827!
1828!===> ... begin here
1829!
1830 do k = 1, nlay
1831 do i = 1, ix
1832 cldtot(i,k) = 0.0
1833 cldcnv(i,k) = 0.0
1834 cwp(i,k) = 0.0
1835 cip(i,k) = 0.0
1836 crp(i,k) = 0.0
1837 csp(i,k) = 0.0
1838 rew(i,k) = re_cloud(i,k)
1839 rei(i,k) = re_ice(i,k)
1840 rer(i,k) = rrain_def ! default rain radius to 1000 micron
1841 res(i,k) = re_snow(i,k)
1842! tem2d (i,k) = min( 1.0, max( 0.0, (con_ttp-tlyr(i,k))*0.05 ) )
1843 clwf(i,k) = 0.0
1844 enddo
1845 enddo
1846!
1847! if ( lcrick ) then
1848! do i = 1, IX
1849! clwf(i,1) = 0.75*clw(i,1) + 0.25*clw(i,2)
1850! clwf(i,nlay) = 0.75*clw(i,nlay) + 0.25*clw(i,nlay-1)
1851! enddo
1852! do k = 2, NLAY-1
1853! do i = 1, IX
1854! clwf(i,K) = 0.25*clw(i,k-1) + 0.5*clw(i,k) + 0.25*clw(i,k+1)
1855! enddo
1856! enddo
1857! else
1858! do k = 1, NLAY
1859! do i = 1, IX
1860! clwf(i,k) = clw(i,k)
1861! enddo
1862! enddo
1863! endif
1864
1865 do k = 1, nlay
1866 do i = 1, ix
1867 clwf(i,k) = clw(i,k,ntcw) + clw(i,k,ntiw)
1868 enddo
1869 enddo
1870
1872
1873 do k = 1, nlay
1874 do i = 1, ix
1875 cwp(i,k) = max(0.0, clw(i,k,ntcw) * gfac * delp(i,k))
1876 cip(i,k) = max(0.0, clw(i,k,ntiw) * gfac * delp(i,k))
1877 crp(i,k) = max(0.0, clw(i,k,ntrw) * gfac * delp(i,k))
1878 csp(i,k) = 0.0
1879 enddo
1880 enddo
1881
1882!mz* if (uni_cld) then ! use unified sgs clouds generated outside
1883!mz* use unified sgs clouds generated outside
1884 if (uni_cld) then
1885 do k = 1, nlay
1886 do i = 1, ix
1887 cldtot(i,k) = cldcov(i,k)
1888 enddo
1889 enddo
1890
1891 else
1892
1894
1895 if (.not. lmfshal) then
1897 & ( ix, nlay, xr_con, xr_exp, plyr, clwf, rhly, qstl, & ! --- inputs
1898 & cldtot ) & ! --- outputs
1899 else
1901 & ( ix, nlay, lmfdeep2, xr_con, xr_exp, plyr, clwf, rhly, &
1902 & qstl, & ! --- inputs
1903 & cldtot ) & ! --- outputs
1904 endif
1905
1906 endif ! if (uni_cld) then
1907
1908 do k = 1, nlay
1909 do i = 1, ix
1910 if (cldtot(i,k) < climit) then
1911 cldtot(i,k) = 0.0
1912 cwp(i,k) = 0.0
1913 cip(i,k) = 0.0
1914 crp(i,k) = 0.0
1915 csp(i,k) = 0.0
1916 endif
1917 enddo
1918 enddo
1919
1920 if ( lcnorm ) then
1921 do k = 1, nlay
1922 do i = 1, ix
1923 if (cldtot(i,k) >= climit) then
1924 tem1 = 1.0 / max(climit2, cldtot(i,k))
1925 cwp(i,k) = cwp(i,k) * tem1
1926 cip(i,k) = cip(i,k) * tem1
1927 crp(i,k) = crp(i,k) * tem1
1928 csp(i,k) = csp(i,k) * tem1
1929 endif
1930 enddo
1931 enddo
1932 endif
1933!
1934 do k = 1, nlay
1935 do i = 1, ix
1936 cld_frac(i,k) = cldtot(i,k)
1937 cld_lwp(i,k) = cwp(i,k)
1938 cld_reliq(i,k) = rew(i,k)
1939 cld_iwp(i,k) = cip(i,k)
1940 cld_reice(i,k) = rei(i,k)
1941 cld_rwp(i,k) = crp(i,k)
1942 cld_rerain(i,k) = rer(i,k)
1943 cld_swp(i,k) = 0.0
1944 cld_resnow(i,k) = 10.0
1945 re_cloud(i,k) = rew(i,k)
1946 re_ice(i,k) = rei(i,k)
1947 re_snow(i,k) = 10.
1948 enddo
1949 enddo
1950!
1951!...................................
1952 end subroutine progcld_fer_hires
1953!...................................
1954
1955
1958 & ( plyr,plvl,tlyr,qlyr,qstl,rhly,clw, & ! --- inputs:
1959 & xlat,xlon,slmsk,dz,delp, &
1960 & ntrac,ntcw,ntiw,ntrw,ntsw,ntgl,con_ttp, &
1961 & xr_cnvcld, ix, nlay, nlp1, xr_con, xr_exp, &
1962 & uni_cld, lmfshal, lmfdeep2, cldcov, cnvw, &
1963 & re_cloud,re_ice,re_snow, &
1964 & lwp_ex, iwp_ex, lwp_fc, iwp_fc, &
1965 & dzlay, cldtot, cldcnv, lcnorm, &
1966 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
1967 & cld_reice,cld_rwp, cld_rerain,cld_swp, cld_resnow &
1968 & )
1969
1970! ================= subprogram documentation block ================ !
1971! !
1972! subprogram: progcld_thompson_wsm6 !
1973! computes cloud related quantities using !
1974! Thompson/WSM6/NSSL cloud microphysics scheme. !
1975! !
1976! abstract: this program computes cloud fractions from cloud !
1977! condensates, !
1978! and computes the low, mid, high, total and boundary layer cloud !
1979! fractions and the vertical indices of low, mid, and high cloud !
1980! top and base. the three vertical cloud domains are set up in the !
1981! initial subroutine "cld_init". !
1982! !
1983! usage: call progcld_thompson_wsm6 !
1984! !
1985! attributes: !
1986! language: fortran 90 !
1987! machine: ibm-sp, sgi !
1988! !
1989! !
1990! ==================== definition of variables ==================== !
1991! !
1992! !
1993! input variables: !
1994! plyr (IX,NLAY) : model layer mean pressure in mb (100Pa) !
1995! plvl (IX,NLP1) : model level pressure in mb (100Pa) !
1996! tlyr (IX,NLAY) : model layer mean temperature in k !
1997! tvly (IX,NLAY) : model layer virtual temperature in k !
1998! qlyr (IX,NLAY) : layer specific humidity in gm/gm !
1999! qstl (IX,NLAY) : layer saturate humidity in gm/gm !
2000! rhly (IX,NLAY) : layer relative humidity (=qlyr/qstl) !
2001! clw (IX,NLAY,ntrac) : layer cloud condensate amount !
2002! xlat (IX) : grid latitude in radians, default to pi/2 -> -pi/2!
2003! range, otherwise see in-line comment !
2004! xlon (IX) : grid longitude in radians (not used) !
2005! slmsk (IX) : sea/land mask array (sea:0,land:1,sea-ice:2) !
2006! dz (ix,nlay) : layer thickness (km) !
2007! delp (ix,nlay) : model layer pressure thickness in mb (100Pa) !
2008! IX : horizontal dimention !
2009! NLAY,NLP1 : vertical layer/level dimensions !
2010! uni_cld : logical - true for cloud fraction from shoc !
2011! lmfshal : logical - true for mass flux shallow convection !
2012! lmfdeep2 : logical - true for mass flux deep convection !
2013! cldcov : layer cloud fraction (used when uni_cld=.true. !
2014! lmfshal : mass-flux shallow conv scheme flag !
2015! lmfdeep2 : scale-aware mass-flux deep conv scheme flag !
2016! lcrick : control flag for eliminating CRICK !
2017! =t: apply layer smoothing to eliminate CRICK !
2018! =f: do not apply layer smoothing !
2019! lcnorm : control flag for in-cld condensate !
2020! =t: normalize cloud condensate !
2021! =f: not normalize cloud condensate !
2022! !
2023! output variables: !
2024! cloud profiles: !
2025! cld_frac (:,:) - layer total cloud fraction !
2026! cld_lwp (:,:) - layer cloud liq water path (g/m**2) !
2027! cld_reliq (:,:) - mean eff radius for liq cloud (micron) !
2028! cld_iwp (:,:) - layer cloud ice water path (g/m**2) !
2029! cld_reice (:,:) - mean eff radius for ice cloud (micron) !
2030! cld_rwp (:,:) - layer rain drop water path not assigned !
2031! cld_rerain(:,:) - mean eff radius for rain drop (micron) !
2032! *** cld_swp (:,:) - layer snow flake water path not assigned !
2033! cld_resnow(:,:) - mean eff radius for snow flake (micron) !
2034! *** fu's scheme need to be normalized by snow density (g/m**3/1.0e6) !
2035! clds (IX,5) : fraction of clouds for low, mid, hi, tot, bl !
2036! mtop (IX,3) : vertical indices for low, mid, hi cloud tops !
2037! mbot (IX,3) : vertical indices for low, mid, hi cloud bases !
2038! de_lgth(ix) : clouds decorrelation length (km) !
2039! !
2040! ==================== end of description ===================== !
2041!
2042 implicit none
2043
2044! --- inputs
2045 integer, intent(in) :: ix, nlay, nlp1
2046 integer, intent(in) :: ntrac, ntcw, ntiw, ntrw, ntsw, ntgl
2047
2048 logical, intent(in) :: uni_cld, lmfshal, lmfdeep2, lcnorm, &
2049 & xr_cnvcld
2050
2051 real (kind=kind_phys), dimension(:,:), intent(in) :: plvl, plyr, &
2052 & tlyr, qlyr, qstl, rhly, cldcov, delp, dz, dzlay, &
2053 & re_cloud, re_ice, re_snow, cnvw
2054 real (kind=kind_phys), dimension(:), intent(inout) :: &
2055 & lwp_ex, iwp_ex, lwp_fc, iwp_fc
2056
2057 real (kind=kind_phys), dimension(:,:,:), intent(in) :: clw
2058
2059 real (kind=kind_phys), dimension(:), intent(in) :: xlat, xlon, &
2060 & slmsk
2061 real (kind=kind_phys), intent(in) :: con_ttp, xr_con, xr_exp
2062! --- inputs/outputs
2063
2064 real (kind=kind_phys), dimension(:,:), intent(inout) :: &
2065 & cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, &
2066 & cld_rwp, cld_rerain, cld_swp, cld_resnow
2067
2068! --- local variables:
2069 real (kind=kind_phys), dimension(IX,NLAY) :: cldtot, cldcnv, &
2070 & cwp, cip, crp, csp, rew, rei, res, rer, tem2d, clwf
2071
2072 real (kind=kind_phys) :: clwmin, clwm, clwt, onemrh, value, &
2073 & tem1, tem2, tem3
2074
2075 integer :: i, k, id, nf
2076
2077! --- constant values
2078 real (kind=kind_phys), parameter :: snow2ice = 0.25
2079 real (kind=kind_phys), parameter :: coef_t = 0.025
2080!
2081!===> ... begin here
2082
2083 do k = 1, nlay
2084 do i = 1, ix
2085 cldtot(i,k) = 0.0
2086 cldcnv(i,k) = 0.0
2087 cwp(i,k) = 0.0
2088 cip(i,k) = 0.0
2089 crp(i,k) = 0.0
2090 csp(i,k) = 0.0
2091 rew(i,k) = re_cloud(i,k)
2092 rei(i,k) = re_ice(i,k)
2093 rer(i,k) = rrain_def ! default rain radius to 1000 micron
2094 res(i,k) = re_snow(i,k)
2095 tem2d(i,k) = min(1.0, max( 0.0, (con_ttp-tlyr(i,k))*coef_t))
2096 clwf(i,k) = 0.0
2097 enddo
2098 enddo
2099!
2100!
2101! if ( lcrick ) then
2102! do i = 1, IX
2103! clwf(i,1) = 0.75*clw(i,1) + 0.25*clw(i,2)
2104! clwf(i,nlay) = 0.75*clw(i,nlay) + 0.25*clw(i,nlay-1)
2105! enddo
2106! do k = 2, NLAY-1
2107! do i = 1, IX
2108! clwf(i,K) = 0.25*clw(i,k-1) + 0.5*clw(i,k) + 0.25*clw(i,k+1)
2109! enddo
2110! enddo
2111! else
2112! do k = 1, NLAY
2113! do i = 1, IX
2114! clwf(i,k) = clw(i,k)
2115! enddo
2116! enddo
2117! endif
2118
2121
2122 if(xr_cnvcld)then
2123 do k = 1, nlay
2124 do i = 1, ix
2125 clwf(i,k) = clw(i,k,ntcw) + clw(i,k,ntiw) + clw(i,k,ntsw)
2126 & + clw(i,k,ntrw) + cnvw(i,k)
2127 enddo
2128 enddo
2129 else
2130 do k = 1, nlay
2131 do i = 1, ix
2132 clwf(i,k) = clw(i,k,ntcw) + clw(i,k,ntiw) + clw(i,k,ntsw)
2133 & + clw(i,k,ntrw)
2134 enddo
2135 enddo
2136 endif
2137
2140 do k = 1, nlay-1
2141 do i = 1, ix
2142 if(xr_cnvcld)then
2143 tem1 = cnvw(i,k)*(1.-tem2d(i,k))
2144 else
2145 tem1 = 0.
2146 endif
2147 cwp(i,k) = max(0.0, (clw(i,k,ntcw)+tem1) *
2148 & gfac * delp(i,k))
2149 if(tem1 > 1.e-12 .and. clw(i,k,ntcw) < 1.e-12)
2150 & rew(i,k)=reliq_def
2151 if(xr_cnvcld)then
2152 tem2 = cnvw(i,k)*tem2d(i,k)
2153 else
2154 tem2 = 0.
2155 endif
2156 cip(i,k) = max(0.0, (clw(i,k,ntiw) +
2157 & snow2ice*clw(i,k,ntsw) + tem2) *
2158 & gfac * delp(i,k))
2159 if(tem2 > 1.e-12 .and. clw(i,k,ntiw) < 1.e-12) then
2160 if(nint(slmsk(i))==1) then
2161 rei(i,k)=creice_def
2162 else
2163 rei(i,k)=reice_def
2164 endif
2165 endif
2166 crp(i,k) = max(0.0, clw(i,k,ntrw) * gfac * delp(i,k))
2167 csp(i,k) = max(0.0, (1.-snow2ice)*clw(i,k,ntsw) *
2168 & gfac * delp(i,k))
2169 enddo
2170 enddo
2171
2173
2174 do i = 1, ix
2175 lwp_ex(i) = 0.0
2176 iwp_ex(i) = 0.0
2177 lwp_fc(i) = 0.0
2178 iwp_fc(i) = 0.0
2179 do k = 1, nlay-1
2180 lwp_ex(i) = lwp_ex(i) + cwp(i,k)
2181 iwp_ex(i) = iwp_ex(i) + cip(i,k) + csp(i,k)
2182 enddo
2183 lwp_ex(i) = lwp_ex(i)*1.e-3
2184 iwp_ex(i) = iwp_ex(i)*1.e-3
2185 enddo
2186
2187 if (uni_cld) then ! use unified sgs clouds generated outside
2188 do k = 1, nlay
2189 do i = 1, ix
2190 cldtot(i,k) = cldcov(i,k)
2191 enddo
2192 enddo
2193
2194 else
2195
2197
2198 if (.not. lmfshal) then
2200 & ( ix, nlay, xr_con, xr_exp, plyr, clwf, rhly, qstl, & ! --- inputs
2201 & cldtot ) & ! --- outputs
2202 else
2204 & ( ix, nlay, lmfdeep2, xr_con, xr_exp, plyr, clwf, rhly, &
2205 & qstl, & ! --- inputs
2206 & cldtot )
2207 endif
2208
2209 endif ! if (uni_cld) then
2210
2211 do k = 1, nlay
2212 do i = 1, ix
2213 if (cldtot(i,k) < climit) then
2214 cldtot(i,k) = 0.0
2215 cwp(i,k) = 0.0
2216 cip(i,k) = 0.0
2217 crp(i,k) = 0.0
2218 csp(i,k) = 0.0
2219 endif
2220 enddo
2221 enddo
2222
2223 ! What portion of water and ice contents is associated with the
2224 ! partly cloudy boxes
2225 do i = 1, ix
2226 do k = 1, nlay-1
2227 if (cldtot(i,k).ge.climit .and. cldtot(i,k).lt.ovcst) then
2228 lwp_fc(i) = lwp_fc(i) + cwp(i,k)
2229 iwp_fc(i) = iwp_fc(i) + cip(i,k) + csp(i,k)
2230 endif
2231 enddo
2232 lwp_fc(i) = lwp_fc(i)*1.e-3
2233 iwp_fc(i) = iwp_fc(i)*1.e-3
2234 enddo
2235
2236 if ( lcnorm ) then
2237 do k = 1, nlay
2238 do i = 1, ix
2239 if (cldtot(i,k) >= climit) then
2240 tem1 = 1.0 / max(climit2, cldtot(i,k))
2241 cwp(i,k) = cwp(i,k) * tem1
2242 cip(i,k) = cip(i,k) * tem1
2243 crp(i,k) = crp(i,k) * tem1
2244 csp(i,k) = csp(i,k) * tem1
2245 endif
2246 enddo
2247 enddo
2248 endif
2249
2250 do k = 1, nlay
2251 do i = 1, ix
2252 cld_frac(i,k) = cldtot(i,k)
2253 cld_lwp(i,k) = cwp(i,k)
2254 cld_reliq(i,k) = rew(i,k)
2255 cld_iwp(i,k) = cip(i,k)
2256 cld_reice(i,k) = rei(i,k)
2257 cld_rwp(i,k) = crp(i,k) ! added for Thompson
2258 cld_rerain(i,k) = rer(i,k)
2259 cld_swp(i,k) = csp(i,k) ! added for Thompson
2260 cld_resnow(i,k) = res(i,k)
2261 enddo
2262 enddo
2263
2264!............................................
2265 end subroutine progcld_thompson_wsm6
2266!............................................
2267
2268
2278 subroutine progcld_thompson &
2279 & ( plyr,plvl,tlyr,qlyr,qstl,rhly,clw, & ! --- inputs:
2280 & xlat,xlon,slmsk,dz,delp, &
2281 & ntrac,ntcw,ntiw,ntrw,ntsw,ntgl, &
2282 & ix, nlay, nlp1, &
2283 & uni_cld, lmfshal, lmfdeep2, cldcov, &
2284 & re_cloud,re_ice,re_snow, &
2285 & lwp_ex, iwp_ex, lwp_fc, iwp_fc, &
2286 & dzlay, gridkm, top_at_1, cldtot, cldcnv, &
2287 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
2288 & cld_reice,cld_rwp, cld_rerain,cld_swp, cld_resnow &
2289 & )
2290
2291! ================= subprogram documentation block ================ !
2292! !
2293! subprogram: progcld_thompson computes cloud related quantities !
2294! using Thompson cloud microphysics scheme. !
2295! !
2296! abstract: this program computes cloud fractions from cloud !
2297! condensates, !
2298! and computes the low, mid, high, total and boundary layer cloud !
2299! fractions and the vertical indices of low, mid, and high cloud !
2300! top and base. the three vertical cloud domains are set up in the !
2301! initial subroutine "cld_init". !
2302! !
2303! usage: call progcld_thompson !
2304! !
2305! attributes: !
2306! language: fortran 90 !
2307! machine: ibm-sp, sgi !
2308! !
2309! !
2310! ==================== definition of variables ==================== !
2311! !
2312! !
2313! input variables: !
2314! plyr (IX,NLAY) : model layer mean pressure in mb (100Pa) !
2315! plvl (IX,NLP1) : model level pressure in mb (100Pa) !
2316! tlyr (IX,NLAY) : model layer mean temperature in k !
2317! tvly (IX,NLAY) : model layer virtual temperature in k !
2318! qlyr (IX,NLAY) : layer specific humidity in gm/gm !
2319! qstl (IX,NLAY) : layer saturate humidity in gm/gm !
2320! rhly (IX,NLAY) : layer relative humidity (=qlyr/qstl) !
2321! clw (IX,NLAY,ntrac) : layer cloud condensate amount !
2322! xlat (IX) : grid latitude in radians, default to pi/2 -> -pi/2!
2323! range, otherwise see in-line comment !
2324! xlon (IX) : grid longitude in radians (not used) !
2325! slmsk (IX) : sea/land mask array (sea:0,land:1,sea-ice:2) !
2326! dz (ix,nlay) : layer thickness (km) !
2327! delp (ix,nlay) : model layer pressure thickness in mb (100Pa) !
2328! gridkm (IX) : grid length in km !
2329! IX : horizontal dimention !
2330! NLAY,NLP1 : vertical layer/level dimensions !
2331! uni_cld : logical - true for cloud fraction from shoc !
2332! lmfshal : logical - true for mass flux shallow convection !
2333! lmfdeep2 : logical - true for mass flux deep convection !
2334! top_at_1 : logical - true if vertical ordereing is toa-2-sfc !
2335! cldcov : layer cloud fraction (used when uni_cld=.true. !
2336! lmfshal : mass-flux shallow conv scheme flag !
2337! lmfdeep2 : scale-aware mass-flux deep conv scheme flag !
2338! lcrick : control flag for eliminating CRICK !
2339! =t: apply layer smoothing to eliminate CRICK !
2340! =f: do not apply layer smoothing !
2341! lcnorm : control flag for in-cld condensate !
2342! =t: normalize cloud condensate !
2343! =f: not normalize cloud condensate !
2344! !
2345! output variables: !
2346! cloud profiles: !
2347! cld_frac (:,:) - layer total cloud fraction !
2348! cld_lwp (:,:) - layer cloud liq water path (g/m**2) !
2349! cld_reliq (:,:) - mean eff radius for liq cloud (micron) !
2350! cld_iwp (:,:) - layer cloud ice water path (g/m**2) !
2351! cld_reice (:,:) - mean eff radius for ice cloud (micron) !
2352! cld_rwp (:,:) - layer rain drop water path not assigned !
2353! cld_rerain(:,:) - mean eff radius for rain drop (micron) !
2354! *** cld_swp (:,:) - layer snow flake water path not assigned !
2355! cld_resnow(:,:) - mean eff radius for snow flake (micron) !
2356! !
2357! ==================== end of description ===================== !
2358!
2359 implicit none
2360
2361! --- inputs
2362 integer, intent(in) :: ix, nlay, nlp1
2363 integer, intent(in) :: ntrac, ntcw, ntiw, ntrw, ntsw, ntgl
2364
2365 logical, intent(in) :: uni_cld, lmfshal, lmfdeep2, top_at_1
2366
2367 real (kind=kind_phys), dimension(:,:), intent(in) :: plvl, plyr, &
2368 & tlyr, qlyr, qstl, rhly, cldcov, delp, dz, dzlay, &
2369 & re_cloud, re_ice, re_snow
2370 real (kind=kind_phys), dimension(:), intent(inout) :: &
2371 & lwp_ex, iwp_ex, lwp_fc, iwp_fc
2372
2373 real (kind=kind_phys), dimension(:,:,:), intent(in) :: clw
2374
2375 real (kind=kind_phys), dimension(:), intent(in) :: xlat, xlon, &
2376 & slmsk
2377 real(kind=kind_phys), dimension(:), intent(in) :: gridkm
2378
2379! --- inputs/outputs
2380
2381 real (kind=kind_phys), dimension(:,:), intent(inout) :: &
2382 & cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, &
2383 & cld_rwp, cld_rerain, cld_swp, cld_resnow
2384
2385! --- local variables:
2386 real (kind=kind_phys), dimension(IX,NLAY) :: cldtot, cldcnv, &
2387 & cwp, cip, crp, csp, rew, rei, res, rer
2388
2389 real (kind=kind_phys), dimension(NLAY) :: cldfra1d, qv1d, &
2390 & qc1d, qi1d, qs1d, dz1d, p1d, t1d
2391
2392 real (kind=kind_phys) :: clwmin, tem1
2393 real (kind=kind_phys) :: corr, xland, snow_mass_factor
2394 real (kind=kind_phys), parameter :: max_relh = 1.5
2395 real (kind=kind_phys), parameter :: snow_max_radius = 130.0
2396
2397 integer :: i, k, k2, id, nf, idx_rei
2398!
2399!===> ... begin here
2400!
2401
2402 clwmin = 1.0e-9
2403
2404 do k = 1, nlay
2405 do i = 1, ix
2406 cldtot(i,k) = 0.0
2407 cldcnv(i,k) = 0.0
2408 cwp(i,k) = 0.0
2409 cip(i,k) = 0.0
2410 crp(i,k) = 0.0
2411 csp(i,k) = 0.0
2412 rew(i,k) = re_cloud(i,k)
2413 rei(i,k) = re_ice(i,k)
2414 rer(i,k) = rrain_def ! default rain radius to 1000 micron
2415 res(i,k) = re_snow(i,k)
2416 enddo
2417 enddo
2418
2422
2423 do k = 1, nlay-1
2424 do i = 1, ix
2425 cwp(i,k) = max(0.0, clw(i,k,ntcw) * dz(i,k)*1.e6)
2426 crp(i,k) = 0.0
2427 snow_mass_factor = 0.90
2428 cip(i,k) = max(0.0, (clw(i,k,ntiw) &
2429 & + (1.0-snow_mass_factor)*clw(i,k,ntsw))*dz(i,k)*1.e6)
2430 if (re_snow(i,k) .gt. snow_max_radius)then
2431 snow_mass_factor = min(snow_mass_factor, &
2432 & (snow_max_radius/re_snow(i,k)) &
2433 & *(snow_max_radius/re_snow(i,k)))
2434 res(i,k) = snow_max_radius
2435 endif
2436 csp(i,k) = max(0.,snow_mass_factor*clw(i,k,ntsw)*dz(i,k)*1.e6)
2437 enddo
2438 enddo
2439
2441
2442 do i = 1, ix
2443 lwp_ex(i) = 0.0
2444 iwp_ex(i) = 0.0
2445 do k = 1, nlay-1
2446 lwp_ex(i) = lwp_ex(i) + cwp(i,k)
2447 iwp_ex(i) = iwp_ex(i) + cip(i,k) + csp(i,k)
2448 enddo
2449 enddo
2450
2459
2460 do i = 1, ix
2461 if (slmsk(i)-0.5 .gt. 0.5 .and. slmsk(i)+0.5 .lt. 1.5) then
2462 xland = 1.0
2463 else
2464 xland = 2.0
2465 endif
2466
2467 cldfra1d(:) = 0.0
2468
2469 if (.not. top_at_1) then
2470 do k = 1, nlay
2471 qv1d(k) = qlyr(i,k)
2472 qc1d(k) = max(0.0, clw(i,k,ntcw))
2473 qi1d(k) = max(0.0, clw(i,k,ntiw))
2474 qs1d(k) = max(0.0, clw(i,k,ntsw))
2475 dz1d(k) = dz(i,k)*1.e3
2476 p1d(k) = plyr(i,k)*100.0
2477 t1d(k) = tlyr(i,k)
2478 enddo
2479 else
2480 do k = nlay, 1, -1
2481 k2 = nlay - k + 1
2482 qv1d(k2) = qlyr(i,k)
2483 qc1d(k2) = max(0.0, clw(i,k,ntcw))
2484 qi1d(k2) = max(0.0, clw(i,k,ntiw))
2485 qs1d(k2) = max(0.0, clw(i,k,ntsw))
2486 dz1d(k2) = dz(i,k)*1.e3
2487 p1d(k2) = plyr(i,k)*100.0
2488 t1d(k2) = tlyr(i,k)
2489 enddo
2490 endif
2491
2492 call cal_cldfra3(cldfra1d, qv1d, qc1d, qi1d, qs1d, dz1d, &
2493 & p1d, t1d, xland, gridkm(i), &
2494 & .false., max_relh, 1, nlay, .false.)
2495
2496 do k = 1, nlay
2497 cldtot(i,k) = cldfra1d(k)
2498 if (qc1d(k).gt.clwmin .and. cldfra1d(k).lt.ovcst) then
2499 cwp(i,k) = qc1d(k) * dz1d(k)*1000.
2500 if ((xland-1.5).GT.0.) then !--- Ocean
2501 rew(i,k) = 9.5
2502 else !--- Land
2503 rew(i,k) = 5.5
2504 endif
2505 endif
2506 if (qi1d(k).gt.clwmin .and. cldfra1d(k).lt.ovcst) then
2507 cip(i,k) = qi1d(k) * dz1d(k)*1000.
2508 idx_rei = int(t1d(k)-179.)
2509 idx_rei = min(max(idx_rei,1),75)
2510 corr = t1d(k) - int(t1d(k))
2511 rei(i,k) = max(5.0, retab(idx_rei)*(1.-corr) + &
2512 & retab(idx_rei+1)*corr)
2513 endif
2514 enddo
2515 enddo
2516
2517 do k = 1, nlay
2518 do i = 1, ix
2519 cld_frac(i,k) = cldtot(i,k)
2520 cld_lwp(i,k) = cwp(i,k)
2521 cld_reliq(i,k) = rew(i,k)
2522 cld_iwp(i,k) = cip(i,k)
2523 cld_reice(i,k) = rei(i,k)
2524 cld_rwp(i,k) = crp(i,k) ! added for Thompson
2525 cld_rerain(i,k) = rer(i,k)
2526 cld_swp(i,k) = csp(i,k) ! added for Thompson
2527 cld_resnow(i,k) = res(i,k)
2528 enddo
2529 enddo
2530
2532
2533 do i = 1, ix
2534 lwp_fc(i) = 0.0
2535 iwp_fc(i) = 0.0
2536 do k = 1, nlay
2537 lwp_fc(i) = lwp_fc(i) + cwp(i,k)
2538 iwp_fc(i) = iwp_fc(i) + cip(i,k) + csp(i,k)
2539 enddo
2540 lwp_fc(i) = max(0.0, lwp_fc(i) - lwp_ex(i))
2541 iwp_fc(i) = max(0.0, iwp_fc(i) - iwp_ex(i))
2542 lwp_fc(i) = lwp_fc(i)*1.e-3
2543 iwp_fc(i) = iwp_fc(i)*1.e-3
2544 lwp_ex(i) = lwp_ex(i)*1.e-3
2545 iwp_ex(i) = iwp_ex(i)*1.e-3
2546 enddo
2547!
2548!............................................
2549 end subroutine progcld_thompson
2550!............................................
2551!mz
2552
2553
2557 subroutine progclduni &
2558 & ( plyr,plvl,tlyr,tvly,ccnd,ncnd, & ! --- inputs:
2559 & xlat,xlon,slmsk,dz,delp, ix, nlay, nlp1, cldtot, &
2560 & effrl,effri,effrr,effrs,effr_in, &
2561 & dzlay, cldtot1, cldcnv, lcrick, lcnorm, con_ttp, &
2562 & cld_frac, cld_lwp, cld_reliq, cld_iwp, & ! --- outputs
2563 & cld_reice,cld_rwp, cld_rerain,cld_swp, cld_resnow &
2564 & )
2565
2566! ================= subprogram documentation block ================ !
2567! !
2568! subprogram: progclduni computes cloud related quantities using !
2569! for unified cloud microphysics scheme. !
2570! !
2571! abstract: this program computes cloud fractions from cloud !
2572! condensates, calculates liquid/ice cloud droplet effective radius, !
2573! and computes the low, mid, high, total and boundary layer cloud !
2574! fractions and the vertical indices of low, mid, and high cloud !
2575! top and base. the three vertical cloud domains are set up in the !
2576! initial subroutine "cld_init". !
2577! This program is written by Moorthi !
2578! to represent unified cloud across all physics while !
2579! using SHOC+MG2/3+convection (RAS or SAS or CSAW) !
2580! !
2581! usage: call progclduni !
2582! !
2583! attributes: !
2584! language: fortran 90 !
2585! machine: ibm-sp, sgi !
2586! !
2587! !
2588! ==================== definition of variables ==================== !
2589! !
2590! input variables: !
2591! plyr (IX,NLAY) : model layer mean pressure in mb (100Pa) !
2592! plvl (IX,NLP1) : model level pressure in mb (100Pa) !
2593! tlyr (IX,NLAY) : model layer mean temperature in k !
2594! tvly (IX,NLAY) : model layer virtual temperature in k !
2595! ccnd (IX,NLAY,ncnd) : layer cloud condensate amount !
2596! water, ice, rain, snow (+ graupel) !
2597! ncnd : number of layer cloud condensate types (max of 4) !
2598! xlat (IX) : grid latitude in radians, default to pi/2 -> -pi/2!
2599! range, otherwise see in-line comment !
2600! xlon (IX) : grid longitude in radians (not used) !
2601! slmsk (IX) : sea/land mask array (sea:0,land:1,sea-ice:2) !
2602! IX : horizontal dimention !
2603! NLAY,NLP1 : vertical layer/level dimensions !
2604! cldtot : unified cloud fracrion from moist physics !
2605! effrl (ix,nlay) : effective radius for liquid water !
2606! effri (ix,nlay) : effective radius for ice water !
2607! effrr (ix,nlay) : effective radius for rain water !
2608! effrs (ix,nlay) : effective radius for snow water !
2609! effr_in : logical - if .true. use input effective radii !
2610! dz (ix,nlay) : layer thickness (km) !
2611! delp (ix,nlay) : model layer pressure thickness in mb (100Pa) !
2612! dzlay(ix,nlay) : thickness between model layer centers (km) !
2613! lmfshal : mass-flux shallow conv scheme flag !
2614! lmfdeep2 : scale-aware mass-flux deep conv scheme flag !
2615! lcrick : control flag for eliminating CRICK !
2616! =t: apply layer smoothing to eliminate CRICK !
2617! =f: do not apply layer smoothing !
2618! lcnorm : control flag for in-cld condensate !
2619! =t: normalize cloud condensate !
2620! =f: not normalize cloud condensate !
2621! !
2622! output variables: !
2623! cloud profiles: !
2624! cld_frac (:,:) - layer total cloud fraction !
2625! cld_lwp (:,:) - layer cloud liq water path (g/m**2) !
2626! cld_reliq (:,:) - mean eff radius for liq cloud (micron) !
2627! cld_iwp (:,:) - layer cloud ice water path (g/m**2) !
2628! cld_reice (:,:) - mean eff radius for ice cloud (micron) !
2629! cld_rwp (:,:) - layer rain drop water path not assigned !
2630! cld_rerain(:,:) - mean eff radius for rain drop (micron) !
2631! *** cld_swp (:,:) - layer snow flake water path not assigned !
2632! cld_resnow(:,:) - mean eff radius for snow flake (micron) !
2633! *** fu's scheme need to be normalized by snow density (g/m**3/1.0e6) !
2634! clds (IX,5) : fraction of clouds for low, mid, hi, tot, bl !
2635! mtop (IX,3) : vertical indices for low, mid, hi cloud tops !
2636! mbot (IX,3) : vertical indices for low, mid, hi cloud bases !
2637! de_lgth(ix) : clouds decorrelation length (km) !
2638! alpha(ix,nlay) : alpha decorrelation parameter !
2639! !
2640! ==================== end of description ===================== !
2641!
2642 implicit none
2643
2644! --- inputs
2645 integer, intent(in) :: ix, nlay, nlp1, ncnd
2646 logical, intent(in) :: effr_in, lcrick, lcnorm
2647
2648 real (kind=kind_phys), intent(in) :: con_ttp
2649 real (kind=kind_phys), dimension(:,:,:), intent(in) :: ccnd
2650 real (kind=kind_phys), dimension(:,:), intent(in) :: plvl, plyr,&
2651 & tlyr, tvly, cldtot, effrl, effri, effrr, effrs, dz, delp, &
2652 & dzlay
2653
2654 real (kind=kind_phys), dimension(:), intent(in) :: xlat, xlon, &
2655 & slmsk
2656
2657 real (kind=kind_phys), dimension(:,:), intent(inout) :: cldtot1
2658
2659! --- inputs/outputs
2660
2661 real (kind=kind_phys), dimension(:,:), intent(inout) :: &
2662 & cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, &
2663 & cld_rwp, cld_rerain, cld_swp, cld_resnow
2664
2665! --- local variables:
2666 real (kind=kind_phys), dimension(IX,NLAY) :: cldcnv, cwp, cip, &
2667 & crp, csp, rew, rei, res, rer
2668 real (kind=kind_phys), dimension(IX,NLAY,ncnd) :: cndf
2669
2670 real (kind=kind_phys) :: tem1, tem2, tem3
2671
2672 integer :: i, k, id, nf, n
2673
2674!
2675!===> ... begin here
2676!
2677 do k = 1, nlay
2678 do i = 1, ix
2679 cldcnv(i,k) = 0.0
2680 cwp(i,k) = 0.0
2681 cip(i,k) = 0.0
2682 crp(i,k) = 0.0
2683 csp(i,k) = 0.0
2684 enddo
2685 enddo
2686 do n=1,ncnd
2687 do k = 1, nlay
2688 do i = 1, ix
2689 cndf(i,k,n) = ccnd(i,k,n)
2690 enddo
2691 enddo
2692 enddo
2693 if ( lcrick ) then ! vertical smoorthing
2694 do n=1,ncnd
2695 do i = 1, ix
2696 cndf(i,1,n) = 0.75*ccnd(i,1,n) + 0.25*ccnd(i,2,n)
2697 cndf(i,nlay,n) = 0.75*ccnd(i,nlay,n) + 0.25*ccnd(i,nlay-1,n)
2698 enddo
2699 do k = 2, nlay-1
2700 do i = 1, ix
2701 cndf(i,k,n) = 0.25 * (ccnd(i,k-1,n) + ccnd(i,k+1,n)) &
2702 & + 0.5 * ccnd(i,k,n)
2703 enddo
2704 enddo
2705 enddo
2706 endif
2707
2709
2710 if (ncnd == 2) then
2711 do k = 1, nlay
2712 do i = 1, ix
2713 tem1 = gfac * delp(i,k)
2714 cwp(i,k) = cndf(i,k,1) * tem1
2715 cip(i,k) = cndf(i,k,2) * tem1
2716 enddo
2717 enddo
2718 elseif (ncnd == 4) then
2719 do k = 1, nlay
2720 do i = 1, ix
2721 tem1 = gfac * delp(i,k)
2722 cwp(i,k) = cndf(i,k,1) * tem1
2723 cip(i,k) = cndf(i,k,2) * tem1
2724 crp(i,k) = cndf(i,k,3) * tem1
2725 csp(i,k) = cndf(i,k,4) * tem1
2726 enddo
2727 enddo
2728 endif
2729
2730 do k = 1, nlay
2731 do i = 1, ix
2732 if (cldtot(i,k) < climit) then
2733 cwp(i,k) = 0.0
2734 cip(i,k) = 0.0
2735 crp(i,k) = 0.0
2736 csp(i,k) = 0.0
2737 endif
2738 enddo
2739 enddo
2740
2741 if ( lcnorm ) then
2742 do k = 1, nlay
2743 do i = 1, ix
2744 if (cldtot(i,k) >= climit) then
2745 tem1 = 1.0 / max(climit2, cldtot(i,k))
2746 cwp(i,k) = cwp(i,k) * tem1
2747 cip(i,k) = cip(i,k) * tem1
2748 crp(i,k) = crp(i,k) * tem1
2749 csp(i,k) = csp(i,k) * tem1
2750 endif
2751 enddo
2752 enddo
2753 endif
2754
2755! assign/calculate efective radii for cloud water, ice, rain, snow
2756
2757 if (effr_in) then
2758 do k = 1, nlay
2759 do i = 1, ix
2760 rew(i,k) = effrl(i,k)
2761 rei(i,k) = max(10.0, min(150.0,effri(i,k)))
2762 rer(i,k) = effrr(i,k)
2763 res(i,k) = effrs(i,k)
2764 enddo
2765 enddo
2766 else
2767 do k = 1, nlay
2768 do i = 1, ix
2769 rew(i,k) = reliq_def ! default liq radius to 10 micron
2770 rei(i,k) = reice_def ! default ice radius to 50 micron
2771 rer(i,k) = rrain_def ! default rain radius to 1000 micron
2772 res(i,k) = rsnow_def ! default snow radius to 250 micron
2773 enddo
2774 enddo
2776 do i = 1, ix
2777 if (nint(slmsk(i)) == 1) then
2778 do k = 1, nlay
2779 tem1 = min(1.0, max(0.0, (con_ttp-tlyr(i,k))*0.05))
2780 rew(i,k) = 5.0 + 5.0 * tem1
2781 enddo
2782 endif
2783 enddo
2784
2787
2788 do k = 1, nlay
2789 do i = 1, ix
2790 tem2 = tlyr(i,k) - con_ttp
2791
2792 if (cip(i,k) > 0.0) then
2793 tem3 = gord * cip(i,k) * plyr(i,k) / (delp(i,k)*tvly(i,k))
2794
2795 if (tem2 < -50.0) then
2796 rei(i,k) = (1250.0/9.917) * tem3 ** 0.109
2797 elseif (tem2 < -40.0) then
2798 rei(i,k) = (1250.0/9.337) * tem3 ** 0.08
2799 elseif (tem2 < -30.0) then
2800 rei(i,k) = (1250.0/9.208) * tem3 ** 0.055
2801 else
2802 rei(i,k) = (1250.0/9.387) * tem3 ** 0.031
2803 endif
2804! rei(i,k) = max(20.0, min(rei(i,k), 300.0))
2805! rei(i,k) = max(10.0, min(rei(i,k), 100.0))
2806 rei(i,k) = max(10.0, min(rei(i,k), 150.0))
2807! rei(i,k) = max(5.0, min(rei(i,k), 130.0))
2808 endif
2809 enddo
2810 enddo
2811 endif
2812
2813 do k = 1, nlay
2814 do i = 1, ix
2815 cldtot1(i,k) = cldtot(i,k)
2816 enddo
2817 enddo
2818!
2819 do k = 1, nlay
2820 do i = 1, ix
2821 cld_frac(i,k) = cldtot(i,k)
2822 cld_lwp(i,k) = cwp(i,k)
2823 cld_reliq(i,k) = rew(i,k)
2824 cld_iwp(i,k) = cip(i,k)
2825 cld_reice(i,k) = rei(i,k)
2826 cld_rwp(i,k) = crp(i,k) ! added for Thompson
2827 cld_rerain(i,k) = rer(i,k)
2828 cld_swp(i,k) = csp(i,k) ! added for Thompson
2829 cld_resnow(i,k) = res(i,k)
2830 enddo
2831 enddo
2832!
2833!...................................
2834 end subroutine progclduni
2835!-----------------------------------
2836
2863 subroutine gethml &
2864 & ( plyr, ptop1, cldtot, cldcnv, dz, de_lgth, alpha, & ! --- inputs:
2865 & ix, nlay, iovr, iovr_rand, iovr_maxrand, iovr_max, &
2866 & iovr_dcorr, iovr_exp, iovr_exprand, top_at_1, si, &
2867 & clds, mtop, mbot & ! --- outputs:
2868 & )
2869
2870! =================================================================== !
2871! !
2872! abstract: compute high, mid, low, total, and boundary cloud fractions !
2873! and cloud top/bottom layer indices for model diagnostic output. !
2874! the three cloud domain boundaries are defined by ptopc. the cloud !
2875! overlapping method is defined by control flag 'iovr', which is also !
2876! used by lw and sw radiation programs. !
2877! !
2878! usage: call gethml !
2879! !
2880! subprograms called: none !
2881! !
2882! attributes: !
2883! language: fortran 90 !
2884! machine: ibm-sp, sgi !
2885! !
2886! !
2887! ==================== definition of variables ==================== !
2888! !
2889! input variables: !
2890! plyr (IX,NLAY) : model layer mean pressure in mb (100Pa) !
2891! ptop1 (IX,4) : pressure limits of cloud domain interfaces !
2892! (sfc,low,mid,high) in mb (100Pa) !
2893! cldtot(IX,NLAY) : total or straiform cloud profile in fraction !
2894! cldcnv(IX,NLAY) : convective cloud (for diagnostic scheme only) !
2895! dz (ix,nlay) : layer thickness (km) !
2896! de_lgth(ix) : clouds vertical de-correlation length (km) !
2897! alpha(ix,nlay) : alpha decorrelation parameter
2898! IX : horizontal dimention !
2899! NLAY : vertical layer dimensions !
2900! !
2901! output variables: !
2902! clds (IX,5) : fraction of clouds for low, mid, hi, tot, bl !
2903! mtop (IX,3) : vertical indices for low, mid, hi cloud tops !
2904! mbot (IX,3) : vertical indices for low, mid, hi cloud bases ! !
2905! internal module variables: !
2906! iovr : control flag for cloud overlap !
2907! =0 random overlapping clouds !
2908! =1 max/ran overlapping clouds !
2909! =2 maximum overlapping ( for mcica only ) !
2910! =3 decorr-length ovlp ( for mcica only ) !
2911! =4: exponential cloud overlap (AER; mcica only) !
2912! =5: exponential-random overlap (AER; mcica only) !
2913! !
2914! ==================== end of description ===================== !
2915!
2916 implicit none!
2917
2918! --- inputs:
2919 logical, intent(in) :: top_at_1
2920 integer, intent(in) :: ix, nlay
2921 integer, intent(in) :: &
2922 & iovr, !
2923 & iovr_rand, ! Flag for random cloud overlap method
2924 & iovr_maxrand, ! Flag for maximum-random cloud overlap method
2925 & iovr_max, ! Flag for maximum cloud overlap method
2926 & iovr_dcorr, ! Flag for decorrelation-length cloud overlap method
2927 & iovr_exp, ! Flag for exponential cloud overlap method
2928 & iovr_exprand ! Flag for exponential-random cloud overlap method
2929
2930 real (kind=kind_phys), dimension(:,:), intent(in) :: plyr, ptop1, &
2931 & cldtot, cldcnv, dz
2932 real (kind=kind_phys), dimension(:), intent(in) :: de_lgth, si
2933 real (kind=kind_phys), dimension(:,:), intent(in) :: alpha
2934
2935! --- outputs
2936 real (kind=kind_phys), dimension(:,:), intent(out) :: clds
2937
2938 integer, dimension(:,:), intent(out) :: mtop, mbot
2939
2940! --- local variables:
2941 real (kind=kind_phys) :: cl1(ix), cl2(ix), dz1(ix)
2942 real (kind=kind_phys) :: pcur, pnxt, ccur, cnxt, alfa
2943
2944 integer, dimension(IX):: idom, kbt1, kth1, kbt2, kth2
2945 integer :: i, k, id, id1, kstr, kend, kinc,kl
2946
2947!
2948!===> ... begin here
2949!
2953
2954 if (top_at_1) then ! data from toa to sfc
2955 lab_do_k0 : do k = nlay, 2, -1
2956 kl = k
2957 if (si(k) < 0.9e0) exit lab_do_k0
2958 enddo lab_do_k0
2959 llyr = kl
2960 else ! data from sfc to top
2961 lab_do_k1 : do k = 2, nlay
2962 kl = k
2963 if (si(k) < 0.9e0) exit lab_do_k1
2964 enddo lab_do_k1
2965
2966 llyr = kl - 1
2967 endif ! end_if_top_at_1
2968
2969 clds(:,:) = 0.0
2970
2971 do i = 1, ix
2972 cl1(i) = 1.0
2973 cl2(i) = 1.0
2974 enddo
2975
2976! --- total and bl clouds, where cl1, cl2 are fractions of clear-sky view
2977! layer processed from surface and up
2978
2981
2982 if (top_at_1) then ! input data from toa to sfc
2983 kstr = nlay
2984 kend = 1
2985 kinc = -1
2986 else ! input data from sfc to toa
2987 kstr = 1
2988 kend = nlay
2989 kinc = 1
2990 endif ! end_if_top_at_1
2991
2992 if ( iovr == iovr_rand ) then ! random overlap
2993
2994 do k = kstr, kend, kinc
2995 do i = 1, ix
2996 ccur = min( ovcst, max( cldtot(i,k), cldcnv(i,k) ))
2997 if (ccur >= climit) cl1(i) = cl1(i) * (1.0 - ccur)
2998 enddo
2999
3000 if (k == llyr) then
3001 do i = 1, ix
3002 clds(i,5) = 1.0 - cl1(i) ! save bl cloud
3003 enddo
3004 endif
3005 enddo
3006
3007 do i = 1, ix
3008 clds(i,4) = 1.0 - cl1(i) ! save total cloud
3009 enddo
3010
3011 elseif ( iovr == iovr_maxrand ) then ! max/ran overlap
3012
3013 do k = kstr, kend, kinc
3014 do i = 1, ix
3015 ccur = min( ovcst, max( cldtot(i,k), cldcnv(i,k) ))
3016 if (ccur >= climit) then ! cloudy layer
3017 cl2(i) = min( cl2(i), (1.0 - ccur) )
3018 else ! clear layer
3019 cl1(i) = cl1(i) * cl2(i)
3020 cl2(i) = 1.0
3021 endif
3022 enddo
3023
3024 if (k == llyr) then
3025 do i = 1, ix
3026 clds(i,5) = 1.0 - cl1(i) * cl2(i) ! save bl cloud
3027 enddo
3028 endif
3029 enddo
3030
3031 do i = 1, ix
3032 clds(i,4) = 1.0 - cl1(i) * cl2(i) ! save total cloud
3033 enddo
3034
3035 elseif ( iovr == iovr_max ) then ! maximum overlap all levels
3036
3037 cl1(:) = 0.0
3038
3039 do k = kstr, kend, kinc
3040 do i = 1, ix
3041 ccur = min( ovcst, max( cldtot(i,k), cldcnv(i,k) ))
3042 if (ccur >= climit) cl1(i) = max( cl1(i), ccur )
3043 enddo
3044
3045 if (k == llyr) then
3046 do i = 1, ix
3047 clds(i,5) = cl1(i) ! save bl cloud
3048 enddo
3049 endif
3050 enddo
3051
3052 do i = 1, ix
3053 clds(i,4) = cl1(i) ! save total cloud
3054 enddo
3055
3056 elseif ( iovr == iovr_dcorr ) then ! random if clear-layer divided,
3057 ! otherwise de-corrlength method
3058 do i = 1, ix
3059 dz1(i) = - dz(i,kstr)
3060 enddo
3061
3062 do k = kstr, kend, kinc
3063 do i = 1, ix
3064 ccur = min( ovcst, max( cldtot(i,k), cldcnv(i,k) ))
3065 if (ccur >= climit) then ! cloudy layer
3066 alfa = exp( -0.5*((dz1(i)+dz(i,k)))/de_lgth(i) )
3067 dz1(i) = dz(i,k)
3068 cl2(i) = alfa * min(cl2(i), (1.0 - ccur)) & ! maximum part
3069 & + (1.0 - alfa) * (cl2(i) * (1.0 - ccur)) ! random part
3070 else ! clear layer
3071 cl1(i) = cl1(i) * cl2(i)
3072 cl2(i) = 1.0
3073 if (k /= kend) dz1(i) = -dz(i,k+kinc)
3074 endif
3075 enddo
3076
3077 if (k == llyr) then
3078 do i = 1, ix
3079 clds(i,5) = 1.0 - cl1(i) * cl2(i) ! save bl cloud
3080 enddo
3081 endif
3082 enddo
3083
3084 do i = 1, ix
3085 clds(i,4) = 1.0 - cl1(i) * cl2(i) ! save total cloud
3086 enddo
3087
3088 elseif ( iovr == iovr_exp .or. iovr == iovr_exprand ) then ! exponential overlap (iovr=4), or
3089 ! exponential-random (iovr=5);
3090 ! distinction defined by alpha
3091
3092 do k = kstr, kend, kinc
3093 do i = 1, ix
3094 ccur = min( ovcst, max( cldtot(i,k), cldcnv(i,k) ))
3095 if (ccur >= climit) then ! cloudy layer
3096 cl2(i) = alpha(i,k) * min(cl2(i), (1.0 - ccur)) & ! maximum part
3097 & + (1.0 - alpha(i,k)) * (cl2(i) * (1.0 - ccur)) ! random part
3098 else ! clear layer
3099 cl1(i) = cl1(i) * cl2(i)
3100 cl2(i) = 1.0
3101 endif
3102 enddo
3103
3104 if (k == llyr) then
3105 do i = 1, ix
3106 clds(i,5) = 1.0 - cl1(i) * cl2(i) ! save bl cloud
3107 enddo
3108 endif
3109 enddo
3110
3111 do i = 1, ix
3112 clds(i,4) = 1.0 - cl1(i) * cl2(i) ! save total cloud
3113 enddo
3114
3115 endif ! end_if_iovr
3116
3117! --- high, mid, low clouds, where cl1, cl2 are cloud fractions
3118! layer processed from one layer below llyr and up
3119! --- change! layer processed from surface to top, so low clouds will
3120! contains both bl and low clouds.
3121
3124 if (top_at_1) then ! input data from toa to sfc
3125
3126 do i = 1, ix
3127 cl1(i) = 0.0
3128 cl2(i) = 0.0
3129 kbt1(i) = nlay
3130 kbt2(i) = nlay
3131 kth1(i) = 0
3132 kth2(i) = 0
3133 idom(i) = 1
3134 mbot(i,1) = nlay
3135 mtop(i,1) = nlay
3136 mbot(i,2) = nlay - 1
3137 mtop(i,2) = nlay - 1
3138 mbot(i,3) = nlay - 1
3139 mtop(i,3) = nlay - 1
3140 enddo
3141
3142!org do k = llyr-1, 1, -1
3143 do k = nlay, 1, -1
3144 do i = 1, ix
3145 id = idom(i)
3146 id1= id + 1
3147
3148 pcur = plyr(i,k)
3149 ccur = min( ovcst, max( cldtot(i,k), cldcnv(i,k) ))
3150
3151 if (k > 1) then
3152 pnxt = plyr(i,k-1)
3153 cnxt = min( ovcst, max( cldtot(i,k-1), cldcnv(i,k-1) ))
3154 else
3155 pnxt = -1.0
3156 cnxt = 0.0
3157 endif
3158
3159 if (pcur < ptop1(i,id1)) then
3160 id = id + 1
3161 id1= id1 + 1
3162 idom(i) = id
3163 endif
3164
3165 if (ccur >= climit) then
3166 if (kth2(i) == 0) kbt2(i) = k
3167 kth2(i) = kth2(i) + 1
3168
3169 if ( iovr == iovr_rand ) then
3170 cl2(i) = cl2(i) + ccur - cl2(i)*ccur
3171 else
3172 cl2(i) = max( cl2(i), ccur )
3173 endif
3174
3175 if (cnxt < climit .or. pnxt < ptop1(i,id1)) then
3176 kbt1(i) = nint( (cl1(i)*kbt1(i) + cl2(i)*kbt2(i) ) &
3177 & / (cl1(i) + cl2(i)) )
3178 kth1(i) = nint( (cl1(i)*kth1(i) + cl2(i)*kth2(i) ) &
3179 & / (cl1(i) + cl2(i)) )
3180 cl1(i) = cl1(i) + cl2(i) - cl1(i)*cl2(i)
3181
3182 kbt2(i) = k - 1
3183 kth2(i) = 0
3184 cl2(i) = 0.0
3185 endif ! end_if_cnxt_or_pnxt
3186 endif ! end_if_ccur
3187
3188 if (pnxt < ptop1(i,id1)) then
3189 clds(i,id) = cl1(i)
3190 mtop(i,id) = min( kbt1(i), kbt1(i)-kth1(i)+1 )
3191 mbot(i,id) = kbt1(i)
3192
3193 cl1(i) = 0.0
3194 kbt1(i) = k - 1
3195 kth1(i) = 0
3196
3197 if (id1 <= nk_clds) then
3198 mbot(i,id1) = kbt1(i)
3199 mtop(i,id1) = kbt1(i)
3200 endif
3201 endif ! end_if_pnxt
3202
3203 enddo ! end_do_i_loop
3204 enddo ! end_do_k_loop
3205
3206 else ! input data from sfc to toa
3207
3208 do i = 1, ix
3209 cl1(i) = 0.0
3210 cl2(i) = 0.0
3211 kbt1(i) = 1
3212 kbt2(i) = 1
3213 kth1(i) = 0
3214 kth2(i) = 0
3215 idom(i) = 1
3216 mbot(i,1) = 1
3217 mtop(i,1) = 1
3218 mbot(i,2) = 2
3219 mtop(i,2) = 2
3220 mbot(i,3) = 2
3221 mtop(i,3) = 2
3222 enddo
3223
3224!org do k = llyr+1, NLAY
3225 do k = 1, nlay
3226 do i = 1, ix
3227 id = idom(i)
3228 id1= id + 1
3229
3230 pcur = plyr(i,k)
3231 ccur = min( ovcst, max( cldtot(i,k), cldcnv(i,k) ))
3232
3233 if (k < nlay) then
3234 pnxt = plyr(i,k+1)
3235 cnxt = min( ovcst, max( cldtot(i,k+1), cldcnv(i,k+1) ))
3236 else
3237 pnxt = -1.0
3238 cnxt = 0.0
3239 endif
3240
3241 if (pcur < ptop1(i,id1)) then
3242 id = id + 1
3243 id1= id1 + 1
3244 idom(i) = id
3245 endif
3246
3247 if (ccur >= climit) then
3248 if (kth2(i) == 0) kbt2(i) = k
3249 kth2(i) = kth2(i) + 1
3250
3251 if ( iovr == iovr_rand ) then
3252 cl2(i) = cl2(i) + ccur - cl2(i)*ccur
3253 else
3254 cl2(i) = max( cl2(i), ccur )
3255 endif
3256
3257 if (cnxt < climit .or. pnxt < ptop1(i,id1)) then
3258 kbt1(i) = nint( (cl1(i)*kbt1(i) + cl2(i)*kbt2(i)) &
3259 & / (cl1(i) + cl2(i)) )
3260 kth1(i) = nint( (cl1(i)*kth1(i) + cl2(i)*kth2(i)) &
3261 & / (cl1(i) + cl2(i)) )
3262 cl1(i) = cl1(i) + cl2(i) - cl1(i)*cl2(i)
3263
3264 kbt2(i) = k + 1
3265 kth2(i) = 0
3266 cl2(i) = 0.0
3267 endif ! end_if_cnxt_or_pnxt
3268 endif ! end_if_ccur
3269
3270 if (pnxt < ptop1(i,id1)) then
3271 clds(i,id) = cl1(i)
3272 mtop(i,id) = max( kbt1(i), kbt1(i)+kth1(i)-1 )
3273 mbot(i,id) = kbt1(i)
3274
3275 cl1(i) = 0.0
3276 kbt1(i) = min(k+1, nlay)
3277 kth1(i) = 0
3278
3279 if (id1 <= nk_clds) then
3280 mbot(i,id1) = kbt1(i)
3281 mtop(i,id1) = kbt1(i)
3282 endif
3283 endif ! end_if_pnxt
3284
3285 enddo ! end_do_i_loop
3286 enddo ! end_do_k_loop
3287
3288 endif ! end_if_top_at_1
3289
3290!
3291!...................................
3292 end subroutine gethml
3293!-----------------------------------
3294
3305 SUBROUTINE cal_cldfra3(CLDFRA, qv, qc, qi, qs, dz, &
3306 & p, t, XLAND, gridkm, &
3307 & modify_qvapor, max_relh, &
3308 & kts,kte, debug_flag)
3309!
3310 USE module_mp_thompson , ONLY : rsif, rslf
3311 IMPLICIT NONE
3312!
3313 INTEGER, INTENT(IN):: kts, kte
3314 LOGICAL, INTENT(IN):: modify_qvapor
3315 REAL, DIMENSION(kts:kte), INTENT(INOUT):: qv, qc, qi, cldfra
3316 REAL, DIMENSION(kts:kte), INTENT(IN):: p, t, dz, qs
3317 REAL, INTENT(IN):: gridkm, xland, max_relh
3318 LOGICAL, INTENT(IN):: debug_flag
3319
3320!..Local vars.
3321 REAL:: rh_00l, rh_00o, rh_00
3322 REAL:: entrmnt=0.5
3323 INTEGER:: k
3324 REAL:: tc, qvsi, qvsw, rhum, delz, var_temp
3325 REAL, DIMENSION(kts:kte):: qvs, rh, rhoa
3326 integer:: ndebug = 0
3327
3328 character*512 dbg_msg
3329
3330!+---+
3331
3332!..Initialize cloud fraction, compute RH, and rho-air.
3333
3334 DO k = kts,kte
3335 cldfra(k) = 0.0
3336
3337 qvsw = rslf(p(k), t(k))
3338 qvsi = rsif(p(k), t(k))
3339
3340 tc = t(k) - 273.15
3341 if (tc .ge. -12.0) then
3342 qvs(k) = qvsw
3343 elseif (tc .lt. -35.0) then
3344 qvs(k) = qvsi
3345 else
3346 qvs(k) = qvsw - (qvsw-qvsi)*(-12.0-tc)/(-12.0+35.)
3347 endif
3348
3349 if (modify_qvapor) then
3350 if (qc(k).gt.1.e-8) then
3351 qv(k) = max(qv(k), qvsw)
3352 qvs(k) = qvsw
3353 endif
3354 if (qc(k).le.1.e-8 .and. qi(k).ge.1.e-9) then
3355 qv(k) = max(qv(k), qvsi*1.005) !..To ensure a tiny bit ice supersaturated
3356 qvs(k) = qvsi
3357 endif
3358 endif
3359
3360 rh(k) = max(0.01, qv(k)/qvs(k))
3361 rhoa(k) = p(k)/(287.0*t(k))
3362
3363 ENDDO
3364
3365
3366!..First cut scale-aware. Higher resolution should require closer to
3367!.. saturated grid box for higher cloud fraction. Simple functions
3368!.. chosen based on Mocko and Cotton (1995) starting point and desire
3369!.. to get near 100% RH as grid spacing moves toward 1.0km, but higher
3370!.. RH over ocean required as compared to over land.
3371
3372 DO k = kts,kte
3373
3374 delz = max(100., dz(k))
3375 rh_00l = 0.77+min(0.22,sqrt(1./(50.0+gridkm*gridkm*delz*0.01)))
3376 rh_00o = 0.85+min(0.14,sqrt(1./(50.0+gridkm*gridkm*delz*0.01)))
3377 rhum = rh(k)
3378
3379 if (qc(k).ge.1.e-6 .or. qi(k).ge.1.e-7 &
3380 & .or. (qs(k).gt.1.e-6 .and. t(k).lt.273.)) then
3381 cldfra(k) = 1.0
3382 elseif (((qc(k)+qi(k)).gt.1.e-10) .and. &
3383 & ((qc(k)+qi(k)).lt.1.e-6)) then
3384 var_temp = min(0.99, 0.1*(11.0 + log10(qc(k)+qi(k))))
3385 cldfra(k) = var_temp*var_temp
3386 else
3387
3388 IF ((xland-1.5).GT.0.) THEN !--- Ocean
3389 rh_00 = rh_00o
3390 ELSE !--- Land
3391 rh_00 = rh_00l
3392 ENDIF
3393
3394 tc = max(-80.0, t(k) - 273.15)
3395 if (tc .lt. -24.0) rh_00 = rh_00l
3396
3397 if (tc .gt. 20.0) then
3398 cldfra(k) = 0.0
3399 elseif (tc .ge. -24.0) then
3400 rhum = min(rh(k), 1.0)
3401 var_temp = sqrt(sqrt((1.001-rhum)/(1.001-rh_00)))
3402 cldfra(k) = max(0., 1.0-var_temp)
3403 else
3404 if (max_relh.gt.1.12 .or. (.NOT.(modify_qvapor)) ) then
3405!..For HRRR model, the following look OK.
3406 rhum = min(rh(k), 1.45)
3407 rh_00 = rh_00 + (1.45-rh_00)*(-24.0-tc)/(-24.0+85.)
3408 var_temp = sqrt(sqrt((1.46-rhum)/(1.46-rh_00)))
3409 cldfra(k) = max(0., 1.0-var_temp)
3410 else
3411!..but for the GFS model, RH is way lower.
3412 rhum = min(rh(k), 1.05)
3413 rh_00 = rh_00 + (1.05-rh_00)*(-24.0-tc)/(-24.0+85.)
3414 var_temp = sqrt(sqrt((1.06-rhum)/(1.06-rh_00)))
3415 cldfra(k) = max(0., 1.0-var_temp)
3416 endif
3417 endif
3418 if (cldfra(k).gt.0.) cldfra(k)=max(0.01,min(cldfra(k),0.99))
3419
3420 endif
3421 if (cldfra(k).gt.0.0 .and. p(k).lt.7000.0) cldfra(k) = 0.0
3422 ENDDO
3423
3424 call find_cloudlayers(qvs, cldfra, t, p, dz, entrmnt, &
3425 & debug_flag, qc, qi, qs, kts,kte)
3426
3427!..Do a final total column adjustment since we may have added more than 1 mm
3428!.. LWP/IWP for multiple cloud decks.
3429
3430 call adjust_cloudfinal(cldfra, qc, qi, rhoa, dz, kts,kte)
3431
3432!..Intended for cold start model runs, we use modify_qvapor to ensure that cloudy
3433!.. areas are actually saturated such that the inserted clouds do not evaporate a
3434!.. timestep later.
3435
3436 if (modify_qvapor) then
3437 DO k = kts,kte
3438 if (cldfra(k).gt.0.2 .and. cldfra(k).lt.1.0) then
3439 qv(k) = max(qv(k),qvs(k))
3440 endif
3441 ENDDO
3442 endif
3443
3444
3445 END SUBROUTINE cal_cldfra3
3446
3450 SUBROUTINE find_cloudlayers(qvs1d, cfr1d, T1d, P1d, Dz1d, entrmnt,&
3451 & debugfl, qc1d, qi1d, qs1d, kts,kte)
3452!
3453 IMPLICIT NONE
3454!
3455 INTEGER, INTENT(IN):: kts, kte
3456 LOGICAL, INTENT(IN):: debugfl
3457 REAL, INTENT(IN):: entrmnt
3458 REAL, DIMENSION(kts:kte), INTENT(IN):: qs1d,qvs1d,t1d,p1d,dz1d
3459 REAL, DIMENSION(kts:kte), INTENT(INOUT):: cfr1d, qc1d, qi1d
3460
3461!..Local vars.
3462 REAL, DIMENSION(kts:kte):: theta
3463 REAL:: theta1, theta2, delz
3464 INTEGER:: k, k2, k_tropo, k_m12c, k_cldb, k_cldt, kbot, k_p200
3465 LOGICAL:: in_cloud
3466 character*512 dbg_msg
3467
3468!+---+
3469
3470 k_m12c = 0
3471 k_p200 = 0
3472 DO k = kte, kts, -1
3473 theta(k) = t1d(k)*((100000.0/p1d(k))**(287.05/1004.))
3474 if (t1d(k)-273.16 .gt. -12.0 .and. p1d(k).gt.10100.0) &
3475 & k_m12c = max(k_m12c, k)
3476 if (p1d(k).gt.19999.0 .and. k_p200.eq.0) k_p200 = k
3477 ENDDO
3478 if (k_m12c .le. kts) k_m12c = kts
3479
3480!..Find tropopause height, best surrogate, because we would not really
3481!.. wish to put fake clouds into the stratosphere. The 10/1500 ratio
3482!.. d(Theta)/d(Z) approximates a vertical line on typical SkewT chart
3483!.. near typical (mid-latitude) tropopause height. Since messy data
3484!.. could give us a false signal of such a transition, do the check over
3485!.. three K-level change, not just a level-to-level check. This method
3486!.. has potential failure in arctic-like conditions with extremely low
3487!.. tropopause height, as would any other diagnostic, so ensure resulting
3488!.. k_tropo level is above 700hPa.
3489
3490 if ( (kte-k_p200) .lt. 3) k_p200 = kte-3
3491 DO k = k_p200-2, kts, -1
3492 theta1 = theta(k)
3493 theta2 = theta(k+2)
3494 delz = 0.5*dz1d(k) + dz1d(k+1) + 0.5*dz1d(k+2)
3495 if ( (((theta2-theta1)/delz).lt.10./1500.) .OR. &
3496 & p1d(k).gt.70000.) EXIT
3497 ENDDO
3498 k_tropo = max(kts+2, min(k+2, kte-1))
3499
3500 if (k_tropo .gt. k_p200) then
3501 DO k = kte-3, k_p200-2, -1
3502 theta1 = theta(k)
3503 theta2 = theta(k+2)
3504 delz = 0.5*dz1d(k) + dz1d(k+1) + 0.5*dz1d(k+2)
3505 if ( (((theta2-theta1)/delz).lt.10./1500.) .AND. &
3506 & p1d(k).gt.9000.) EXIT
3507 ENDDO
3508 k_tropo = max(k_p200-1, min(k+2, kte-1))
3509 endif
3510
3511!..Eliminate possible fractional clouds above supposed tropopause.
3512 DO k = k_tropo+1, kte
3513 if (cfr1d(k).gt.0.0 .and. cfr1d(k).lt.1.0) then
3514 cfr1d(k) = 0.
3515 endif
3516 ENDDO
3517
3518!..Be a bit more conservative with lower cloud fraction in scenario with
3519!.. well-mixed convective boundary layer below LCL.
3520
3521 kbot = kts+1
3522 DO k = kbot, k_m12c
3523 if ( (theta(k)-theta(k-1)) .gt. 0.010e-3*dz1d(k)) EXIT
3524 ENDDO
3525 kbot = max(kts+1, k-2)
3526 DO k = kts, kbot
3527 if (cfr1d(k).gt.0.0 .and. cfr1d(k).lt.1.0) &
3528 & cfr1d(k) = max(0.01,0.33*cfr1d(k))
3529 ENDDO
3530 DO k = kts,k_tropo
3531 if (cfr1d(k).gt.0.0) kbot = min(k,kbot)
3532 ENDDO
3533
3534!..Starting below tropo height, if cloud fraction greater than 1 percent,
3535!.. compute an approximate total layer depth of cloud, determine a total
3536!.. liquid water/ice path (LWP/IWP), then reduce that amount with tuning
3537!.. parameter to represent entrainment factor, then divide up LWP/IWP
3538!.. into delta-Z weighted amounts for individual levels per cloud layer.
3539
3540 k_cldb = k_tropo
3541 in_cloud = .false.
3542 k = k_tropo
3543 DO WHILE (.not. in_cloud .AND. k.gt.k_m12c+1)
3544 k_cldt = 0
3545 if (cfr1d(k).ge.0.01) then
3546 in_cloud = .true.
3547 k_cldt = max(k_cldt, k)
3548 endif
3549 if (in_cloud) then
3550 DO k2 = k_cldt-1, k_m12c, -1
3551 if (cfr1d(k2).lt.0.01 .or. k2.eq.k_m12c) then
3552 k_cldb = k2+1
3553 goto 87
3554 endif
3555 ENDDO
3556 87 continue
3557 in_cloud = .false.
3558 endif
3559 if ((k_cldt - k_cldb + 1) .ge. 2) then
3560 call adjust_cloudice(cfr1d, qi1d, qs1d, qvs1d, t1d, dz1d, &
3561 & entrmnt, k_cldb,k_cldt,kts,kte)
3562 k = k_cldb
3563 elseif ((k_cldt - k_cldb + 1) .eq. 1) then
3564 if (cfr1d(k_cldb).gt.0.and.cfr1d(k_cldb).lt.1.) &
3565 & qi1d(k_cldb)=qi1d(k_cldb)+0.05*qvs1d(k_cldb)*cfr1d(k_cldb)
3566 k = k_cldb
3567 endif
3568 k = k - 1
3569 ENDDO
3570
3571 k_cldb = k_m12c + 3
3572 in_cloud = .false.
3573 k = k_m12c + 2
3574 DO WHILE (.not. in_cloud .AND. k.gt.kbot)
3575 k_cldt = 0
3576 if (cfr1d(k).ge.0.01) then
3577 in_cloud = .true.
3578 k_cldt = max(k_cldt, k)
3579 endif
3580 if (in_cloud) then
3581 DO k2 = k_cldt-1, kbot, -1
3582 if (cfr1d(k2).lt.0.01 .or. k2.eq.kbot) then
3583 k_cldb = k2+1
3584 goto 88
3585 endif
3586 ENDDO
3587 88 continue
3588 in_cloud = .false.
3589 endif
3590 if ((k_cldt - k_cldb + 1) .ge. 2) then
3591 call adjust_cloudh2o(cfr1d, qc1d, qvs1d, t1d, dz1d, &
3592 & entrmnt, k_cldb,k_cldt,kts,kte)
3593 k = k_cldb
3594 elseif ((k_cldt - k_cldb + 1) .eq. 1) then
3595 if (cfr1d(k_cldb).gt.0.and.cfr1d(k_cldb).lt.1.) &
3596 & qc1d(k_cldb)=qc1d(k_cldb)+0.05*qvs1d(k_cldb)*cfr1d(k_cldb)
3597 k = k_cldb
3598 endif
3599 k = k - 1
3600 ENDDO
3601
3602
3603 END SUBROUTINE find_cloudlayers
3604
3605!+---+-----------------------------------------------------------------+
3606
3608 SUBROUTINE adjust_cloudice(cfr,qi,qs,qvs,T,dz,entr, k1,k2,kts,kte)
3609!
3610 IMPLICIT NONE
3611!
3612 INTEGER, INTENT(IN):: k1,k2, kts,kte
3613 REAL, INTENT(IN):: entr
3614 REAL, DIMENSION(kts:kte), INTENT(IN):: cfr, qs, qvs, t, dz
3615 REAL, DIMENSION(kts:kte), INTENT(INOUT):: qi
3616 REAL:: iwc, max_iwc, tdz, this_iwc, this_dz
3617 INTEGER:: k
3618
3619 tdz = 0.
3620 do k = k1, k2
3621 tdz = tdz + dz(k)
3622 enddo
3623! max_iwc = ABS(qvs(k2)-qvs(k1))
3624 max_iwc = max(0.0, qvs(k1)-qvs(k2))
3625
3626 do k = k1, k2
3627 max_iwc = max(1.e-6, max_iwc - (qi(k)+qs(k)))
3628 enddo
3629 max_iwc = min(1.e-4, max_iwc)
3630
3631 this_dz = 0.0
3632 do k = k1, k2
3633 if (k.eq.k1) then
3634 this_dz = this_dz + 0.5*dz(k)
3635 else
3636 this_dz = this_dz + dz(k)
3637 endif
3638 this_iwc = max_iwc*this_dz/tdz
3639 iwc = max(1.e-6, this_iwc*(1.-entr))
3640 if (cfr(k).gt.0.0.and.cfr(k).lt.1.0.and.t(k).ge.203.16) then
3641 qi(k) = qi(k) + cfr(k)*cfr(k)*iwc
3642 endif
3643 enddo
3644
3645 END SUBROUTINE adjust_cloudice
3646
3647!+---+-----------------------------------------------------------------+
3648
3650 SUBROUTINE adjust_cloudh2o(cfr, qc, qvs,T,dz,entr, k1,k2,kts,kte)
3651!
3652 IMPLICIT NONE
3653!
3654 INTEGER, INTENT(IN):: k1,k2, kts,kte
3655 REAL, INTENT(IN):: entr
3656 REAL, DIMENSION(kts:kte), INTENT(IN):: cfr, qvs, t, dz
3657 REAL, DIMENSION(kts:kte), INTENT(INOUT):: qc
3658 REAL:: lwc, max_lwc, tdz, this_lwc, this_dz
3659 INTEGER:: k
3660
3661 tdz = 0.
3662 do k = k1, k2
3663 tdz = tdz + dz(k)
3664 enddo
3665! max_lwc = ABS(qvs(k2)-qvs(k1))
3666 max_lwc = max(0.0, qvs(k1)-qvs(k2))
3667! print*, ' max_lwc = ', max_lwc, ' over DZ=',tdz
3668
3669 do k = k1, k2
3670 max_lwc = max(1.e-6, max_lwc - qc(k))
3671 enddo
3672 max_lwc = min(1.e-4, max_lwc)
3673
3674 this_dz = 0.0
3675 do k = k1, k2
3676 if (k.eq.k1) then
3677 this_dz = this_dz + 0.5*dz(k)
3678 else
3679 this_dz = this_dz + dz(k)
3680 endif
3681 this_lwc = max_lwc*this_dz/tdz
3682 lwc = max(1.e-6, this_lwc*(1.-entr))
3683 if (cfr(k).gt.0.0.and.cfr(k).lt.1.0.and.t(k).ge.258.16) then
3684 qc(k) = qc(k) + cfr(k)*cfr(k)*lwc
3685 endif
3686 enddo
3687
3688 END SUBROUTINE adjust_cloudh2o
3689
3690!+---+-----------------------------------------------------------------+
3691
3694 SUBROUTINE adjust_cloudfinal(cfr, qc, qi, Rho,dz, kts,kte)
3695!
3696 IMPLICIT NONE
3697!
3698 INTEGER, INTENT(IN):: kts,kte
3699 REAL, DIMENSION(kts:kte), INTENT(IN):: cfr, rho, dz
3700 REAL, DIMENSION(kts:kte), INTENT(INOUT):: qc, qi
3701 REAL:: lwp, iwp, xfac
3702 INTEGER:: k
3703
3704 lwp = 0.
3705 iwp = 0.
3706 do k = kts, kte
3707 if (cfr(k).gt.0.0 .and. cfr(k).lt.1.0) then
3708 lwp = lwp + qc(k)*rho(k)*dz(k)
3709 iwp = iwp + qi(k)*rho(k)*dz(k)
3710 endif
3711 enddo
3712
3713 if (lwp .gt. 1.0) then
3714 xfac = 1.0/lwp
3715 do k = kts, kte
3716 if (cfr(k).gt.0.0 .and. cfr(k).lt.1.0) then
3717 qc(k) = qc(k)*xfac
3718 endif
3719 enddo
3720 endif
3721
3722 if (iwp .gt. 1.0) then
3723 xfac = 1.0/iwp
3724 do k = kts, kte
3725 if (cfr(k).gt.0.0 .and. cfr(k).lt.1.0) then
3726 qi(k) = qi(k)*xfac
3727 endif
3728 enddo
3729 endif
3730
3731 END SUBROUTINE adjust_cloudfinal
3732
3735 & ( ix, nlay, xr_con, xr_exp, plyr, clwf, rhly, qstl, & ! --- inputs
3736 & cldtot ) & ! --- outputs
3737
3738! --- inputs:
3739 integer, intent(in) :: IX, NLAY
3740 real (kind=kind_phys), intent(in) :: xr_con, xr_exp
3741 real (kind=kind_phys), dimension(:,:), intent(in) :: plyr, clwf, &
3742 & rhly, qstl
3743
3744! --- outputs
3745 real (kind=kind_phys), dimension(:,:), intent(inout) :: cldtot
3746
3747! --- local variables:
3748
3749 real (kind=kind_phys) :: clwmin, clwm, clwt, onemrh, value, &
3750 & tem1, tem2
3751 integer :: i, k
3752
3754
3755 clwmin = 0.0
3756 do k = 1, nlay
3757 do i = 1, ix
3758 clwt = 1.0e-6 * (plyr(i,k)*0.001)
3759
3760 if (clwf(i,k) > clwt) then
3761
3762 onemrh= max( 1.e-10, 1.0-rhly(i,k) )
3763 clwm = clwmin / max( 0.01, plyr(i,k)*0.001 )
3764
3765 tem1 = min(max((onemrh*qstl(i,k))**xr_exp,0.0001),1.0)
3766 tem1 = xr_con / tem1
3767
3768 value = max( min( tem1*(clwf(i,k)-clwm), 50.0 ), 0.0 )
3769 tem2 = sqrt(sqrt(rhly(i,k)))
3770
3771 cldtot(i,k) = max( tem2*(1.0-exp(-value)), 0.0 )
3772 endif
3773 enddo
3774 enddo
3775
3776 end subroutine cloud_fraction_xurandall
3777
3780 & ( ix, nlay, lmfdeep2, xrc3, xr_exp, plyr, clwf, rhly, qstl, & ! --- inputs
3781 & cldtot ) & ! --- outputs
3782
3783! --- inputs:
3784 integer, intent(in) :: IX, NLAY
3785 real (kind=kind_phys), intent(in) :: xrc3, xr_exp
3786 real (kind=kind_phys), dimension(:,:), intent(in) :: plyr, clwf, &
3787 & rhly, qstl
3788 logical, intent(in) :: lmfdeep2
3789
3790! --- outputs
3791 real (kind=kind_phys), dimension(:,:), intent(inout) :: cldtot
3792
3793! --- local variables:
3794
3795 real (kind=kind_phys) :: clwmin, clwm, clwt, onemrh, value, &
3796 & tem1, tem2
3797 integer :: i, k
3798
3800
3801 clwmin = 0.0
3802 do k = 1, nlay
3803 do i = 1, ix
3804 clwt = 1.0e-6 * (plyr(i,k)*0.001)
3805! clwt = 2.0e-6 * (plyr(i,k)*0.001)
3806
3807 if (clwf(i,k) > clwt) then
3808 onemrh= max( 1.e-10, 1.0-rhly(i,k) )
3809 clwm = clwmin / max( 0.01, plyr(i,k)*0.001 )
3810!
3811 tem1 = min(max((onemrh*qstl(i,k))**xr_exp,0.0001),1.0) !jhan
3812 if (lmfdeep2) then
3813 tem1 = xrc3 / tem1
3814 else
3815 tem1 = 100.0 / tem1
3816 endif
3817!
3818 value = max( min( tem1*(clwf(i,k)-clwm), 50.0 ), 0.0 )
3819 tem2 = sqrt( sqrt(rhly(i,k)) )
3820
3821 cldtot(i,k) = max( tem2*(1.0-exp(-value)), 0.0 )
3822 endif
3823 enddo
3824 enddo
3825
3826 end subroutine cloud_fraction_mass_flx_1
3827
3830 & ( ix, nlay, lmfdeep2, xrc3, xr_exp, plyr, clwf, rhly, qstl, & ! --- inputs
3831 & cldtot ) & ! --- outputs
3832
3833! --- inputs:
3834 integer, intent(in) :: IX, NLAY
3835 real (kind=kind_phys), intent(in) :: xrc3, xr_exp
3836 real (kind=kind_phys), dimension(:,:), intent(in) :: plyr, clwf, &
3837 & rhly, qstl
3838 logical, intent(in) :: lmfdeep2
3839
3840! --- outputs
3841 real (kind=kind_phys), dimension(:,:), intent(inout) :: cldtot
3842
3843! --- local variables:
3844
3845 real (kind=kind_phys) :: clwmin, clwm, clwt, onemrh, value, &
3846 & tem1, tem2
3847 integer :: i, k
3848
3850
3851 clwmin = 0.0
3852 do k = 1, nlay-1
3853 do i = 1, ix
3854 clwt = 1.0e-6 * (plyr(i,k)*0.001)
3855
3856 if (clwf(i,k) > clwt) then
3857 if(rhly(i,k) > 0.99) then
3858 cldtot(i,k) = 1.
3859 else
3860 onemrh= max( 1.e-10, 1.0-rhly(i,k) )
3861 clwm = clwmin / max( 0.01, plyr(i,k)*0.001 )
3862
3863 tem1 = min(max((onemrh*qstl(i,k))**xr_exp,0.0001),1.0)
3864 if (lmfdeep2) then
3865 tem1 = xrc3 / tem1
3866 else
3867 tem1 = 100.0 / tem1
3868 endif
3869
3870 value = max( min( tem1*(clwf(i,k)-clwm), 50.0 ), 0.0 )
3871 tem2 = sqrt( sqrt(rhly(i,k)) )
3872
3873 cldtot(i,k) = max( tem2*(1.0-exp(-value)), 0.0 )
3874 endif
3875 else
3876 cldtot(i,k) = 0.0
3877 endif
3878 enddo
3879 enddo
3880
3881 end subroutine cloud_fraction_mass_flx_2
3882!........................................!
3883 end module module_radiation_clouds
real function rsif(p, t)
THIS FUNCTION CALCULATES THE ICE SATURATION VAPOR MIXING RATIO AS A FUNCTION OF TEMPERATURE AND PRESS...
real function rslf(p, t)
THIS FUNCTION CALCULATES THE LIQUID SATURATION VAPOR MIXING RATIO AS A FUNCTION OF TEMPERATURE AND PR...
real(kind=kind_phys), dimension(nk_clds+1, 2), save ptopc
pressure limits of cloud domain interfaces (low, mid, high) in mb (0.1kPa)
subroutine, public progcld_thompson(plyr, plvl, tlyr, qlyr, qstl, rhly, clw, xlat, xlon, slmsk, dz, delp, ntrac, ntcw, ntiw, ntrw, ntsw, ntgl, ix, nlay, nlp1, uni_cld, lmfshal, lmfdeep2, cldcov, re_cloud, re_ice, re_snow, lwp_ex, iwp_ex, lwp_fc, iwp_fc, dzlay, gridkm, top_at_1, cldtot, cldcnv, cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, cld_rwp, cld_rerain, cld_swp, cld_resnow)
This subroutine added by G. Thompson specifically to account for explicit (microphysics-produced) clo...
subroutine cloud_fraction_xurandall(ix, nlay, xr_con, xr_exp, plyr, clwf, rhly, qstl, cldtot)
This subroutine computes the Xu-Randall cloud fraction scheme.
real(kind=kind_phys), parameter reliq_def
default liq radius to 10 micron
subroutine, public progcld_thompson_wsm6(plyr, plvl, tlyr, qlyr, qstl, rhly, clw, xlat, xlon, slmsk, dz, delp, ntrac, ntcw, ntiw, ntrw, ntsw, ntgl, con_ttp, xr_cnvcld, ix, nlay, nlp1, xr_con, xr_exp, uni_cld, lmfshal, lmfdeep2, cldcov, cnvw, re_cloud, re_ice, re_snow, lwp_ex, iwp_ex, lwp_fc, iwp_fc, dzlay, cldtot, cldcnv, lcnorm, cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, cld_rwp, cld_rerain, cld_swp, cld_resnow)
This subroutine is used by Thompson/WSM6/NSSL cloud microphysics (EMC)
integer, parameter, public nk_clds
number of cloud vertical domains
real(kind=kind_phys), parameter rsnow_def
default snow radius to 250 micron
real(kind=kind_phys), parameter cldssa_def
default cld single scat albedo
real(kind=kind_phys), parameter climit2
integer llyr
upper limit of boundary layer clouds
real(kind=kind_phys), dimension(95), parameter retab
real(kind=kind_phys) gord
subroutine, public adjust_cloudice(cfr, qi, qs, qvs, t, dz, entr, k1, k2, kts, kte)
real(kind=kind_phys), parameter ovcst
real(kind=kind_phys) gfac
subroutine cloud_fraction_mass_flx_2(ix, nlay, lmfdeep2, xrc3, xr_exp, plyr, clwf, rhly, qstl, cldtot)
subroutine, public gethml(plyr, ptop1, cldtot, cldcnv, dz, de_lgth, alpha, ix, nlay, iovr, iovr_rand, iovr_maxrand, iovr_max, iovr_dcorr, iovr_exp, iovr_exprand, top_at_1, si, clds, mtop, mbot)
This subroutine computes high, mid, low, total, and boundary cloud fractions and cloud top/bottom lay...
subroutine, public progcld_zhao_carr_pdf(plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, clw, cnvw, cnvc, xlat, xlon, slmsk, dz, delp, ix, nlay, nlp1, deltaq, sup, kdt, me, dzlay, cldtot, cldcnv, lcrick, lcnorm, con_thgni, con_ttp, cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, cld_rwp, cld_rerain, cld_swp, cld_resnow)
This subroutine computes cloud related quantities using zhao/moorthi's prognostic cloud microphysics ...
integer, parameter, public nf_clds
number of fields in cloud array
subroutine, public progcld_gfdl_lin(plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, clw, cnvw, cnvc, xlat, xlon, slmsk, cldtot, dz, delp, ix, nlay, nlp1, dzlay, cldtot1, cldcnv, lcrick, lcnorm, con_ttp, cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, cld_rwp, cld_rerain, cld_swp, cld_resnow)
This subroutine computes cloud related quantities using GFDL Lin MP prognostic cloud microphysics sch...
real(kind=kind_phys), parameter reice_def
default ice radius to 50 micron
subroutine, public adjust_cloudh2o(cfr, qc, qvs, t, dz, entr, k1, k2, kts, kte)
real(kind=kind_phys), parameter rrain_def
default rain radius to 1000 micron
subroutine, public progcld_fer_hires(plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, clw, xlat, xlon, slmsk, dz, delp, ntrac, ntcw, ntiw, ntrw, ix, nlay, nlp1, icloud, xr_con, xr_exp, uni_cld, lmfshal, lmfdeep2, cldcov, re_cloud, re_ice, re_snow, dzlay, cldtot, cldcnv, lcnorm, cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, cld_rwp, cld_rerain, cld_swp, cld_resnow)
This subroutine computes cloud related quantities using Ferrier-Aligo cloud microphysics scheme.
subroutine, public adjust_cloudfinal(cfr, qc, qi, rho, dz, kts, kte)
Do not alter any grid-explicitly resolved hydrometeors, rather only the supposed amounts due to the c...
real(kind=kind_phys), parameter cldasy_def
default cld asymmetry factor
subroutine, public cld_init(si, nlay, imp_physics, me, con_g, con_rd, errflg, errmsg)
This subroutine is an initialization program for cloud-radiation calculations and sets up boundary la...
subroutine, public progcld_zhao_carr(plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, clw, xlat, xlon, slmsk, dz, delp, ix, nlay, nlp1, uni_cld, lmfshal, lmfdeep2, xr_con, xr_exp, cldcov, effrl, effri, effrr, effrs, effr_in, dzlay, cldtot, cldcnv, lcrick, lcnorm, con_ttp, cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, cld_rwp, cld_rerain, cld_swp, cld_resnow)
This subroutine computes cloud related quantities using zhao/moorthi's prognostic cloud microphysics ...
subroutine, public progclduni(plyr, plvl, tlyr, tvly, ccnd, ncnd, xlat, xlon, slmsk, dz, delp, ix, nlay, nlp1, cldtot, effrl, effri, effrr, effrs, effr_in, dzlay, cldtot1, cldcnv, lcrick, lcnorm, con_ttp, cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, cld_rwp, cld_rerain, cld_swp, cld_resnow)
This subroutine computes cloud related quantities using for unified cloud microphysics scheme.
real(kind=kind_phys), parameter creice_def
default convective ice radius to 25 micron overland
real(kind=kind_phys), parameter climit
subroutine cloud_fraction_mass_flx_1(ix, nlay, lmfdeep2, xrc3, xr_exp, plyr, clwf, rhly, qstl, cldtot)
subroutine, public radiation_clouds_prop(plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, ccnd, ncndl, cnvw, cnvc, tracer1, xlat, xlon, slmsk, dz, delp, ix, lm, nlay, nlp1, deltaq, sup, dcorr_con, me, icloud, kdt, ntrac, ntcw, ntiw, ntrw, ntsw, ntgl, ntclamt, imp_physics, imp_physics_nssl, imp_physics_fer_hires, imp_physics_gfdl, imp_physics_thompson, imp_physics_wsm6, imp_physics_tempo, imp_physics_zhao_carr, imp_physics_zhao_carr_pdf, imp_physics_mg, iovr, iovr_rand, iovr_maxrand, iovr_max, iovr_dcorr, iovr_exp, iovr_exprand, idcor, idcor_con, idcor_hogan, idcor_oreopoulos, lcrick, lcnorm, imfdeepcnv, imfdeepcnv_gf, imfdeepcnv_c3, do_mynnedmf, lgfdlmprad, xr_cnvcld, uni_cld, lmfshal, lmfdeep2, cldcov, clouds1, effrl, effri, effrr, effrs, effr_in, effrl_inout, effri_inout, effrs_inout, lwp_ex, iwp_ex, lwp_fc, iwp_fc, dzlay, latdeg, julian, yearlen, gridkm, top_at_1, si, xr_con, xr_exp, con_ttp, con_pi, con_g, con_rd, con_thgni, cld_frac, cld_lwp, cld_reliq, cld_iwp, cld_reice, cld_rwp, cld_rerain, cld_swp, cld_resnow, clds, mtop, mbot, de_lgth, alpha)
Subroutine radiation_clouds_prop computes cloud related quantities for different cloud microphysics s...
subroutine, public cal_cldfra3(cldfra, qv, qc, qi, qs, dz, p, t, xland, gridkm, modify_qvapor, max_relh, kts, kte, debug_flag)
Cloud fraction scheme by G. Thompson (NCAR-RAL), not intended for combining with any cumulus or shall...
subroutine, public find_cloudlayers(qvs1d, cfr1d, t1d, p1d, dz1d, entrmnt, debugfl, qc1d, qi1d, qs1d, kts, kte)
From cloud fraction array, find clouds of multi-level depth and compute a reasonable value of LWP or ...
this module defines fortran unit numbers for input/output data files for the ncep gfs model.
Definition iounitdef.f:53
This module computes the moisture tendencies of water vapor, cloud droplets, rain,...
This module contains the calculation of cloud overlap parameters for both RRTMG and RRTMGP.
This module computes cloud related quantities for radiation computations.