CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
rrtmgp_lw_cloud_optics.F90
1
3
10 use machine, only: kind_phys
11 use mo_rte_kind, only: wl
12 use mo_cloud_optics, only: ty_cloud_optics
13 use rrtmgp_lw_gas_optics, only: lw_gas_props
14 use radiation_tools, only: check_error_msg
15 use netcdf
16#ifdef MPI
17 use mpi_f08
18#endif
19
20 implicit none
21
22 type(ty_cloud_optics) :: lw_cloud_props
23 integer :: &
24 nrghice_fromfilelw, nbandlw, nsize_liqlw, nsize_icelw, nsizereglw, &
25 ncoeff_extlw, ncoeff_ssa_glw, nboundlw, npairslw
26 real(kind_phys) :: &
27 radliq_faclw, & !< Factor for calculating LUT interpolation indices for liquid
28 radice_faclw
29 real(kind_phys), dimension(:,:), allocatable :: &
30 lut_extliqlw, & !< LUT shortwave liquid extinction coefficient
31 lut_ssaliqlw, & !< LUT shortwave liquid single scattering albedo
32 lut_asyliqlw, & !< LUT shortwave liquid asymmetry parameter
33 band_limscldlw
34 real(kind_phys), dimension(:,:,:), allocatable :: &
35 lut_exticelw, & !< LUT shortwave ice extinction coefficient
36 lut_ssaicelw, & !< LUT shortwave ice single scattering albedo
37 lut_asyicelw
38 real(kind_phys), dimension(:), allocatable :: &
39 pade_sizereg_extliqlw, & !< Particle size regime boundaries for shortwave liquid extinction
40
41 pade_sizereg_ssaliqlw, &
43 pade_sizereg_asyliqlw, &
45 pade_sizereg_exticelw, &
47 pade_sizereg_ssaicelw, &
49 pade_sizereg_asyicelw
51 real(kind_phys), dimension(:,:,:), allocatable :: &
52 pade_extliqlw, & !< PADE coefficients for shortwave liquid extinction
53 pade_ssaliqlw, & !< PADE coefficients for shortwave liquid single scattering albedo
54 pade_asyliqlw
55 real(kind_phys), dimension(:,:,:,:), allocatable :: &
56 pade_exticelw, & !< PADE coefficients for shortwave ice extinction
57 pade_ssaicelw, & !< PADE coefficients for shortwave ice single scattering albedo
58 pade_asyicelw
59
60 ! Parameters used for rain and snow(+groupel) RRTMGP cloud-optics
61 real(kind_phys), parameter :: &
62 absrain = 0.33e-3, &
63 abssnow0 = 1.5, &
64 abssnow1 = 2.34e-3
65 real(kind_phys) :: &
66 radliq_lwrlw, & !< Liquid particle size lower bound for LUT interpolation
67 radliq_uprlw, & !< Liquid particle size upper bound for LUT interpolation
68 radice_lwrlw, & !< Ice particle size upper bound for LUT interpolation
69 radice_uprlw
70
71contains
72
73 ! ######################################################################################
74 ! SUBROUTINE rrtmgp_lw_cloud_optics_init()
75 ! ######################################################################################
77 subroutine rrtmgp_lw_cloud_optics_init(rrtmgp_root_dir, rrtmgp_lw_file_clouds, &
78 doGP_cldoptics_PADE, doGP_cldoptics_LUT, nrghice, mpicomm, mpirank, mpiroot, &
79 errmsg, errflg)
80
81 ! Inputs
82 character(len=128),intent(in) :: &
83 rrtmgp_root_dir, & !< RTE-RRTMGP root directory
84 rrtmgp_lw_file_clouds
85
86 logical, intent(in) :: &
87 doGP_cldoptics_PADE,& !< Use RRTMGP cloud-optics: PADE approximation?
88 doGP_cldoptics_LUT
89 integer, intent(inout) :: &
90 nrghice
91 type(mpi_comm), intent(in) :: &
92 mpicomm
93 integer, intent(in) :: &
94 mpirank, & !< Current MPI rank
95 mpiroot
96
97 ! Outputs
98 character(len=*), intent(out) :: &
99 errmsg
100 integer, intent(out) :: &
101 errflg
102
103 ! Local variables
104 integer :: dimID,varID,status,ncid,mpierr
105 character(len=264) :: lw_cloud_props_file
106
107 ! Initialize
108 errmsg = ''
109 errflg = 0
110
111 ! Filenames are set in the physics_nml
112 lw_cloud_props_file = trim(rrtmgp_root_dir)//trim(rrtmgp_lw_file_clouds)
113
114 ! #######################################################################################
115 !
116 ! Read dimensions for longwave cloud-optics fields...
117 ! (ONLY master processor(0), if MPI enabled)
118 !
119 ! #######################################################################################
120#ifdef MPI
121 if (mpirank .eq. mpiroot) then
122#endif
123 write (*,*) 'Reading RRTMGP longwave cloud-optics metadata ... '
124
125 ! Open file
126 status = nf90_open(trim(lw_cloud_props_file), nf90_nowrite, ncid)
127
128 ! Read dimensions
129 status = nf90_inq_dimid(ncid, 'nband', dimid)
130 status = nf90_inquire_dimension(ncid, dimid, len=nbandlw)
131 status = nf90_inq_dimid(ncid, 'nrghice', dimid)
132 status = nf90_inquire_dimension(ncid, dimid, len=nrghice_fromfilelw)
133 status = nf90_inq_dimid(ncid, 'nsize_liq', dimid)
134 status = nf90_inquire_dimension(ncid, dimid, len=nsize_liqlw)
135 status = nf90_inq_dimid(ncid, 'nsize_ice', dimid)
136 status = nf90_inquire_dimension(ncid, dimid, len=nsize_icelw)
137 status = nf90_inq_dimid(ncid, 'nsizereg', dimid)
138 status = nf90_inquire_dimension(ncid, dimid, len=nsizereglw)
139 status = nf90_inq_dimid(ncid, 'ncoeff_ext', dimid)
140 status = nf90_inquire_dimension(ncid, dimid, len=ncoeff_extlw)
141 status = nf90_inq_dimid(ncid, 'ncoeff_ssa_g', dimid)
142 status = nf90_inquire_dimension(ncid, dimid, len=ncoeff_ssa_glw)
143 status = nf90_inq_dimid(ncid, 'nbound', dimid)
144 status = nf90_inquire_dimension(ncid, dimid, len=nboundlw)
145 status = nf90_inq_dimid(ncid, 'pair', dimid)
146 status = nf90_inquire_dimension(ncid, dimid, len=npairslw)
147
148#ifdef MPI
149 endif ! On master processor
150
151 ! Other processors waiting...
152 call mpi_barrier(mpicomm, mpierr)
153
154 ! #######################################################################################
155 !
156 ! Broadcast dimensions...
157 ! (ALL processors)
158 !
159 ! #######################################################################################
160 call mpi_bcast(nbandlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
161 call mpi_bcast(nsize_liqlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
162 call mpi_bcast(nsize_icelw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
163 call mpi_bcast(nsizereglw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
164 call mpi_bcast(ncoeff_extlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
165 call mpi_bcast(ncoeff_ssa_glw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
166 call mpi_bcast(nboundlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
167 call mpi_bcast(npairslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
168#endif
169
170 ! Has the number of ice-roughnesses to use been provided from the namelist?
171 ! If so, override nrghice from cloud-optics file
172 if (nrghice .ne. 0) nrghice_fromfilelw = nrghice
173#ifdef MPI
174 call mpi_bcast(nrghice_fromfilelw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
175#endif
176
177 ! #######################################################################################
178 !
179 ! Allocate space for arrays...
180 ! (ALL processors)
181 !
182 ! #######################################################################################
183 if (dogp_cldoptics_lut) then
184 allocate(lut_extliqlw(nsize_liqlw, nbandlw))
185 allocate(lut_ssaliqlw(nsize_liqlw, nbandlw))
186 allocate(lut_asyliqlw(nsize_liqlw, nbandlw))
187 allocate(lut_exticelw(nsize_icelw, nbandlw, nrghice_fromfilelw))
188 allocate(lut_ssaicelw(nsize_icelw, nbandlw, nrghice_fromfilelw))
189 allocate(lut_asyicelw(nsize_icelw, nbandlw, nrghice_fromfilelw))
190 endif
191 if (dogp_cldoptics_pade) then
192 allocate(pade_extliqlw(nbandlw, nsizereglw, ncoeff_extlw ))
193 allocate(pade_ssaliqlw(nbandlw, nsizereglw, ncoeff_ssa_glw))
194 allocate(pade_asyliqlw(nbandlw, nsizereglw, ncoeff_ssa_glw))
195 allocate(pade_exticelw(nbandlw, nsizereglw, ncoeff_extlw, nrghice_fromfilelw))
196 allocate(pade_ssaicelw(nbandlw, nsizereglw, ncoeff_ssa_glw, nrghice_fromfilelw))
197 allocate(pade_asyicelw(nbandlw, nsizereglw, ncoeff_ssa_glw, nrghice_fromfilelw))
198 allocate(pade_sizereg_extliqlw(nboundlw))
199 allocate(pade_sizereg_ssaliqlw(nboundlw))
200 allocate(pade_sizereg_asyliqlw(nboundlw))
201 allocate(pade_sizereg_exticelw(nboundlw))
202 allocate(pade_sizereg_ssaicelw(nboundlw))
203 allocate(pade_sizereg_asyicelw(nboundlw))
204 endif
205 allocate(band_limscldlw(2,nbandlw))
206
207 ! #######################################################################################
208 !
209 ! Read in data ...
210 ! (ONLY master processor(0), if MPI enabled)
211 !
212 ! #######################################################################################
213#ifdef MPI
214 if (mpirank .eq. mpiroot) then
215#endif
216 ! Read in fields from file
217 if (dogp_cldoptics_lut) then
218 write (*,*) 'Reading RRTMGP longwave cloud data (LUT) ... '
219 status = nf90_inq_varid(ncid,'radliq_lwr',varid)
220 status = nf90_get_var(ncid,varid,radliq_lwrlw)
221 status = nf90_inq_varid(ncid,'radliq_upr',varid)
222 status = nf90_get_var(ncid,varid,radliq_uprlw)
223 status = nf90_inq_varid(ncid,'radliq_fac',varid)
224 status = nf90_get_var(ncid,varid,radliq_faclw)
225 status = nf90_inq_varid(ncid,'radice_lwr',varid)
226 status = nf90_get_var(ncid,varid,radice_lwrlw)
227 status = nf90_inq_varid(ncid,'radice_upr',varid)
228 status = nf90_get_var(ncid,varid,radice_uprlw)
229 status = nf90_inq_varid(ncid,'radice_fac',varid)
230 status = nf90_get_var(ncid,varid,radice_faclw)
231 status = nf90_inq_varid(ncid,'lut_extliq',varid)
232 status = nf90_get_var(ncid,varid,lut_extliqlw)
233 status = nf90_inq_varid(ncid,'lut_ssaliq',varid)
234 status = nf90_get_var(ncid,varid,lut_ssaliqlw)
235 status = nf90_inq_varid(ncid,'lut_asyliq',varid)
236 status = nf90_get_var(ncid,varid,lut_asyliqlw)
237 status = nf90_inq_varid(ncid,'lut_extice',varid)
238 status = nf90_get_var(ncid,varid,lut_exticelw)
239 status = nf90_inq_varid(ncid,'lut_ssaice',varid)
240 status = nf90_get_var(ncid,varid,lut_ssaicelw)
241 status = nf90_inq_varid(ncid,'lut_asyice',varid)
242 status = nf90_get_var(ncid,varid,lut_asyicelw)
243 status = nf90_inq_varid(ncid,'bnd_limits_wavenumber',varid)
244 status = nf90_get_var(ncid,varid,band_limscldlw)
245 endif
246 if (dogp_cldoptics_pade) then
247 write (*,*) 'Reading RRTMGP longwave cloud data (PADE) ... '
248 status = nf90_inq_varid(ncid,'radliq_lwr',varid)
249 status = nf90_get_var(ncid,varid,radliq_lwrlw)
250 status = nf90_inq_varid(ncid,'radliq_upr',varid)
251 status = nf90_get_var(ncid,varid,radliq_uprlw)
252 status = nf90_inq_varid(ncid,'radliq_fac',varid)
253 status = nf90_get_var(ncid,varid,radliq_faclw)
254 status = nf90_inq_varid(ncid,'radice_lwr',varid)
255 status = nf90_get_var(ncid,varid,radice_lwrlw)
256 status = nf90_inq_varid(ncid,'radice_upr',varid)
257 status = nf90_get_var(ncid,varid,radice_uprlw)
258 status = nf90_inq_varid(ncid,'radice_fac',varid)
259 status = nf90_get_var(ncid,varid,radice_faclw)
260 status = nf90_inq_varid(ncid,'pade_extliq',varid)
261 status = nf90_get_var(ncid,varid,pade_extliqlw)
262 status = nf90_inq_varid(ncid,'pade_ssaliq',varid)
263 status = nf90_get_var(ncid,varid,pade_ssaliqlw)
264 status = nf90_inq_varid(ncid,'pade_asyliq',varid)
265 status = nf90_get_var(ncid,varid,pade_asyliqlw)
266 status = nf90_inq_varid(ncid,'pade_extice',varid)
267 status = nf90_get_var(ncid,varid,pade_exticelw)
268 status = nf90_inq_varid(ncid,'pade_ssaice',varid)
269 status = nf90_get_var(ncid,varid,pade_ssaicelw)
270 status = nf90_inq_varid(ncid,'pade_asyice',varid)
271 status = nf90_get_var(ncid,varid,pade_asyicelw)
272 status = nf90_inq_varid(ncid,'pade_sizreg_extliq',varid)
273 status = nf90_get_var(ncid,varid,pade_sizereg_extliqlw)
274 status = nf90_inq_varid(ncid,'pade_sizreg_ssaliq',varid)
275 status = nf90_get_var(ncid,varid,pade_sizereg_ssaliqlw)
276 status = nf90_inq_varid(ncid,'pade_sizreg_asyliq',varid)
277 status = nf90_get_var(ncid,varid,pade_sizereg_asyliqlw)
278 status = nf90_inq_varid(ncid,'pade_sizreg_extice',varid)
279 status = nf90_get_var(ncid,varid,pade_sizereg_exticelw)
280 status = nf90_inq_varid(ncid,'pade_sizreg_ssaice',varid)
281 status = nf90_get_var(ncid,varid,pade_sizereg_ssaicelw)
282 status = nf90_inq_varid(ncid,'pade_sizreg_asyice',varid)
283 status = nf90_get_var(ncid,varid,pade_sizereg_asyicelw)
284 status = nf90_inq_varid(ncid,'bnd_limits_wavenumber',varid)
285 status = nf90_get_var(ncid,varid,band_limscldlw)
286 endif
287
288 ! Close file
289 status = nf90_close(ncid)
290#ifdef MPI
291 endif ! Master process
292
293 ! Other processors waiting...
294 call mpi_barrier(mpicomm, mpierr)
295
296 ! #######################################################################################
297 !
298 ! Broadcast data...
299 ! (ALL processors)
300 !
301 ! #######################################################################################
302
303 ! Real scalars
304 call mpi_bcast(radliq_faclw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
305 call mpi_bcast(radice_faclw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
306 call mpi_bcast(radliq_lwrlw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
307 call mpi_bcast(radliq_uprlw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
308 call mpi_bcast(radice_lwrlw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
309 call mpi_bcast(radice_uprlw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
310
311 ! Real arrays
312 call mpi_bcast(band_limscldlw, size(band_limscldlw), &
313 mpi_double_precision, mpiroot, mpicomm, mpierr)
314 if (dogp_cldoptics_lut) then
315 call mpi_bcast(lut_extliqlw, size(lut_extliqlw), &
316 mpi_double_precision, mpiroot, mpicomm, mpierr)
317 call mpi_bcast(lut_ssaliqlw, size(lut_ssaliqlw), &
318 mpi_double_precision, mpiroot, mpicomm, mpierr)
319 call mpi_bcast(lut_asyliqlw, size(lut_asyliqlw), &
320 mpi_double_precision, mpiroot, mpicomm, mpierr)
321 call mpi_bcast(lut_exticelw, size(lut_exticelw), &
322 mpi_double_precision, mpiroot, mpicomm, mpierr)
323 call mpi_bcast(lut_ssaicelw, size(lut_ssaicelw), &
324 mpi_double_precision, mpiroot, mpicomm, mpierr)
325 call mpi_bcast(lut_asyicelw, size(lut_asyicelw), &
326 mpi_double_precision, mpiroot, mpicomm, mpierr)
327 endif
328 if (dogp_cldoptics_pade) then
329 call mpi_bcast(pade_extliqlw, size(pade_extliqlw), &
330 mpi_double_precision, mpiroot, mpicomm, mpierr)
331 call mpi_bcast(pade_ssaliqlw, size(pade_ssaliqlw), &
332 mpi_double_precision, mpiroot, mpicomm, mpierr)
333 call mpi_bcast(pade_asyliqlw, size(pade_asyliqlw), &
334 mpi_double_precision, mpiroot, mpicomm, mpierr)
335 call mpi_bcast(pade_exticelw, size(pade_exticelw), &
336 mpi_double_precision, mpiroot, mpicomm, mpierr)
337 call mpi_bcast(pade_ssaicelw, size(pade_ssaicelw), &
338 mpi_double_precision, mpiroot, mpicomm, mpierr)
339 call mpi_bcast(pade_asyicelw, size(pade_asyicelw), &
340 mpi_double_precision, mpiroot, mpicomm, mpierr)
341 call mpi_bcast(pade_sizereg_extliqlw, size(pade_sizereg_extliqlw), &
342 mpi_double_precision, mpiroot, mpicomm, mpierr)
343 call mpi_bcast(pade_sizereg_ssaliqlw, size(pade_sizereg_ssaliqlw), &
344 mpi_double_precision, mpiroot, mpicomm, mpierr)
345 call mpi_bcast(pade_sizereg_asyliqlw, size(pade_sizereg_asyliqlw), &
346 mpi_double_precision, mpiroot, mpicomm, mpierr)
347 call mpi_bcast(pade_sizereg_exticelw, size(pade_sizereg_exticelw), &
348 mpi_double_precision, mpiroot, mpicomm, mpierr)
349 call mpi_bcast(pade_sizereg_ssaicelw, size(pade_sizereg_ssaicelw), &
350 mpi_double_precision, mpiroot, mpicomm, mpierr)
351 call mpi_bcast(pade_sizereg_asyicelw, size(pade_sizereg_asyicelw), &
352 mpi_double_precision, mpiroot, mpicomm, mpierr)
353 endif
354#endif
355
356 ! #######################################################################################
357 !
358 ! Initialize RRTMGP DDT's...
359 !
360 ! #######################################################################################
361 if (dogp_cldoptics_lut) then
362 call check_error_msg('lw_cloud_optics_init',lw_cloud_props%load(band_limscldlw, &
363 radliq_lwrlw, radliq_uprlw, radliq_faclw, radice_lwrlw, radice_uprlw, &
364 radice_faclw, lut_extliqlw, lut_ssaliqlw, lut_asyliqlw, lut_exticelw, &
365 lut_ssaicelw, lut_asyicelw))
366 endif
367
368 if (dogp_cldoptics_pade) then
369 call check_error_msg('lw_cloud_optics_init', lw_cloud_props%load(band_limscldlw, &
370 pade_extliqlw, pade_ssaliqlw, pade_asyliqlw, pade_exticelw, pade_ssaicelw, &
371 pade_asyicelw, pade_sizereg_extliqlw, pade_sizereg_ssaliqlw, &
372 pade_sizereg_asyliqlw, pade_sizereg_exticelw, pade_sizereg_ssaicelw, &
373 pade_sizereg_asyicelw))
374 endif
375
376 call check_error_msg('lw_cloud_optics_init',lw_cloud_props%set_ice_roughness(nrghice))
377
378 end subroutine rrtmgp_lw_cloud_optics_init
379end module rrtmgp_lw_cloud_optics
This module contains tools for radiation.
This module contains two routines: The first initializes data and functions needed to compute the lon...
This module contains two routines: One to initialize the k-distribution data and functions needed to ...