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