CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
rrtmgp_lw_gas_optics.F90
1
8 use machine, only: kind_phys
9 use mo_rte_kind, only: wl
10 use mo_gas_optics_rrtmgp, only: ty_gas_optics_rrtmgp
11 use mo_gas_concentrations, only: ty_gas_concs
12 use radiation_tools, only: check_error_msg
13 use netcdf
14#ifdef MPI
15 use mpi_f08
16#endif
17
18 implicit none
19
20 type(ty_gas_optics_rrtmgp) :: lw_gas_props
21 integer :: &
22 ntempslw, npresslw, ngptslw, nabsorberslw, nextrabsorberslw, nminorabsorberslw,&
23 nmixingfracslw, nlayerslw, nbndslw, npairslw, ninternalsourcetempslw, &
24 nminor_absorber_intervals_lowerlw, nminor_absorber_intervals_upperlw, &
25 ncontributors_lowerlw, ncontributors_upperlw, nfit_coeffslw
26 integer, dimension(:), allocatable :: &
27 kminor_start_lowerlw, & !< Starting index in the [1, nContributors] vector for a contributor
28
29 kminor_start_upperlw
31 integer, dimension(:,:), allocatable :: &
32 band2gptlw, & !< Beginning and ending gpoint for each band
33 minor_limits_gpt_lowerlw, & !< Beginning and ending gpoint for each minor interval in lower atmosphere
34 minor_limits_gpt_upperlw
35 integer, dimension(:,:,:), allocatable :: &
36 key_specieslw
37 real(kind_phys) :: &
38 press_ref_troplw, & !< Reference pressure separating the lower and upper atmosphere [Pa]
39 temp_ref_plw, & !< Standard spectroscopic reference pressure [Pa]
40 temp_ref_tlw
41 real(kind_phys), dimension(:), allocatable :: &
42 press_reflw, & !< Pressures for reference atmosphere; press_ref(# reference layers) [Pa]
43 temp_reflw
44 real(kind_phys), dimension(:,:), allocatable :: &
45 band_limslw, & !< Beginning and ending wavenumber [cm -1] for each band
46 totplnklw, & !< Integrated Planck function by band
47 optimal_angle_fitlw
48 real(kind_phys), dimension(:,:,:), allocatable :: &
49 vmr_reflw, & !< volume mixing ratios for reference atmospherer
50 kminor_lowerlw, & !< (transformed from [nTemp x nEta x nGpt x nAbsorbers] array to
51
52 kminor_upperlw, &
54 rayl_lowerlw, &
55 rayl_upperlw
56 real(kind_phys), dimension(:,:,:,:), allocatable :: &
57 kmajorlw, & !< Stored absorption coefficients due to major absorbing gases
58 planck_fraclw
59 character(len=32), dimension(:), allocatable :: &
60 gas_nameslw, & !< Names of absorbing gases
61 gas_minorlw, & !< Name of absorbing minor gas
62 identifier_minorlw, & !< Unique string identifying minor gas
63 minor_gases_lowerlw, & !< Names of minor absorbing gases in lower atmosphere
64 minor_gases_upperlw, & !< Names of minor absorbing gases in upper atmosphere
65 scaling_gas_lowerlw, & !< Absorption also depends on the concentration of this gas
66 scaling_gas_upperlw
67 logical(wl), dimension(:), allocatable :: &
68 minor_scales_with_density_lowerlw, & !< Density scaling is applied to minor absorption coefficients
69 minor_scales_with_density_upperlw, & !< Density scaling is applied to minor absorption coefficients
70 scale_by_complement_lowerlw, & !< Absorption is scaled by concentration of scaling_gas (F) or its complement (T)
71 scale_by_complement_upperlw
72
73contains
74
76 subroutine rrtmgp_lw_gas_optics_init(rrtmgp_root_dir, rrtmgp_lw_file_gas, &
77 active_gases_array, mpicomm, mpirank, mpiroot, errmsg, errflg)
78
79 ! Inputs
80 character(len=128),intent(in) :: &
81 rrtmgp_root_dir, & !< RTE-RRTMGP root directory
82 rrtmgp_lw_file_gas
83 character(len=*), dimension(:), intent(in) :: &
84 active_gases_array
85 type(mpi_comm),intent(in) :: &
86 mpicomm
87 integer,intent(in) :: &
88 mpirank, & !< Current MPI rank
89 mpiroot
90
91 ! Outputs
92 character(len=*), intent(out) :: &
93 errmsg
94 integer, intent(out) :: &
95 errflg
96
97 ! Local variables
98 integer :: ncid, dimID, varID, status, ii, mpierr, iChar
99 integer,dimension(:),allocatable :: temp1, temp2, temp3, temp4
100 character(len=264) :: lw_gas_props_file
101 type(ty_gas_concs) :: gas_concs ! RRTMGP DDT: trace gas concentrations (vmr)
102
103 ! Initialize
104 errmsg = ''
105 errflg = 0
106
107 ! Filenames are set in the physics_nml
108 lw_gas_props_file = trim(rrtmgp_root_dir)//trim(rrtmgp_lw_file_gas)
109
110 ! #######################################################################################
111 !
112 ! Read dimensions for k-distribution fields...
113 ! (ONLY master processor(0), if MPI enabled)
114 !
115 ! #######################################################################################
116#ifdef MPI
117 if (mpirank .eq. mpiroot) then
118#endif
119 write (*,*) 'Reading RRTMGP longwave k-distribution metadata ... '
120
121 ! Open file
122 status = nf90_open(trim(lw_gas_props_file), nf90_nowrite, ncid)
123
124 ! Read dimensions for k-distribution fields
125 status = nf90_inq_dimid(ncid, 'temperature', dimid)
126 status = nf90_inquire_dimension(ncid, dimid, len = ntempslw)
127 status = nf90_inq_dimid(ncid, 'pressure', dimid)
128 status = nf90_inquire_dimension(ncid, dimid, len = npresslw)
129 status = nf90_inq_dimid(ncid, 'absorber', dimid)
130 status = nf90_inquire_dimension(ncid, dimid, len = nabsorberslw)
131 status = nf90_inq_dimid(ncid, 'minor_absorber', dimid)
132 status = nf90_inquire_dimension(ncid, dimid, len = nminorabsorberslw)
133 status = nf90_inq_dimid(ncid, 'absorber_ext', dimid)
134 status = nf90_inquire_dimension(ncid, dimid, len = nextrabsorberslw)
135 status = nf90_inq_dimid(ncid, 'mixing_fraction', dimid)
136 status = nf90_inquire_dimension(ncid, dimid, len = nmixingfracslw)
137 status = nf90_inq_dimid(ncid, 'atmos_layer', dimid)
138 status = nf90_inquire_dimension(ncid, dimid, len = nlayerslw)
139 status = nf90_inq_dimid(ncid, 'bnd', dimid)
140 status = nf90_inquire_dimension(ncid, dimid, len = nbndslw)
141 status = nf90_inq_dimid(ncid, 'gpt', dimid)
142 status = nf90_inquire_dimension(ncid, dimid, len = ngptslw)
143 status = nf90_inq_dimid(ncid, 'pair', dimid)
144 status = nf90_inquire_dimension(ncid, dimid, len = npairslw)
145 status = nf90_inq_dimid(ncid, 'contributors_lower', dimid)
146 status = nf90_inquire_dimension(ncid, dimid, len = ncontributors_lowerlw)
147 status = nf90_inq_dimid(ncid, 'contributors_upper', dimid)
148 status = nf90_inquire_dimension(ncid, dimid, len = ncontributors_upperlw)
149 status = nf90_inq_dimid(ncid, 'fit_coeffs', dimid)
150 status = nf90_inquire_dimension(ncid, dimid, len = nfit_coeffslw)
151 status = nf90_inq_dimid(ncid, 'minor_absorber_intervals_lower', dimid)
152 status = nf90_inquire_dimension(ncid, dimid, len = nminor_absorber_intervals_lowerlw)
153 status = nf90_inq_dimid(ncid, 'minor_absorber_intervals_upper', dimid)
154 status = nf90_inquire_dimension(ncid, dimid, len = nminor_absorber_intervals_upperlw)
155 status = nf90_inq_dimid(ncid, 'temperature_Planck', dimid)
156 status = nf90_inquire_dimension(ncid, dimid, len = ninternalsourcetempslw)
157#ifdef MPI
158 endif ! On master processor
159
160 ! Other processors waiting...
161 call mpi_barrier(mpicomm, mpierr)
162
163 ! #######################################################################################
164 !
165 ! Broadcast dimensions...
166 ! (ALL processors)
167 !
168 ! #######################################################################################
169 call mpi_bcast(ntempslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
170 call mpi_bcast(npresslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
171 call mpi_bcast(ngptslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
172 call mpi_bcast(nabsorberslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
173 call mpi_bcast(nextrabsorberslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
174 call mpi_bcast(nminorabsorberslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
175 call mpi_bcast(nmixingfracslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
176 call mpi_bcast(nlayerslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
177 call mpi_bcast(nbndslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
178 call mpi_bcast(npairslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
179 call mpi_bcast(ninternalsourcetempslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
180 call mpi_bcast(nminor_absorber_intervals_lowerlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
181 call mpi_bcast(nminor_absorber_intervals_upperlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
182 call mpi_bcast(ncontributors_lowerlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
183 call mpi_bcast(ncontributors_upperlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
184 call mpi_bcast(nfit_coeffslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
185#endif
186
187 ! Allocate space for arrays
188 if (.not. allocated(gas_nameslw)) &
189 allocate(gas_nameslw(nabsorberslw))
190 if (.not. allocated(scaling_gas_lowerlw)) &
191 allocate(scaling_gas_lowerlw(nminor_absorber_intervals_lowerlw))
192 if (.not. allocated(scaling_gas_upperlw)) &
193 allocate(scaling_gas_upperlw(nminor_absorber_intervals_upperlw))
194 if (.not. allocated(gas_minorlw)) &
195 allocate(gas_minorlw(nminorabsorberslw))
196 if (.not. allocated(identifier_minorlw)) &
197 allocate(identifier_minorlw(nminorabsorberslw))
198 if (.not. allocated(minor_gases_lowerlw)) &
199 allocate(minor_gases_lowerlw(nminor_absorber_intervals_lowerlw))
200 if (.not. allocated(minor_gases_upperlw)) &
201 allocate(minor_gases_upperlw(nminor_absorber_intervals_upperlw))
202 if (.not. allocated(minor_limits_gpt_lowerlw)) &
203 allocate(minor_limits_gpt_lowerlw(npairslw, nminor_absorber_intervals_lowerlw))
204 if (.not. allocated(minor_limits_gpt_upperlw)) &
205 allocate(minor_limits_gpt_upperlw(npairslw, nminor_absorber_intervals_upperlw))
206 if (.not. allocated(band2gptlw)) &
207 allocate(band2gptlw(2, nbndslw))
208 if (.not. allocated(key_specieslw)) &
209 allocate(key_specieslw(2, nlayerslw, nbndslw))
210 if (.not. allocated(band_limslw)) &
211 allocate(band_limslw(2, nbndslw))
212 if (.not. allocated(press_reflw)) &
213 allocate(press_reflw(npresslw))
214 if (.not. allocated(temp_reflw)) &
215 allocate(temp_reflw(ntempslw))
216 if (.not. allocated(vmr_reflw)) &
217 allocate(vmr_reflw(nlayerslw, nextrabsorberslw, ntempslw))
218 if (.not. allocated(kminor_lowerlw)) &
219 allocate(kminor_lowerlw(ncontributors_lowerlw, nmixingfracslw, ntempslw))
220 if (.not. allocated(kmajorlw)) &
221 allocate(kmajorlw(ngptslw, nmixingfracslw, npresslw+1, ntempslw))
222 if (.not. allocated(kminor_start_lowerlw)) &
223 allocate(kminor_start_lowerlw(nminor_absorber_intervals_lowerlw))
224 if (.not. allocated(kminor_upperlw)) &
225 allocate(kminor_upperlw(ncontributors_upperlw, nmixingfracslw, ntempslw))
226 if (.not. allocated(kminor_start_upperlw)) &
227 allocate(kminor_start_upperlw(nminor_absorber_intervals_upperlw))
228 if (.not. allocated(optimal_angle_fitlw)) &
229 allocate(optimal_angle_fitlw(nfit_coeffslw, nbndslw))
230 if (.not. allocated(minor_scales_with_density_lowerlw)) &
231 allocate(minor_scales_with_density_lowerlw(nminor_absorber_intervals_lowerlw))
232 if (.not. allocated(minor_scales_with_density_upperlw)) &
233 allocate(minor_scales_with_density_upperlw(nminor_absorber_intervals_upperlw))
234 if (.not. allocated(scale_by_complement_lowerlw)) &
235 allocate(scale_by_complement_lowerlw(nminor_absorber_intervals_lowerlw))
236 if (.not. allocated(scale_by_complement_upperlw)) &
237 allocate(scale_by_complement_upperlw(nminor_absorber_intervals_upperlw))
238 if (.not. allocated(temp1)) &
239 allocate(temp1(nminor_absorber_intervals_lowerlw))
240 if (.not. allocated(temp2)) &
241 allocate(temp2(nminor_absorber_intervals_upperlw))
242 if (.not. allocated(temp3)) &
243 allocate(temp3(nminor_absorber_intervals_lowerlw))
244 if (.not. allocated(temp4)) &
245 allocate(temp4(nminor_absorber_intervals_upperlw))
246 if (.not. allocated(totplnklw)) &
247 allocate(totplnklw(ninternalsourcetempslw, nbndslw))
248 if (.not. allocated(planck_fraclw)) &
249 allocate(planck_fraclw(ngptslw, nmixingfracslw, npresslw+1, ntempslw))
250
251 ! #######################################################################################
252 !
253 ! Read in data ...
254 ! (ONLY master processor(0), if MPI enabled)
255 !
256 ! #######################################################################################
257#ifdef MPI
258 if (mpirank .eq. mpiroot) then
259#endif
260 write (*,*) 'Reading RRTMGP longwave k-distribution data ... '
261 status = nf90_inq_varid(ncid, 'gas_names', varid)
262 status = nf90_get_var( ncid, varid, gas_nameslw)
263 status = nf90_inq_varid(ncid, 'scaling_gas_lower', varid)
264 status = nf90_get_var( ncid, varid, scaling_gas_lowerlw)
265 status = nf90_inq_varid(ncid, 'scaling_gas_upper', varid)
266 status = nf90_get_var( ncid, varid, scaling_gas_upperlw)
267 status = nf90_inq_varid(ncid, 'gas_minor', varid)
268 status = nf90_get_var( ncid, varid, gas_minorlw)
269 status = nf90_inq_varid(ncid, 'identifier_minor', varid)
270 status = nf90_get_var( ncid, varid, identifier_minorlw)
271 status = nf90_inq_varid(ncid, 'minor_gases_lower', varid)
272 status = nf90_get_var( ncid, varid, minor_gases_lowerlw)
273 status = nf90_inq_varid(ncid, 'minor_gases_upper', varid)
274 status = nf90_get_var( ncid, varid, minor_gases_upperlw)
275 status = nf90_inq_varid(ncid, 'minor_limits_gpt_lower', varid)
276 status = nf90_get_var( ncid, varid, minor_limits_gpt_lowerlw)
277 status = nf90_inq_varid(ncid, 'minor_limits_gpt_upper', varid)
278 status = nf90_get_var( ncid, varid, minor_limits_gpt_upperlw)
279 status = nf90_inq_varid(ncid, 'bnd_limits_gpt', varid)
280 status = nf90_get_var( ncid, varid, band2gptlw)
281 status = nf90_inq_varid(ncid, 'key_species', varid)
282 status = nf90_get_var( ncid, varid, key_specieslw)
283 status = nf90_inq_varid(ncid, 'bnd_limits_wavenumber', varid)
284 status = nf90_get_var( ncid, varid, band_limslw)
285 status = nf90_inq_varid(ncid, 'press_ref', varid)
286 status = nf90_get_var( ncid, varid, press_reflw)
287 status = nf90_inq_varid(ncid, 'temp_ref', varid)
288 status = nf90_get_var( ncid, varid, temp_reflw)
289 status = nf90_inq_varid(ncid, 'absorption_coefficient_ref_P', varid)
290 status = nf90_get_var( ncid, varid, temp_ref_plw)
291 status = nf90_inq_varid(ncid, 'absorption_coefficient_ref_T', varid)
292 status = nf90_get_var( ncid, varid, temp_ref_tlw)
293 status = nf90_inq_varid(ncid, 'press_ref_trop', varid)
294 status = nf90_get_var( ncid, varid, press_ref_troplw)
295 status = nf90_inq_varid(ncid, 'kminor_lower', varid)
296 status = nf90_get_var( ncid, varid, kminor_lowerlw)
297 status = nf90_inq_varid(ncid, 'kminor_upper', varid)
298 status = nf90_get_var( ncid, varid, kminor_upperlw)
299 status = nf90_inq_varid(ncid, 'vmr_ref', varid)
300 status = nf90_get_var( ncid, varid, vmr_reflw)
301 status = nf90_inq_varid(ncid, 'optimal_angle_fit',varid)
302 status = nf90_get_var( ncid, varid, optimal_angle_fitlw)
303 status = nf90_inq_varid(ncid, 'kmajor', varid)
304 status = nf90_get_var( ncid, varid, kmajorlw)
305 status = nf90_inq_varid(ncid, 'kminor_start_lower', varid)
306 status = nf90_get_var( ncid, varid, kminor_start_lowerlw)
307 status = nf90_inq_varid(ncid, 'kminor_start_upper', varid)
308 status = nf90_get_var( ncid, varid, kminor_start_upperlw)
309 status = nf90_inq_varid(ncid, 'totplnk', varid)
310 status = nf90_get_var( ncid, varid, totplnklw)
311 status = nf90_inq_varid(ncid, 'plank_fraction', varid)
312 status = nf90_get_var( ncid, varid, planck_fraclw)
313
314 ! Logical fields are read in as integers and then converted to logicals.
315 status = nf90_inq_varid(ncid,'minor_scales_with_density_lower', varid)
316 status = nf90_get_var( ncid, varid,temp1)
317 status = nf90_inq_varid(ncid,'minor_scales_with_density_upper', varid)
318 status = nf90_get_var( ncid, varid,temp2)
319 status = nf90_inq_varid(ncid,'scale_by_complement_lower', varid)
320 status = nf90_get_var( ncid, varid,temp3)
321 status = nf90_inq_varid(ncid,'scale_by_complement_upper', varid)
322 status = nf90_get_var( ncid, varid,temp4)
323 status = nf90_close(ncid)
324
325 do ii=1,nminor_absorber_intervals_lowerlw
326 if (temp1(ii) .eq. 0) minor_scales_with_density_lowerlw(ii) = .false.
327 if (temp1(ii) .eq. 1) minor_scales_with_density_lowerlw(ii) = .true.
328 if (temp3(ii) .eq. 0) scale_by_complement_lowerlw(ii) = .false.
329 if (temp3(ii) .eq. 1) scale_by_complement_lowerlw(ii) = .true.
330 enddo
331 do ii=1,nminor_absorber_intervals_upperlw
332 if (temp2(ii) .eq. 0) minor_scales_with_density_upperlw(ii) = .false.
333 if (temp2(ii) .eq. 1) minor_scales_with_density_upperlw(ii) = .true.
334 if (temp4(ii) .eq. 0) scale_by_complement_upperlw(ii) = .false.
335 if (temp4(ii) .eq. 1) scale_by_complement_upperlw(ii) = .true.
336 enddo
337#ifdef MPI
338 endif ! Master process
339
340 ! Other processors waiting...
341 call mpi_barrier(mpicomm, mpierr)
342
343 ! #######################################################################################
344 !
345 ! Broadcast data...
346 ! (ALL processors)
347 !
348 ! #######################################################################################
349
350 ! Real scalars
351 call mpi_bcast(press_ref_troplw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
352 call mpi_bcast(temp_ref_plw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
353 call mpi_bcast(temp_ref_tlw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
354
355 ! Integer arrays
356 call mpi_bcast(kminor_start_lowerlw, &
357 size(kminor_start_lowerlw), mpi_integer, mpiroot, mpicomm, mpierr)
358 call mpi_bcast(kminor_start_upperlw, &
359 size(kminor_start_upperlw), mpi_integer, mpiroot, mpicomm, mpierr)
360 call mpi_bcast(band2gptlw, &
361 size(band2gptlw), mpi_integer, mpiroot, mpicomm, mpierr)
362 call mpi_bcast(minor_limits_gpt_lowerlw, &
363 size(minor_limits_gpt_lowerlw), mpi_integer, mpiroot, mpicomm, mpierr)
364 call mpi_bcast(minor_limits_gpt_upperlw, &
365 size(minor_limits_gpt_upperlw), mpi_integer, mpiroot, mpicomm, mpierr)
366 call mpi_bcast(key_specieslw, &
367 size(key_specieslw), mpi_integer, mpiroot, mpicomm, mpierr)
368
369 ! Real arrays
370 call mpi_bcast(press_reflw, &
371 size(press_reflw), mpi_double_precision, mpiroot, mpicomm, mpierr)
372 call mpi_bcast(temp_reflw, &
373 size(temp_reflw), mpi_double_precision, mpiroot, mpicomm, mpierr)
374 call mpi_bcast(band_limslw, &
375 size(band_limslw), mpi_double_precision, mpiroot, mpicomm, mpierr)
376 call mpi_bcast(totplnklw, &
377 size(totplnklw), mpi_double_precision, mpiroot, mpicomm, mpierr)
378 call mpi_bcast(optimal_angle_fitlw, &
379 size(optimal_angle_fitlw), mpi_double_precision, mpiroot, mpicomm, mpierr)
380 call mpi_bcast(vmr_reflw, &
381 size(vmr_reflw), mpi_double_precision, mpiroot, mpicomm, mpierr)
382 call mpi_bcast(kminor_lowerlw, &
383 size(kminor_lowerlw), mpi_double_precision, mpiroot, mpicomm, mpierr)
384 call mpi_bcast(kminor_upperlw, &
385 size(kminor_upperlw), mpi_double_precision, mpiroot, mpicomm, mpierr)
386 call mpi_bcast(kmajorlw, &
387 size(kmajorlw), mpi_double_precision, mpiroot, mpicomm, mpierr)
388 call mpi_bcast(planck_fraclw, &
389 size(planck_fraclw), mpi_double_precision, mpiroot, mpicomm, mpierr)
390
391
392 ! Characters
393 do ichar=1,nabsorberslw
394 call mpi_bcast(gas_nameslw(ichar), &
395 len(gas_nameslw(ichar)), mpi_character, mpiroot, mpicomm, mpierr)
396 enddo
397 do ichar=1,nminorabsorberslw
398 call mpi_bcast(gas_minorlw(ichar), &
399 len(gas_minorlw(ichar)), mpi_character, mpiroot, mpicomm, mpierr)
400 call mpi_bcast(identifier_minorlw(ichar), &
401 len(identifier_minorlw(ichar)), mpi_character, mpiroot, mpicomm, mpierr)
402 enddo
403 do ichar=1,nminor_absorber_intervals_lowerlw
404 call mpi_bcast(minor_gases_lowerlw(ichar), &
405 len(minor_gases_lowerlw(ichar)), mpi_character, mpiroot, mpicomm, mpierr)
406 call mpi_bcast(scaling_gas_lowerlw(ichar), &
407 len(scaling_gas_lowerlw(ichar)), mpi_character, mpiroot, mpicomm, mpierr)
408 enddo
409 do ichar=1,nminor_absorber_intervals_upperlw
410 call mpi_bcast(minor_gases_upperlw(ichar), &
411 len(minor_gases_upperlw(ichar)), mpi_character, mpiroot, mpicomm, mpierr)
412 call mpi_bcast(scaling_gas_upperlw(ichar), &
413 len(scaling_gas_upperlw(ichar)), mpi_character, mpiroot, mpicomm, mpierr)
414 enddo
415
416 ! Logicals
417 call mpi_bcast(minor_scales_with_density_lowerlw, &
418 size(minor_scales_with_density_lowerlw), mpi_logical, mpiroot, mpicomm, mpierr)
419 call mpi_bcast(minor_scales_with_density_upperlw, &
420 size(minor_scales_with_density_upperlw), mpi_logical, mpiroot, mpicomm, mpierr)
421 call mpi_bcast(scale_by_complement_lowerlw, &
422 size(scale_by_complement_lowerlw), mpi_logical, mpiroot, mpicomm, mpierr)
423 call mpi_bcast(scale_by_complement_upperlw, &
424 size(scale_by_complement_upperlw), mpi_logical, mpiroot, mpicomm, mpierr)
425
426 call mpi_barrier(mpicomm, mpierr)
427#endif
428
429 ! #######################################################################################
430 !
431 ! Initialize RRTMGP DDT's...
432 !
433 ! #######################################################################################
434 call check_error_msg('rrtmgp_lw_gas_optics_init_gas_concs',gas_concs%init(active_gases_array))
435 call check_error_msg('rrtmgp_lw_gas_optics_init_load',lw_gas_props%load(gas_concs, &
436 gas_nameslw, key_specieslw, band2gptlw, band_limslw, press_reflw, press_ref_troplw,&
437 temp_reflw, temp_ref_plw, temp_ref_tlw, vmr_reflw, kmajorlw, kminor_lowerlw, &
438 kminor_upperlw, gas_minorlw, identifier_minorlw, minor_gases_lowerlw, &
439 minor_gases_upperlw, minor_limits_gpt_lowerlw, minor_limits_gpt_upperlw, &
440 minor_scales_with_density_lowerlw, minor_scales_with_density_upperlw, &
441 scaling_gas_lowerlw, scaling_gas_upperlw, scale_by_complement_lowerlw, &
442 scale_by_complement_upperlw, kminor_start_lowerlw, kminor_start_upperlw, totplnklw,&
443 planck_fraclw, rayl_lowerlw, rayl_upperlw, optimal_angle_fitlw))
444
445 end subroutine rrtmgp_lw_gas_optics_init
446
447end module rrtmgp_lw_gas_optics