Line data Source code
1 : module cloud_fraction
2 :
3 : ! Cloud fraction parameterization.
4 :
5 :
6 : use shr_kind_mod, only: r8 => shr_kind_r8
7 : use ppgrid, only: pcols, pver, pverp
8 : use ref_pres, only: pref_mid
9 : use spmd_utils, only: masterproc
10 : use cam_logfile, only: iulog
11 : use cam_abortutils, only: endrun
12 : use ref_pres, only: trop_cloud_top_lev
13 :
14 : implicit none
15 : private
16 : save
17 :
18 : ! Public interfaces
19 : public &
20 : cldfrc_readnl, &! read cldfrc_nl namelist
21 : cldfrc_register, &! add fields to pbuf
22 : cldfrc_init, &! Inititialization of cloud_fraction run-time parameters
23 : cldfrc_getparams, &! public access of tuning parameters
24 : dp1, &! parameter for deep convection cloud fraction needed in clubb_intr
25 : dp2 ! parameter for deep convection cloud fraction needed in clubb_intr
26 :
27 : ! Private data
28 : real(r8), parameter :: unset_r8 = huge(1.0_r8)
29 :
30 : ! Top level
31 : integer :: top_lev = 1
32 :
33 : ! Physics buffer indices
34 : integer :: sh_frac_idx = 0
35 : integer :: dp_frac_idx = 0
36 :
37 : ! Namelist variables
38 : logical :: cldfrc_freeze_dry ! switch for Vavrus correction
39 : logical :: cldfrc_ice ! switch to compute ice cloud fraction
40 : real(r8) :: cldfrc_rhminl = unset_r8 ! minimum rh for low stable clouds
41 : real(r8) :: cldfrc_rhminl_adj_land = unset_r8 ! rhminl adjustment for snowfree land
42 : real(r8) :: cldfrc_rhminh = unset_r8 ! minimum rh for high stable clouds
43 : real(r8) :: cldfrc_sh1 = unset_r8 ! parameter for shallow convection cloud fraction
44 : real(r8) :: cldfrc_sh2 = unset_r8 ! parameter for shallow convection cloud fraction
45 : real(r8) :: cldfrc_dp1 = unset_r8 ! parameter for deep convection cloud fraction
46 : real(r8) :: cldfrc_dp2 = unset_r8 ! parameter for deep convection cloud fraction
47 : real(r8) :: cldfrc_premit = unset_r8 ! top pressure bound for mid level cloud
48 : real(r8) :: cldfrc_premib = unset_r8 ! bottom pressure bound for mid level cloud
49 : integer :: cldfrc_iceopt ! option for ice cloud closure
50 : ! 1=wang & sassen 2=schiller (iciwc)
51 : ! 3=wood & field, 4=Wilson (based on smith)
52 : real(r8) :: cldfrc_icecrit = unset_r8 ! Critical RH for ice clouds in Wilson & Ballard closure (smaller = more ice clouds)
53 :
54 : real(r8) :: rhminl ! set from namelist input cldfrc_rhminl
55 : real(r8) :: rhminl_adj_land ! set from namelist input cldfrc_rhminl_adj_land
56 : real(r8) :: rhminh ! set from namelist input cldfrc_rhminh
57 : real(r8) :: sh1, sh2 ! set from namelist input cldfrc_sh1, cldfrc_sh2
58 : real(r8) :: dp1,dp2 ! set from namelist input cldfrc_dp1, cldfrc_dp2
59 : real(r8) :: premit ! set from namelist input cldfrc_premit
60 : real(r8) :: premib ! set from namelist input cldfrc_premib
61 : integer :: iceopt ! set from namelist input cldfrc_iceopt
62 : real(r8) :: icecrit ! set from namelist input cldfrc_icecrit
63 :
64 : ! constants
65 : real(r8), parameter :: pnot = 1.e5_r8 ! reference pressure
66 : real(r8), parameter :: lapse = 6.5e-3_r8 ! U.S. Standard Atmosphere lapse rate
67 : real(r8), parameter :: pretop = 1.0e2_r8 ! pressure bounding high cloud
68 :
69 : integer count
70 :
71 : logical :: inversion_cld_off ! Turns off stratification-based cld frc
72 :
73 : integer :: k700 ! model level nearest 700 mb
74 :
75 : !================================================================================================
76 : contains
77 : !================================================================================================
78 :
79 2048 : subroutine cldfrc_readnl(nlfile)
80 :
81 : use namelist_utils, only: find_group_name
82 : use units, only: getunit, freeunit
83 : use mpishorthand
84 :
85 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
86 :
87 : ! Local variables
88 : integer :: unitn, ierr
89 : character(len=*), parameter :: subname = 'cldfrc_readnl'
90 :
91 : namelist /cldfrc_nl/ cldfrc_freeze_dry, cldfrc_ice, cldfrc_rhminl, &
92 : cldfrc_rhminl_adj_land, cldfrc_rhminh, cldfrc_sh1, &
93 : cldfrc_sh2, cldfrc_dp1, cldfrc_dp2, &
94 : cldfrc_premit, cldfrc_premib, cldfrc_iceopt, &
95 : cldfrc_icecrit
96 : !-----------------------------------------------------------------------------
97 :
98 1026 : if (masterproc) then
99 2 : unitn = getunit()
100 2 : open( unitn, file=trim(nlfile), status='old' )
101 2 : call find_group_name(unitn, 'cldfrc_nl', status=ierr)
102 2 : if (ierr == 0) then
103 2 : read(unitn, cldfrc_nl, iostat=ierr)
104 2 : if (ierr /= 0) then
105 0 : call endrun(subname // ':: ERROR reading namelist')
106 : end if
107 : end if
108 2 : close(unitn)
109 2 : call freeunit(unitn)
110 :
111 : ! set local variables
112 2 : rhminl = cldfrc_rhminl
113 2 : rhminl_adj_land = cldfrc_rhminl_adj_land
114 2 : rhminh = cldfrc_rhminh
115 2 : sh1 = cldfrc_sh1
116 2 : sh2 = cldfrc_sh2
117 2 : dp1 = cldfrc_dp1
118 2 : dp2 = cldfrc_dp2
119 2 : premit = cldfrc_premit
120 2 : premib = cldfrc_premib
121 2 : iceopt = cldfrc_iceopt
122 2 : icecrit = cldfrc_icecrit
123 :
124 : end if
125 :
126 : #ifdef SPMD
127 : ! Broadcast namelist variables
128 1024 : call mpibcast(cldfrc_freeze_dry, 1, mpilog, 0, mpicom)
129 1024 : call mpibcast(cldfrc_ice, 1, mpilog, 0, mpicom)
130 1024 : call mpibcast(rhminl, 1, mpir8, 0, mpicom)
131 1024 : call mpibcast(rhminl_adj_land, 1, mpir8, 0, mpicom)
132 1024 : call mpibcast(rhminh, 1, mpir8, 0, mpicom)
133 1024 : call mpibcast(sh1 , 1, mpir8, 0, mpicom)
134 1024 : call mpibcast(sh2 , 1, mpir8, 0, mpicom)
135 1024 : call mpibcast(dp1 , 1, mpir8, 0, mpicom)
136 1024 : call mpibcast(dp2 , 1, mpir8, 0, mpicom)
137 1024 : call mpibcast(premit, 1, mpir8, 0, mpicom)
138 1024 : call mpibcast(premib, 1, mpir8, 0, mpicom)
139 1024 : call mpibcast(iceopt, 1, mpiint, 0, mpicom)
140 1024 : call mpibcast(icecrit, 1, mpir8, 0, mpicom)
141 : #endif
142 :
143 1024 : end subroutine cldfrc_readnl
144 :
145 : !================================================================================================
146 :
147 1024 : subroutine cldfrc_register
148 :
149 : ! Register fields in the physics buffer.
150 :
151 : use physics_buffer, only : pbuf_add_field, dtype_r8
152 :
153 : !-----------------------------------------------------------------------
154 :
155 1024 : call pbuf_add_field('SH_FRAC', 'physpkg', dtype_r8, (/pcols,pver/), sh_frac_idx)
156 1024 : call pbuf_add_field('DP_FRAC', 'physpkg', dtype_r8, (/pcols,pver/), dp_frac_idx)
157 :
158 1024 : end subroutine cldfrc_register
159 :
160 : !================================================================================================
161 :
162 1024 : subroutine cldfrc_getparams(rhminl_out, rhminl_adj_land_out, rhminh_out, premit_out, &
163 : premib_out, iceopt_out, icecrit_out)
164 : !-----------------------------------------------------------------------
165 : ! Purpose: Return cldfrc tuning parameters
166 : !-----------------------------------------------------------------------
167 :
168 : real(r8), intent(out), optional :: rhminl_out
169 : real(r8), intent(out), optional :: rhminl_adj_land_out
170 : real(r8), intent(out), optional :: rhminh_out
171 : real(r8), intent(out), optional :: premit_out
172 : real(r8), intent(out), optional :: premib_out
173 : integer, intent(out), optional :: iceopt_out
174 : real(r8), intent(out), optional :: icecrit_out
175 :
176 1024 : if ( present(rhminl_out) ) rhminl_out = rhminl
177 1024 : if ( present(rhminl_adj_land_out) ) rhminl_adj_land_out = rhminl_adj_land
178 1024 : if ( present(rhminh_out) ) rhminh_out = rhminh
179 1024 : if ( present(premit_out) ) premit_out = premit
180 1024 : if ( present(premib_out) ) premib_out = premib
181 1024 : if ( present(iceopt_out) ) iceopt_out = iceopt
182 1024 : if ( present(icecrit_out) ) icecrit_out = icecrit
183 :
184 1024 : end subroutine cldfrc_getparams
185 :
186 : !===============================================================================
187 :
188 1024 : subroutine cldfrc_init
189 :
190 : ! Initialize cloud fraction run-time parameters
191 :
192 : use cam_history, only: addfld
193 : use phys_control, only: phys_getopts
194 :
195 : use convective_cloud_cover, only: convective_cloud_cover_init
196 : use compute_cloud_fraction, only: compute_cloud_fraction_init
197 :
198 : ! query interfaces for scheme settings
199 : character(len=16) :: shallow_scheme, eddy_scheme, macrop_scheme
200 :
201 : integer :: k
202 :
203 : ! outputs from CCPPized scheme
204 : character(len=512) :: errmsg
205 : integer :: errflg
206 :
207 : !-----------------------------------------------------------------------------
208 :
209 : call phys_getopts(shallow_scheme_out = shallow_scheme ,&
210 : eddy_scheme_out = eddy_scheme ,&
211 1024 : macrop_scheme_out = macrop_scheme )
212 :
213 : ! Limit CAM5 cloud physics to below top cloud level.
214 1024 : if ( .not. macrop_scheme == "rk") top_lev = trop_cloud_top_lev
215 :
216 : ! Turn off inversion_cld if any UW PBL scheme is being used
217 1024 : if ( eddy_scheme .eq. 'diag_TKE' .or. shallow_scheme .eq. 'UW' ) then
218 0 : inversion_cld_off = .true.
219 : else
220 1024 : inversion_cld_off = .false.
221 : endif
222 :
223 1024 : if (pref_mid(top_lev) > 7.e4_r8) &
224 0 : call endrun ('cldfrc_init: model levels bracketing 700 mb not found')
225 :
226 : ! Find vertical level nearest 700 mb.
227 27648 : k700 = minloc(abs(pref_mid(top_lev:pver) - 7.e4_r8), 1)
228 :
229 1024 : if (masterproc) then
230 2 : write(iulog,*)'cldfrc_init: model level nearest 700 mb is',k700,'which is',pref_mid(k700),'pascals'
231 : end if
232 :
233 2048 : call addfld ('SH_CLD', (/ 'lev' /), 'A', 'fraction', 'Shallow convective cloud cover' )
234 2048 : call addfld ('DP_CLD', (/ 'lev' /), 'A', 'fraction', 'Deep convective cloud cover' )
235 :
236 : ! populate namelist parameters in CCPPized schemes
237 : call convective_cloud_cover_init( &
238 : amIRoot = masterproc, &
239 : iulog = iulog, &
240 : sh1_in = sh1, &
241 : sh2_in = sh2, &
242 : dp1_in = dp1, &
243 : dp2_in = dp2, &
244 1024 : errmsg = errmsg, errflg = errflg)
245 :
246 : call compute_cloud_fraction_init( &
247 : amIRoot = masterproc, &
248 : iulog = iulog, &
249 : pver = pver, &
250 : pref_mid = pref_mid, &
251 : inversion_cld_off_in = inversion_cld_off, &
252 : cldfrc_freeze_dry_in = cldfrc_freeze_dry, &
253 : cldfrc_ice_in = cldfrc_ice, &
254 : iceopt_in = iceopt, &
255 : rhminl_in = rhminl, &
256 : rhminl_adj_land_in = rhminl_adj_land, &
257 : rhminh_in = rhminh, &
258 : premit_in = premit, &
259 : premib_in = premib, &
260 : icecrit_in = icecrit, &
261 1024 : errmsg = errmsg, errflg = errflg)
262 :
263 1024 : end subroutine cldfrc_init
264 :
265 : end module cloud_fraction
|