source: LMDZ6/branches/Ocean_skin/libf/phylmd/ecrad/radiation_config.F90 @ 5407

Last change on this file since 5407 was 4368, checked in by lguez, 2 years ago

Sync latest trunk changes to Ocean_skin

File size: 88.0 KB
Line 
1! radiation_config.F90 - Derived type to configure the radiation scheme
2!
3! (C) Copyright 2014- ECMWF.
4!
5! This software is licensed under the terms of the Apache Licence Version 2.0
6! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
7!
8! In applying this licence, ECMWF does not waive the privileges and immunities
9! granted to it by virtue of its status as an intergovernmental organisation
10! nor does it submit to any jurisdiction.
11!
12! Author:  Robin Hogan
13! Email:   r.j.hogan@ecmwf.int
14!
15! Modifications
16!   2017-07-22  R. Hogan  Added Yi et al. ice optics model
17!   2017-10-23  R. Hogan  Renamed single-character variables
18!   2018-03-15  R. Hogan  Added logicals controlling surface spectral treatment
19!   2018-08-29  R. Hogan  Added monochromatic single-scattering albedo / asymmetry factor
20!   2018-09-03  R. Hogan  Added min_cloud_effective_size
21!   2018-09-04  R. Hogan  Added encroachment_scaling
22!   2018-09-13  R. Hogan  Added IEncroachmentFractal
23!   2019-01-02  R. Hogan  Added Cloudless solvers
24!   2019-01-14  R. Hogan  Added out_of_bounds_[1,2,3]d for checker routines
25!   2019-01-18  R. Hogan  Added albedo weighting
26!   2019-02-03  R. Hogan  Added ability to fix out-of-physical-bounds inputs
27!   2019-02-10  R. Hogan  Renamed "encroachment" to "entrapment"
28!
29! Note: The aim is for ecRad in the IFS to be as similar as possible
30! to the offline version, so if you make any changes to this or any
31! files in this directory, please inform Robin Hogan.
32!
33
34module radiation_config
35
36  use parkind1,                      only : jprb
37
38  use radiation_cloud_optics_data,   only : cloud_optics_type
39  use radiation_aerosol_optics_data, only : aerosol_optics_type
40  use radiation_pdf_sampler,         only : pdf_sampler_type
41  use radiation_cloud_cover,         only : OverlapName, &
42       & IOverlapMaximumRandom, IOverlapExponentialRandom, IOverlapExponential
43
44  implicit none
45  public
46
47  ! Configuration codes: use C-style enumerators to avoid having to
48  ! remember the numbers
49
50  ! Solvers: can be specified for longwave and shortwave
51  ! independently, except for "Homogeneous", which must be the same
52  ! for both
53  enum, bind(c)
54     enumerator ISolverCloudless, ISolverHomogeneous, ISolverMcICA, &
55          &     ISolverSpartacus, ISolverTripleclouds
56  end enum
57  character(len=*), parameter :: SolverName(0:4) = (/ 'Cloudless   ', &
58       &                                              'Homogeneous ', &
59       &                                              'McICA       ', &
60       &                                              'SPARTACUS   ', &
61       &                                              'Tripleclouds' /)
62
63  ! SPARTACUS shortwave solver can treat the reflection of radiation
64  ! back up into different regions in various ways
65  enum, bind(c)
66     enumerator &
67       & IEntrapmentZero, &     ! No entrapment, as Tripleclouds
68       & IEntrapmentEdgeOnly, & ! Only radiation passed through cloud edge is horizontally homogenized
69       & IEntrapmentExplicit, & ! Estimate horiz migration dist, account for fractal clouds
70       & IEntrapmentExplicitNonFractal, & ! As above but ignore fractal nature of clouds
71       & IEntrapmentMaximum ! Complete horizontal homogenization within regions (old SPARTACUS assumption)
72
73  end enum
74 
75  ! Names available in the radiation namelist for variable
76  ! sw_entrapment_name
77  character(len=*), parameter :: EntrapmentName(0:4)   = [ 'Zero       ', &
78       &                                                   'Edge-only  ', &
79       &                                                   'Explicit   ', &
80       &                                                   'Non-fractal', &
81       &                                                   'Maximum    ' ]
82  ! For backwards compatibility, the radiation namelist also supports
83  ! the equivalent variable sw_encroachment_name with the following
84  ! names
85  character(len=*), parameter :: EncroachmentName(0:4) = [ 'Zero    ', &
86       &                                                   'Minimum ', &
87       &                                                   'Fractal ', &
88       &                                                   'Computed', &
89       &                                                   'Maximum ' ]
90
91  ! Two-stream models
92  ! This is not configurable at run-time
93
94  ! Gas models
95  enum, bind(c)
96     enumerator IGasModelMonochromatic, IGasModelIFSRRTMG
97  end enum
98  character(len=*), parameter :: GasModelName(0:1) = (/ 'Monochromatic', &
99       &                                                'RRTMG-IFS    ' /)
100
101  ! Hydrometeor scattering models
102  enum, bind(c)
103     enumerator ILiquidModelMonochromatic, &
104          &     ILiquidModelSOCRATES, ILiquidModelSlingo
105  end enum
106  character(len=*), parameter :: LiquidModelName(0:2) = (/ 'Monochromatic', &
107       &                                                   'SOCRATES     ', &
108       &                                                   'Slingo       ' /)
109
110  enum, bind(c)
111     enumerator IIceModelMonochromatic, IIceModelFu, &
112          &  IIceModelBaran, IIceModelBaran2016, IIceModelBaran2017,   &
113          &  IIceModelYi
114  end enum
115  character(len=*), parameter :: IceModelName(0:5) = (/ 'Monochromatic', &
116       &                                                'Fu-IFS       ', &
117       &                                                'Baran        ', &
118       &                                                'Baran2016    ', &
119       &                                                'Baran2017    ', &
120       &                                                'Yi           ' /)
121
122  ! Cloud PDF distribution shapes
123  enum, bind(c)
124    enumerator IPdfShapeLognormal, IPdfShapeGamma
125  end enum
126  character(len=*), parameter :: PdfShapeName(0:1) = (/ 'Lognormal', &
127       &                                                'Gamma    ' /)
128
129  ! Maximum number of different aerosol types that can be provided
130  integer, parameter :: NMaxAerosolTypes = 256
131
132  ! Maximum number of shortwave albedo and longwave emissivity
133  ! intervals
134  integer, parameter :: NMaxAlbedoIntervals = 256
135
136  ! Length of string buffer for printing config information
137  integer, parameter :: NPrintStringLen = 60
138
139  !---------------------------------------------------------------------
140  ! Derived type containing all the configuration information needed
141  ! to run the radiation scheme.  The intention is that this is fixed
142  ! for a given model run.  The parameters are to list first those
143  ! quantities that can be set directly by the user, for example using a
144  ! namelist, and second those quantities that are computed afterwards
145  ! from the user-supplied numbers, especially the details of the gas
146  ! optics model.
147  type config_type
148    ! USER-CONFIGURABLE PARAMETERS
149
150    ! Override default solar spectrum
151    logical :: use_spectral_solar_scaling = .false.
152
153    ! Directory in which gas, cloud and aerosol data files are to be
154    ! found
155    character(len=511) :: directory_name = '.'
156
157    ! Cloud is deemed to be present in a layer if cloud fraction
158    ! exceeds this value
159    real(jprb) :: cloud_fraction_threshold = 1.0e-6_jprb
160    ! ...and total cloud water mixing ratio exceeds this value
161    real(jprb) :: cloud_mixing_ratio_threshold = 1.0e-9_jprb
162
163    ! Overlap scheme
164    integer :: i_overlap_scheme = IOverlapExponentialRandom
165
166    ! Use the Shonk et al. (2010) "beta" overlap parameter, rather
167    ! than the "alpha" overlap parameter of Hogan and Illingworth
168    ! (2000)?
169    logical :: use_beta_overlap = .false.
170
171    ! Shape of sub-grid cloud water PDF
172    integer :: i_cloud_pdf_shape = IPdfShapeGamma
173
174    ! The ratio of the overlap decorrelation length for cloud
175    ! inhomogeneities to the overlap decorrelation length for cloud
176    ! boundaries.  Observations suggest this has a value of 0.5
177    ! (e.g. from the decorrelation lengths of Hogan and Illingworth
178    ! 2003 and Hogan and Illingworth 2000).
179    real(jprb) :: cloud_inhom_decorr_scaling = 0.5_jprb
180
181    ! Factor controlling how much of the cloud edge length interfaces
182    ! directly between the clear-sky region (region a) and the
183    ! optically thick cloudy region (region c).  If Lxy is the length
184    ! of the interfaces between regions x and y, and Lab and Lbc have
185    ! been computed already, then
186    !   Lac=clear_to_thick_fraction*min(Lab,Lbc).
187    real(jprb) :: clear_to_thick_fraction = 0.0_jprb
188
189    ! Factor allowing lateral transport when the sun is close to
190    ! overhead; consider atand(overhead_sun_factor) to be the number
191    ! of degrees that the sun angle is perturbed from zenith for the
192    ! purposes of computing lateral transport.  A value of up to 0.1
193    ! seems to be necessary to account for the fact that some forward
194    ! scattered radiation is treated as unscattered by delta-Eddington
195    ! scaling; therefore it ought to have the chance to escape.
196    real(jprb) :: overhead_sun_factor = 0.0_jprb
197
198    ! Minimum gas optical depth in a single layer at any wavelength,
199    ! for stability
200    real(jprb) :: min_gas_od_lw = 1.0e-15_jprb
201    real(jprb) :: min_gas_od_sw = 0.0_jprb
202
203    ! Maximum gas optical depth in a layer before that g-point will
204    ! not be considered for 3D treatment: a limit is required to avoid
205    ! expensive computation of matrix exponentials on matrices with
206    ! large elements
207    real(jprb) :: max_gas_od_3d = 8.0_jprb
208
209    ! Maximum total optical depth of a cloudy region for stability:
210    ! optical depth will be capped at this value in the SPARTACUS
211    ! solvers
212    real(jprb) :: max_cloud_od = 16.0_jprb
213
214    ! How much longwave scattering is included?
215    logical :: do_lw_cloud_scattering = .true.
216    logical :: do_lw_aerosol_scattering = .true.
217
218    ! Number of regions used to describe clouds and clear skies. A
219    ! value of 2 means one clear and one cloudy region, so clouds are
220    ! horizontally homogeneous, while a value of 3 means two cloudy
221    ! regions with different optical depth, thereby representing
222    ! inhomogeneity via the Shonk & Hogan (2008) "Tripleclouds"
223    ! method.
224    integer :: nregions = 3
225
226    ! Code specifying the solver to be used: use the enumerations
227    ! defined above
228    integer :: i_solver_sw = ISolverMcICA
229    integer :: i_solver_lw = ISolverMcICA
230
231    ! Do shortwave delta-Eddington scaling on the cloud-aerosol-gas
232    ! mixture (as in the original IFS scheme), rather than the more
233    ! correct approach of separately scaling the cloud and aerosol
234    ! scattering properties before merging with gases.  Note that
235    ! .true. is not compatible with the SPARTACUS solver.
236    logical :: do_sw_delta_scaling_with_gases = .false.
237
238    ! Codes describing the gas and cloud scattering models to use, the
239    ! latter of which is currently not used
240    integer :: i_gas_model = IGasModelIFSRRTMG
241    !     integer :: i_cloud_model
242
243    ! Optics if i_gas_model==IGasModelMonochromatic.
244    ! The wavelength to use for the Planck function in metres. If this
245    ! is positive then the output longwave fluxes will be in units of
246    ! W m-2 um-1.  If this is zero or negative (the default) then
247    ! sigma*T^4 will be used and the output longwave fluxes will be in
248    ! W m-2.
249    real(jprb) :: mono_lw_wavelength = -1.0_jprb
250    ! Total zenith optical depth of the atmosphere in the longwave and
251    ! shortwave, distributed vertically according to the pressure.
252    ! Default is zero.
253    real(jprb) :: mono_lw_total_od = 0.0_jprb
254    real(jprb) :: mono_sw_total_od = 0.0_jprb
255    ! Single-scattering albedo and asymmetry factor: values typical
256    ! for liquid clouds with effective radius of 10 microns, at (SW)
257    ! 0.55 micron wavelength and (LW) 10.7 microns wavelength
258    real(jprb) :: mono_sw_single_scattering_albedo = 0.999999_jprb
259    real(jprb) :: mono_sw_asymmetry_factor = 0.86_jprb
260    real(jprb) :: mono_lw_single_scattering_albedo = 0.538_jprb
261    real(jprb) :: mono_lw_asymmetry_factor = 0.925_jprb
262
263    ! Codes describing particle scattering models
264    integer :: i_liq_model = ILiquidModelSOCRATES
265    integer :: i_ice_model = IIceModelBaran
266   
267    ! The mapping from albedo/emissivity intervals to SW/LW bands can
268    ! either be done by finding the interval containing the central
269    ! wavenumber of the band (nearest neighbour), or by a weighting
270    ! according to the spectral overlap of each interval with each
271    ! band
272    logical :: do_nearest_spectral_sw_albedo = .true.
273    logical :: do_nearest_spectral_lw_emiss  = .true.
274
275    ! User-defined monotonically increasing wavelength bounds (m)
276    ! between input surface albedo/emissivity intervals. Implicitly
277    ! the first interval starts at zero and the last ends at infinity.
278    real(jprb) :: sw_albedo_wavelength_bound(NMaxAlbedoIntervals-1) = -1.0_jprb
279    real(jprb) :: lw_emiss_wavelength_bound( NMaxAlbedoIntervals-1)  = -1.0_jprb
280
281    ! The index to the surface albedo/emissivity intervals for each of
282    ! the wavelength bounds specified in sw_albedo_wavelength_bound
283    ! and lw_emiss_wavelength_bound
284    integer :: i_sw_albedo_index(NMaxAlbedoIntervals) = 0
285    integer :: i_lw_emiss_index (NMaxAlbedoIntervals)  = 0
286
287    ! Do we compute longwave and/or shortwave radiation?
288    logical :: do_lw = .true.
289    logical :: do_sw = .true.
290
291    ! Do we compute clear-sky fluxes and/or solar direct fluxes?
292    logical :: do_clear = .true.
293    logical :: do_sw_direct = .true.
294
295    ! Do we include 3D effects?
296    logical :: do_3d_effects = .true.
297   
298    ! To what extent do we include "entrapment" effects in the
299    ! SPARTACUS solver? This essentially means that in a situation
300    ! like this
301    !
302    ! 000111
303    ! 222222
304    !
305    ! Radiation downwelling from region 1 may be reflected back into
306    ! region 0 due to some degree of homogenization of the radiation
307    ! in region 2.  Hogan and Shonk (2013) referred to this as
308    ! "anomalous horizontal transport" for a 1D model, although for 3D
309    ! calculations it is desirable to include at least some of it. The
310    ! options are described by the IEntrapment* parameters above.
311    integer :: i_3d_sw_entrapment = IEntrapmentExplicit
312
313    ! In the longwave, the equivalent process it either "on" (like
314    ! maximum entrapment) or "off" (like zero entrapment):
315    logical :: do_3d_lw_multilayer_effects = .false.
316
317    ! Do we account for the effective emissivity of the side of
318    ! clouds?
319    logical :: do_lw_side_emissivity = .true.
320
321    ! The 3D transfer rate "X" is such that if transport out of a
322    ! region was the only process occurring then by the base of a
323    ! layer only exp(-X) of the original flux would remain in that
324    ! region. The transfer rate computed geometrically can be very
325    ! high for the clear-sky regions in layers with high cloud
326    ! fraction.  For stability reasons it is necessary to provide a
327    ! maximum possible 3D transfer rate.
328    real(jprb) :: max_3d_transfer_rate = 10.0_jprb
329
330    ! It has also sometimes been found necessary to set a minimum
331    ! cloud effective size for stability (metres)
332    real(jprb) :: min_cloud_effective_size = 100.0_jprb
333
334    ! Given a horizontal migration distance, there is still
335    ! uncertainty about how much entrapment occurs associated with how
336    ! one assumes cloud boundaries line up in adjacent layers. This
337    ! factor can be varied between 0.0 (the boundaries line up to the
338    ! greatest extent possible given the overlap parameter) and 1.0
339    ! (the boundaries line up to the minimum extent possible). In the
340    ! Hogan et al. entrapment paper it is referred to as the overhang
341    ! factor zeta, and a value of 0 matches the Monte Carlo
342    ! calculations best.
343    real(jprb) :: overhang_factor = 0.0_jprb
344
345    ! By default, the Meador & Weaver (1980) expressions are used
346    ! instead of the matrix exponential whenever 3D effects can be
347    ! neglected (e.g. cloud-free layers or clouds with infinitely
348    ! large effective cloud size), but setting the following to true
349    ! uses the matrix exponential everywhere, enabling the two
350    ! methods to be compared. Note that Meador & Weaver will still be
351    ! used for very optically thick g points where the matrix
352    ! exponential can produce incorrect results.
353    logical :: use_expm_everywhere = .false.
354
355    ! Aerosol descriptors: aerosol_type_mapping must be of length
356    ! n_aerosol_types, and contains 0 if that type is to be ignored,
357    ! positive numbers to map on to the indices of hydrophobic
358    ! aerosols in the aerosol optics configuration file, and negative
359    ! numbers to map on to (the negative of) the indices of
360    ! hydrophilic aerosols in the configuration file.
361    logical :: use_aerosols = .false.
362    integer :: n_aerosol_types = 0
363    integer :: i_aerosol_type_map(NMaxAerosolTypes)
364
365    ! Save the gas and cloud optical properties for each g point in
366    ! "radiative_properties.nc"?
367    logical :: do_save_radiative_properties = .false.
368
369    ! Save the flux profiles in each band?
370    logical :: do_save_spectral_flux = .false.
371
372    ! Save the surface downwelling shortwave fluxes in each band?
373    logical :: do_surface_sw_spectral_flux = .true.
374
375    ! Compute the longwave derivatives needed to apply the approximate
376    ! radiation updates of Hogan and Bozzo (2015)
377    logical :: do_lw_derivatives = .false.
378
379    ! Save the flux profiles in each g-point (overrides
380    ! do_save_spectral_flux if TRUE)?
381    logical :: do_save_gpoint_flux = .false.
382
383    ! In the IFS environment, setting up RRTM has already been done
384    ! so not needed to do it again
385    logical :: do_setup_ifsrrtm = .true.
386
387    ! In the IFS environment the old scheme has a bug in the Fu
388    ! longwave ice optics whereby the single scattering albedo is one
389    ! minus what it should be.  Unfortunately fixing it makes
390    ! forecasts worse. Setting the following to true reproduces the
391    ! bug.
392    logical :: do_fu_lw_ice_optics_bug = .false.
393
394    ! Control verbosity: 0=none (no output to standard output; write
395    ! to standard error only if an error occurs), 1=warning, 2=info,
396    ! 3=progress, 4=detailed, 5=debug.  Separate settings for the
397    ! setup of the scheme and the execution of it.
398    integer :: iverbosesetup = 2
399    integer :: iverbose = 1
400
401    ! Are we doing radiative transfer in complex surface canopies
402    ! (streets/vegetation), in which case tailored downward fluxes are
403    ! needed at the top of the canopy?
404    logical :: do_canopy_fluxes_sw = .false.
405    logical :: do_canopy_fluxes_lw = .false.
406    ! If so, do we use the full spectrum as in the atmosphere, or just
407    ! the reduced spectrum in which the shortwave albedo and longwave
408    ! emissivity are provided?
409    logical :: use_canopy_full_spectrum_sw = .false.
410    logical :: use_canopy_full_spectrum_lw = .false.
411    ! Do we treat gas radiative transfer in streets/vegetation?
412    logical :: do_canopy_gases_sw = .false.
413    logical :: do_canopy_gases_lw = .false.
414
415    ! Optics file names for overriding the ones generated from the
416    ! other options. If these remain empty then the generated names
417    ! will be used (see the "consolidate_config" routine below). If
418    ! the user assigns one of these and it starts with a '/' character
419    ! then that will be used instead. If the user assigns one and it
420    ! doesn't start with a '/' character then it will be prepended by
421    ! the contents of directory_name.
422    character(len=511) :: ice_optics_override_file_name = ''
423    character(len=511) :: liq_optics_override_file_name = ''
424    character(len=511) :: aerosol_optics_override_file_name = ''
425
426    ! Optionally override the look-up table file for the cloud-water
427    ! PDF used by the McICA solver
428    character(len=511) :: cloud_pdf_override_file_name = ''
429
430    ! Has "consolidate" been called? 
431    logical :: is_consolidated = .false.
432
433    ! COMPUTED PARAMETERS
434    ! Users of this library should not edit these parameters directly;
435    ! they are set by the "consolidate" routine
436
437    ! Wavenumber range for each band, in cm-1, which will be allocated
438    ! to be of length n_bands_sw or n_bands_lw
439    real(jprb), allocatable, dimension(:) :: wavenumber1_sw
440    real(jprb), allocatable, dimension(:) :: wavenumber2_sw
441    real(jprb), allocatable, dimension(:) :: wavenumber1_lw
442    real(jprb), allocatable, dimension(:) :: wavenumber2_lw
443
444    ! If the nearest surface albedo/emissivity interval is to be used
445    ! for each SW/LW band then the following arrays will be allocated
446    ! to the length of the number of bands and contain the index to
447    ! the relevant interval
448    integer, allocatable, dimension(:) :: i_albedo_from_band_sw
449    integer, allocatable, dimension(:) :: i_emiss_from_band_lw
450
451    ! ...alternatively, this matrix dimensioned
452    ! (n_albedo_intervals,n_bands_sw) providing the weights needed for
453    ! computing the albedo in each ecRad band from the albedo in each
454    ! native albedo band - see radiation_single_level.F90
455    real(jprb), allocatable, dimension(:,:) :: sw_albedo_weights
456    ! ...and similarly in the longwave, dimensioned
457    ! (n_emiss_intervals,n_bands_lw)
458    real(jprb), allocatable, dimension(:,:) :: lw_emiss_weights
459
460    ! Arrays of length the number of g-points that convert from
461    ! g-point to the band index
462    integer, allocatable, dimension(:) :: i_band_from_g_lw
463    integer, allocatable, dimension(:) :: i_band_from_g_sw
464
465    ! We allow for the possibility for g-points to be ordered in terms
466    ! of likely absorption (weakest to strongest) across the shortwave
467    ! or longwave spectrum, in order that in SPARTACUS we select only
468    ! the first n g-points that will not have too large an absorption,
469    ! and therefore matrix exponentials that are both finite and not
470    ! too expensive to compute.  The following two arrays map the
471    ! reordered g-points to the original ones.
472    integer, allocatable, dimension(:) :: i_g_from_reordered_g_lw
473    integer, allocatable, dimension(:) :: i_g_from_reordered_g_sw
474
475    ! The following map the reordered g-points to the bands
476    integer, allocatable, dimension(:) :: i_band_from_reordered_g_lw
477    integer, allocatable, dimension(:) :: i_band_from_reordered_g_sw
478
479    ! The following map the reordered g-points to the spectral
480    ! information being saved: if do_save_gpoint_flux==TRUE then this
481    ! will map on to the original g points, but if only
482    ! do_save_spectral_flux==TRUE then this will map on to the bands
483    integer, pointer, dimension(:) :: i_spec_from_reordered_g_lw
484    integer, pointer, dimension(:) :: i_spec_from_reordered_g_sw
485
486    ! Number of spectral intervals used for the canopy radiative
487    ! transfer calculation; they are either equal to
488    ! n_albedo_intervals/n_emiss_intervals or n_g_sw/n_g_lw
489    integer :: n_canopy_bands_sw = 1
490    integer :: n_canopy_bands_lw = 1
491
492    ! Data structure containing cloud scattering data
493    type(cloud_optics_type)      :: cloud_optics
494
495    ! Data structure containing aerosol scattering data
496    type(aerosol_optics_type)    :: aerosol_optics
497
498    ! Object for sampling from a gamma or lognormal distribution
499    type(pdf_sampler_type)       :: pdf_sampler
500
501    ! Optics file names
502    character(len=511) :: ice_optics_file_name, &
503         &                liq_optics_file_name, &
504         &                aerosol_optics_file_name
505   
506    ! McICA PDF look-up table file name
507    character(len=511) :: cloud_pdf_file_name
508
509    ! Number of gpoints and bands in the shortwave and longwave - set
510    ! to zero as will be set properly later
511    integer :: n_g_sw = 0, n_g_lw = 0
512    integer :: n_bands_sw = 0, n_bands_lw = 0
513
514    ! Number of spectral points to save (equal either to the number of
515    ! g points or the number of bands
516    integer :: n_spec_sw = 0, n_spec_lw = 0
517
518    ! Dimensions to store variables that are only needed if longwave
519    ! scattering is included. "n_g_lw_if_scattering" is equal to
520    ! "n_g_lw" if aerosols are allowed to scatter in the longwave,
521    ! and zero otherwise. "n_bands_lw_if_scattering" is equal to
522    ! "n_bands_lw" if clouds are allowed to scatter in the longwave,
523    ! and zero otherwise.
524    integer :: n_g_lw_if_scattering = 0, n_bands_lw_if_scattering = 0
525
526    ! Treat clouds as horizontally homogeneous within the gribox
527    logical :: is_homogeneous = .false.
528
529    ! If the solvers are both "Cloudless" then we don't need to do any
530    ! cloud processing
531    logical :: do_clouds = .true.
532
533   contains
534     procedure :: read => read_config_from_namelist
535     procedure :: consolidate => consolidate_config
536     procedure :: set  => set_config
537     procedure :: print => print_config
538     procedure :: get_sw_weights
539     procedure :: define_sw_albedo_intervals
540     procedure :: define_lw_emiss_intervals
541     procedure :: consolidate_intervals
542
543  end type config_type
544
545!  procedure, private :: print_logical, print_real, print_int
546
547contains
548
549
550  !---------------------------------------------------------------------
551  ! This subroutine reads configuration data from a namelist file, and
552  ! anything that is not in the namelists will be set to default
553  ! values. If optional output argument "is_success" is present, then
554  ! on error (e.g. missing file) it will be set to .false.; if this
555  ! argument is missing then on error the program will be aborted. You
556  ! may either specify the file_name or the unit of an open file to
557  ! read, but not both.
558  subroutine read_config_from_namelist(this, file_name, unit, is_success)
559
560    use yomhook,      only : lhook, dr_hook
561    use radiation_io, only : nulout, nulerr, nulrad, radiation_abort
562
563    class(config_type), intent(inout)         :: this
564    character(*),       intent(in),  optional :: file_name
565    integer,            intent(in),  optional :: unit
566    logical,            intent(out), optional :: is_success
567
568    integer :: iosopen, iosread ! Status after calling open and read
569
570    ! The following variables are read from the namelists and map
571    ! directly onto members of the config_type derived type
572
573    ! To be read from the radiation_config namelist
574    logical :: do_sw, do_lw, do_clear, do_sw_direct
575    logical :: do_3d_effects, use_expm_everywhere, use_aerosols
576    logical :: do_lw_side_emissivity
577    logical :: do_3d_lw_multilayer_effects, do_fu_lw_ice_optics_bug
578    logical :: do_lw_aerosol_scattering, do_lw_cloud_scattering
579    logical :: do_save_radiative_properties, do_save_spectral_flux
580    logical :: do_save_gpoint_flux, do_surface_sw_spectral_flux
581    logical :: use_beta_overlap, do_lw_derivatives
582    logical :: do_sw_delta_scaling_with_gases
583    logical :: do_canopy_fluxes_sw, do_canopy_fluxes_lw
584    logical :: use_canopy_full_spectrum_sw, use_canopy_full_spectrum_lw
585    logical :: do_canopy_gases_sw, do_canopy_gases_lw
586    integer :: n_regions, iverbose, iverbosesetup, n_aerosol_types
587    real(jprb):: mono_lw_wavelength, mono_lw_total_od, mono_sw_total_od
588    real(jprb):: mono_lw_single_scattering_albedo, mono_sw_single_scattering_albedo
589    real(jprb):: mono_lw_asymmetry_factor, mono_sw_asymmetry_factor
590    real(jprb):: cloud_inhom_decorr_scaling, cloud_fraction_threshold
591    real(jprb):: clear_to_thick_fraction, max_gas_od_3d, max_cloud_od
592    real(jprb):: cloud_mixing_ratio_threshold, overhead_sun_factor
593    real(jprb):: max_3d_transfer_rate, min_cloud_effective_size
594    real(jprb):: overhang_factor, encroachment_scaling
595    character(511) :: directory_name, aerosol_optics_override_file_name
596    character(511) :: liq_optics_override_file_name, ice_optics_override_file_name
597    character(511) :: cloud_pdf_override_file_name
598    character(63)  :: liquid_model_name, ice_model_name, gas_model_name
599    character(63)  :: sw_solver_name, lw_solver_name, overlap_scheme_name
600    character(63)  :: sw_entrapment_name, sw_encroachment_name, cloud_pdf_shape_name
601    integer :: i_aerosol_type_map(NMaxAerosolTypes) ! More than 256 is an error
602
603    logical :: do_nearest_spectral_sw_albedo = .true.
604    logical :: do_nearest_spectral_lw_emiss  = .true.
605    real(jprb) :: sw_albedo_wavelength_bound(NMaxAlbedoIntervals-1)
606    real(jprb) :: lw_emiss_wavelength_bound( NMaxAlbedoIntervals-1)
607    integer :: i_sw_albedo_index(NMaxAlbedoIntervals)
608    integer :: i_lw_emiss_index (NMaxAlbedoIntervals)
609
610    integer :: iunit ! Unit number of namelist file
611
612    logical :: lldeb_conf = .false.
613
614    namelist /radiation/ do_sw, do_lw, do_sw_direct, &
615         &  do_3d_effects, do_lw_side_emissivity, do_clear, &
616         &  do_save_radiative_properties, sw_entrapment_name, sw_encroachment_name, &
617         &  do_3d_lw_multilayer_effects, do_fu_lw_ice_optics_bug, &
618         &  do_save_spectral_flux, do_save_gpoint_flux, &
619         &  do_surface_sw_spectral_flux, do_lw_derivatives, &
620         &  do_lw_aerosol_scattering, do_lw_cloud_scattering, &
621         &  n_regions, directory_name, gas_model_name, &
622         &  ice_optics_override_file_name, liq_optics_override_file_name, &
623         &  aerosol_optics_override_file_name, cloud_pdf_override_file_name, &
624         &  liquid_model_name, ice_model_name, max_3d_transfer_rate, &
625         &  min_cloud_effective_size, overhang_factor, encroachment_scaling, &
626         &  use_canopy_full_spectrum_sw, use_canopy_full_spectrum_lw, &
627         &  do_canopy_fluxes_sw, do_canopy_fluxes_lw, &
628         &  do_canopy_gases_sw, do_canopy_gases_lw, &
629         &  do_sw_delta_scaling_with_gases, overlap_scheme_name, &
630         &  sw_solver_name, lw_solver_name, use_beta_overlap, &
631         &  use_expm_everywhere, iverbose, iverbosesetup, &
632         &  cloud_inhom_decorr_scaling, cloud_fraction_threshold, &
633         &  clear_to_thick_fraction, max_gas_od_3d, max_cloud_od, &
634         &  cloud_mixing_ratio_threshold, overhead_sun_factor, &
635         &  n_aerosol_types, i_aerosol_type_map, use_aerosols, &
636         &  mono_lw_wavelength, mono_lw_total_od, mono_sw_total_od, &
637         &  mono_lw_single_scattering_albedo, mono_sw_single_scattering_albedo, &
638         &  mono_lw_asymmetry_factor, mono_sw_asymmetry_factor, &
639         &  cloud_pdf_shape_name, &
640         &  do_nearest_spectral_sw_albedo, do_nearest_spectral_lw_emiss, &
641         &  sw_albedo_wavelength_bound, lw_emiss_wavelength_bound, &
642         &  i_sw_albedo_index, i_lw_emiss_index
643
644    real(jprb) :: hook_handle
645
646    if (lhook) call dr_hook('radiation_config:read',0,hook_handle)
647
648    ! Copy default values from the original structure
649    do_sw = this%do_sw
650    do_lw = this%do_lw
651    do_sw_direct = this%do_sw_direct
652    do_3d_effects = this%do_3d_effects
653    do_3d_lw_multilayer_effects = this%do_3d_lw_multilayer_effects
654    do_lw_side_emissivity = this%do_lw_side_emissivity
655    do_clear = this%do_clear
656    do_lw_aerosol_scattering = this%do_lw_aerosol_scattering
657    do_lw_cloud_scattering = this%do_lw_cloud_scattering
658    do_sw_delta_scaling_with_gases = this%do_sw_delta_scaling_with_gases
659    do_fu_lw_ice_optics_bug = this%do_fu_lw_ice_optics_bug
660    do_canopy_fluxes_sw = this%do_canopy_fluxes_sw
661    do_canopy_fluxes_lw = this%do_canopy_fluxes_lw
662    use_canopy_full_spectrum_sw = this%use_canopy_full_spectrum_sw
663    use_canopy_full_spectrum_lw = this%use_canopy_full_spectrum_lw
664    do_canopy_gases_sw = this%do_canopy_gases_sw
665    do_canopy_gases_lw = this%do_canopy_gases_lw
666    n_regions = this%nregions
667    directory_name = this%directory_name
668    cloud_pdf_override_file_name = this%cloud_pdf_override_file_name
669    liq_optics_override_file_name = this%liq_optics_override_file_name
670    ice_optics_override_file_name = this%ice_optics_override_file_name
671    aerosol_optics_override_file_name = this%aerosol_optics_override_file_name
672    use_expm_everywhere = this%use_expm_everywhere
673    use_aerosols = this%use_aerosols
674    do_save_radiative_properties = this%do_save_radiative_properties
675    do_save_spectral_flux = this%do_save_spectral_flux
676    do_save_gpoint_flux = this%do_save_gpoint_flux
677    do_lw_derivatives = this%do_lw_derivatives
678    do_surface_sw_spectral_flux = this%do_surface_sw_spectral_flux
679    iverbose = this%iverbose
680    iverbosesetup = this%iverbosesetup
681    cloud_fraction_threshold = this%cloud_fraction_threshold
682    cloud_mixing_ratio_threshold = this%cloud_mixing_ratio_threshold
683    use_beta_overlap = this%use_beta_overlap
684    cloud_inhom_decorr_scaling = this%cloud_inhom_decorr_scaling
685    clear_to_thick_fraction = this%clear_to_thick_fraction
686    overhead_sun_factor = this%overhead_sun_factor
687    max_gas_od_3d = this%max_gas_od_3d
688    max_cloud_od = this%max_cloud_od
689    max_3d_transfer_rate = this%max_3d_transfer_rate
690    min_cloud_effective_size = this%min_cloud_effective_size
691    overhang_factor = this%overhang_factor
692    encroachment_scaling = -1.0_jprb
693    gas_model_name = '' !DefaultGasModelName
694    liquid_model_name = '' !DefaultLiquidModelName
695    ice_model_name = '' !DefaultIceModelName
696    sw_solver_name = '' !DefaultSwSolverName
697    lw_solver_name = '' !DefaultLwSolverName
698    sw_entrapment_name = ''
699    sw_encroachment_name = ''
700    overlap_scheme_name = ''
701    cloud_pdf_shape_name = ''
702    n_aerosol_types = this%n_aerosol_types
703    mono_lw_wavelength = this%mono_lw_wavelength
704    mono_lw_total_od = this%mono_lw_total_od
705    mono_sw_total_od = this%mono_sw_total_od
706    mono_lw_single_scattering_albedo = this%mono_lw_single_scattering_albedo
707    mono_sw_single_scattering_albedo = this%mono_sw_single_scattering_albedo
708    mono_lw_asymmetry_factor = this%mono_lw_asymmetry_factor
709    mono_sw_asymmetry_factor = this%mono_sw_asymmetry_factor
710    i_aerosol_type_map = this%i_aerosol_type_map
711    do_nearest_spectral_sw_albedo = this%do_nearest_spectral_sw_albedo
712    do_nearest_spectral_lw_emiss  = this%do_nearest_spectral_lw_emiss
713    sw_albedo_wavelength_bound    = this%sw_albedo_wavelength_bound
714    lw_emiss_wavelength_bound     = this%lw_emiss_wavelength_bound
715    i_sw_albedo_index             = this%i_sw_albedo_index
716    i_lw_emiss_index              = this%i_lw_emiss_index
717
718    if (present(file_name) .and. present(unit)) then
719      write(nulerr,'(a)') '*** Error: cannot specify both file_name and unit in call to config_type%read'
720      call radiation_abort('Radiation configuration error')
721    else if (.not. present(file_name) .and. .not. present(unit)) then
722      write(nulerr,'(a)') '*** Error: neither file_name nor unit specified in call to config_type%read'
723      call radiation_abort('Radiation configuration error')
724    end if
725
726    if (present(file_name)) then
727      ! Open the namelist file
728      iunit = nulrad
729      open(unit=iunit, iostat=iosopen, file=trim(file_name))
730    else
731      ! Assume that iunit represents and open file
732      iosopen = 0
733      iunit = unit
734    end if
735
736    if (iosopen /= 0) then
737      ! An error occurred opening the file
738      if (present(is_success)) then
739        is_success = .false.
740        ! We now continue the subroutine so that the default values
741        ! are placed in the config structure
742      else
743        write(nulerr,'(a,a,a)') '*** Error: namelist file "', &
744             &                trim(file_name), '" not found'
745        call radiation_abort('Radiation configuration error')
746      end if
747    else
748      read(unit=iunit, iostat=iosread, nml=radiation)
749      if (iosread /= 0) then
750        ! An error occurred reading the file
751        if (present(is_success)) then
752          is_success = .false.
753          ! We now continue the subroutine so that the default values
754          ! are placed in the config structure
755        else if (present(file_name)) then
756          write(nulerr,'(a,a,a)') '*** Error reading namelist "radiation" from file "', &
757               &      trim(file_name), '"'
758          close(unit=iunit)
759          call radiation_abort('Radiation configuration error')
760        else
761          write(nulerr,'(a,i0)') '*** Error reading namelist "radiation" from unit ', &
762               &      iunit
763          call radiation_abort('Radiation configuration error')
764        end if
765      end if
766
767      if (present(file_name)) then
768        close(unit=iunit)
769      end if
770    end if
771
772    ! Copy namelist data into configuration object
773
774    ! Start with verbosity levels, which should be within limits
775    if (iverbose < 0) then
776      iverbose = 0
777    end if
778    this%iverbose = iverbose
779
780    if (iverbosesetup < 0) then
781      iverbosesetup = 0
782    end if
783    this%iverbosesetup = iverbosesetup
784
785    this%do_lw = do_lw
786    this%do_sw = do_sw
787    this%do_clear = do_clear
788    this%do_sw_direct = do_sw_direct
789    this%do_3d_effects = do_3d_effects
790    this%do_3d_lw_multilayer_effects = do_3d_lw_multilayer_effects
791    this%do_lw_side_emissivity = do_lw_side_emissivity
792    this%use_expm_everywhere = use_expm_everywhere
793    this%use_aerosols = use_aerosols
794    this%do_lw_cloud_scattering = do_lw_cloud_scattering
795    this%do_lw_aerosol_scattering = do_lw_aerosol_scattering
796    this%nregions = n_regions
797    this%do_surface_sw_spectral_flux = do_surface_sw_spectral_flux
798    this%do_sw_delta_scaling_with_gases = do_sw_delta_scaling_with_gases
799    this%do_fu_lw_ice_optics_bug = do_fu_lw_ice_optics_bug
800    this%do_canopy_fluxes_sw = do_canopy_fluxes_sw
801    this%do_canopy_fluxes_lw = do_canopy_fluxes_lw
802    this%use_canopy_full_spectrum_sw = use_canopy_full_spectrum_sw
803    this%use_canopy_full_spectrum_lw = use_canopy_full_spectrum_lw
804    this%do_canopy_gases_sw = do_canopy_gases_sw
805    this%do_canopy_gases_lw = do_canopy_gases_lw
806    this%mono_lw_wavelength = mono_lw_wavelength
807    this%mono_lw_total_od = mono_lw_total_od
808    this%mono_sw_total_od = mono_sw_total_od
809    this%mono_lw_single_scattering_albedo = mono_lw_single_scattering_albedo
810    this%mono_sw_single_scattering_albedo = mono_sw_single_scattering_albedo
811    this%mono_lw_asymmetry_factor = mono_lw_asymmetry_factor
812    this%mono_sw_asymmetry_factor = mono_sw_asymmetry_factor
813    this%use_beta_overlap = use_beta_overlap
814    this%cloud_inhom_decorr_scaling = cloud_inhom_decorr_scaling
815    this%clear_to_thick_fraction = clear_to_thick_fraction
816    this%overhead_sun_factor = overhead_sun_factor
817    this%max_gas_od_3d = max_gas_od_3d
818    this%max_cloud_od = max_cloud_od
819    this%max_3d_transfer_rate = max_3d_transfer_rate
820    this%min_cloud_effective_size = max(1.0e-6_jprb, min_cloud_effective_size)
821    if (encroachment_scaling >= 0.0_jprb) then
822      this%overhang_factor = encroachment_scaling
823      if (iverbose >= 1) then
824        write(nulout, '(a)') 'Warning: radiation namelist parameter "encroachment_scaling" is deprecated: use "overhang_factor"'
825      end if
826    else
827      this%overhang_factor = overhang_factor
828    end if
829    this%directory_name = directory_name
830    this%cloud_pdf_override_file_name = cloud_pdf_override_file_name
831    this%liq_optics_override_file_name = liq_optics_override_file_name
832    this%ice_optics_override_file_name = ice_optics_override_file_name
833    this%aerosol_optics_override_file_name = aerosol_optics_override_file_name
834    this%cloud_fraction_threshold = cloud_fraction_threshold
835    this%cloud_mixing_ratio_threshold = cloud_mixing_ratio_threshold
836    this%n_aerosol_types = n_aerosol_types
837    this%do_save_radiative_properties = do_save_radiative_properties
838    this%do_lw_derivatives = do_lw_derivatives
839    this%do_save_spectral_flux = do_save_spectral_flux
840    this%do_save_gpoint_flux = do_save_gpoint_flux
841    this%do_nearest_spectral_sw_albedo = do_nearest_spectral_sw_albedo
842    this%do_nearest_spectral_lw_emiss  = do_nearest_spectral_lw_emiss
843    this%sw_albedo_wavelength_bound    = sw_albedo_wavelength_bound
844    this%lw_emiss_wavelength_bound     = lw_emiss_wavelength_bound
845    this%i_sw_albedo_index             = i_sw_albedo_index
846    this%i_lw_emiss_index              = i_lw_emiss_index
847
848! AI mars 2022
849if (lldeb_conf) then
850print*,'**************PARAMETRES DE CONFIGURATION OFFLINE*******************'
851print*,'config%iverbosesetup   = ', iverbosesetup
852print*,'config%do_lw   = ', do_lw
853print*,'config%do_sw   = ', do_sw
854print*,'config%do_clear   = ', do_clear
855print*,'config%do_sw_direct   = ', do_sw_direct
856print*,'config%do_3d_effects   = ', do_3d_effects
857print*,'config%do_3d_lw_multilayer_effects   = ', do_3d_lw_multilayer_effects
858print*,'config%do_lw_side_emissivity   = ', do_lw_side_emissivity
859print*,'config%use_expm_everywhere   = ', use_expm_everywhere
860print*,'config%use_aerosols   = ', use_aerosols
861print*,'config%do_lw_cloud_scattering   = ', do_lw_cloud_scattering
862print*,'config%do_lw_aerosol_scattering   = ', do_lw_aerosol_scattering
863print*,'config%nregions   = ', n_regions
864print*,'config%do_surface_sw_spectral_flux   = ', do_surface_sw_spectral_flux
865print*,'config%do_sw_delta_scaling_with_gases   = ', &
866do_sw_delta_scaling_with_gases
867print*,'config%do_fu_lw_ice_optics_bug   = ', do_fu_lw_ice_optics_bug
868print*,'config%do_canopy_fluxes_sw   = ', do_canopy_fluxes_sw
869print*,'config%do_canopy_fluxes_lw   = ', do_canopy_fluxes_lw
870print*,'config%use_canopy_full_spectrum_sw   = ', use_canopy_full_spectrum_sw
871print*,'config%use_canopy_full_spectrum_lw   = ', use_canopy_full_spectrum_lw
872print*,'config%do_canopy_gases_sw   = ', do_canopy_gases_sw
873print*,'config%do_canopy_gases_lw   = ', do_canopy_gases_lw
874print*,'config%mono_lw_wavelength   = ', mono_lw_wavelength
875print*,'config%mono_lw_total_od   = ', mono_lw_total_od
876print*,'config%mono_sw_total_od   = ', mono_sw_total_od
877print*,'config%mono_lw_single_scattering_albedo   = ', &
878mono_lw_single_scattering_albedo
879print*,'config%mono_sw_single_scattering_albedo   = ', &
880mono_sw_single_scattering_albedo
881print*,'config%mono_lw_asymmetry_factor   = ', mono_lw_asymmetry_factor
882print*,'config%mono_sw_asymmetry_factor   = ', mono_sw_asymmetry_factor
883print*,'config%use_beta_overlap   = ', use_beta_overlap
884print*,'config%cloud_inhom_decorr_scaling   = ', cloud_inhom_decorr_scaling
885print*,'config%clear_to_thick_fraction   = ', clear_to_thick_fraction
886print*,'config%overhead_sun_factor   = ', overhead_sun_factor
887print*,'config%max_gas_od_3d   = ', max_gas_od_3d
888print*,'config%max_cloud_od   = ', max_cloud_od
889print*,'config%max_3d_transfer_rate   = ', max_3d_transfer_rate
890print*,'config%min_cloud_effective_size   = ', &
891max(1.0e-6_jprb,min_cloud_effective_size)
892print*,'config%overhang_factor   = ', encroachment_scaling
893
894print*,'config%directory_name  = ',directory_name
895print*,'config%cloud_pdf_override_file_name  = ',cloud_pdf_override_file_name
896print*,'config%liq_optics_override_file_name  = ',liq_optics_override_file_name
897print*,'config%ice_optics_override_file_name  = ',ice_optics_override_file_name
898print*,'config%aerosol_optics_override_file_name  = ', &
899aerosol_optics_override_file_name
900print*,'config%cloud_fraction_threshold  = ',cloud_fraction_threshold
901print*,'config%cloud_mixing_ratio_threshold  = ',cloud_mixing_ratio_threshold
902print*,'config%n_aerosol_types  = ',n_aerosol_types
903print*,'config%do_save_radiative_properties  = ',do_save_radiative_properties
904print*,'config%do_lw_derivatives  = ',do_lw_derivatives
905print*,'config%do_save_spectral_flux  = ',do_save_spectral_flux
906print*,'config%do_save_gpoint_flux  = ',do_save_gpoint_flux
907print*,'config%do_nearest_spectral_sw_albedo  = ',do_nearest_spectral_sw_albedo
908print*,'config%do_nearest_spectral_lw_emiss   = ',do_nearest_spectral_lw_emiss
909print*,'config%sw_albedo_wavelength_bound     = ',sw_albedo_wavelength_bound
910print*,'config%lw_emiss_wavelength_bound      = ',lw_emiss_wavelength_bound
911print*,'config%i_sw_albedo_index              = ',i_sw_albedo_index
912print*,'config%i_lw_emiss_index               = ',i_lw_emiss_index
913print*,'************************************************************************'
914endif
915    if (do_save_gpoint_flux) then
916      ! Saving the fluxes every g-point overrides saving as averaged
917      ! in a band, but this%do_save_spectral_flux needs to be TRUE as
918      ! it is tested inside the solver routines to decide whether to
919      ! save anything
920      this%do_save_spectral_flux = .true.
921      print*,'config%do_save_spectral_flux = .true.'
922    end if
923
924    ! Determine liquid optics model
925    call get_enum_code(liquid_model_name, LiquidModelName, &
926         &            'liquid_model_name', this%i_liq_model)
927    print*,'config%i_liq_model =', this%i_liq_model
928
929    ! Determine ice optics model
930    call get_enum_code(ice_model_name, IceModelName, &
931         &            'ice_model_name', this%i_ice_model)
932    print*,'config%i_ice_model =', this%i_ice_model
933    ! Determine gas optics model
934    call get_enum_code(gas_model_name, GasModelName, &
935         &            'gas_model_name', this%i_gas_model)
936    print*,'config%%i_gas_model = ', this%i_gas_model
937
938    ! Determine solvers
939    call get_enum_code(sw_solver_name, SolverName, &
940         &            'sw_solver_name', this%i_solver_sw)
941    print*,'config%i_solver_sw = ', this%i_solver_sw
942    call get_enum_code(lw_solver_name, SolverName, &
943         &            'lw_solver_name', this%i_solver_lw)
944    print*,'config%i_solver_lw = ', this%i_solver_lw
945    if (len_trim(sw_encroachment_name) > 1) then
946      call get_enum_code(sw_encroachment_name, EncroachmentName, &
947           &             'sw_encroachment_name', this%i_3d_sw_entrapment)
948      write(nulout, '(a)') 'Warning: radiation namelist string "sw_encroachment_name" is deprecated: use "sw_entrapment_name"'
949    else
950      call get_enum_code(sw_entrapment_name, EntrapmentName, &
951           &             'sw_entrapment_name', this%i_3d_sw_entrapment)
952      print*,'config%i_3d_sw_entrapment = ', this%i_3d_sw_entrapment
953    end if
954
955    ! Determine overlap scheme
956    call get_enum_code(overlap_scheme_name, OverlapName, &
957         &             'overlap_scheme_name', this%i_overlap_scheme)
958    print*,'config%i_overlap_scheme = ', this%i_overlap_scheme
959    ! Determine cloud PDF shape
960    call get_enum_code(cloud_pdf_shape_name, PdfShapeName, &
961         &             'cloud_pdf_shape_name', this%i_cloud_pdf_shape)
962    print*,'config%i_cloud_pdf_shape = ', this%i_cloud_pdf_shape
963    this%i_aerosol_type_map = 0
964    if (this%use_aerosols) then
965      this%i_aerosol_type_map(1:n_aerosol_types) &
966           &  = i_aerosol_type_map(1:n_aerosol_types)
967      print*,'config%i_aerosol_type_map = ', this%i_aerosol_type_map
968    end if
969
970    ! Will clouds be used at all?
971    if ((this%do_sw .and. this%i_solver_sw /= ISolverCloudless) &
972         &  .or. (this%do_lw .and. this%i_solver_lw /= ISolverCloudless)) then
973      this%do_clouds = .true.
974    else
975      this%do_clouds = .false.
976    end if
977    print*,'config%do_clouds = ', this%do_clouds
978
979    ! Normal subroutine exit
980    if (present(is_success)) then
981      is_success = .true.
982    end if
983
984    if (lhook) call dr_hook('radiation_config:read',1,hook_handle)
985
986  end subroutine read_config_from_namelist
987
988
989  !---------------------------------------------------------------------
990  ! This routine is called by radiation_interface:setup_radiation and
991  ! it converts the user specified options into some more specific
992  ! data such as data file names
993  subroutine consolidate_config(this)
994
995    use yomhook,      only : lhook, dr_hook
996    use radiation_io, only : nulout, nulerr, radiation_abort
997
998    class(config_type), intent(inout)         :: this
999
1000    real(jprb) :: hook_handle
1001
1002    if (lhook) call dr_hook('radiation_config:consolidate',0,hook_handle)
1003
1004    ! Check consistency of models
1005    if (this%do_canopy_fluxes_sw .and. .not. this%do_surface_sw_spectral_flux) then
1006      if (this%iverbosesetup >= 1) then
1007        write(nulout,'(a)') 'Warning: turning on do_surface_sw_spectral_flux as required by do_canopy_fluxes_sw'
1008      end if
1009      this%do_surface_sw_spectral_flux = .true.
1010    end if
1011
1012    ! Will clouds be used at all?
1013    if ((this%do_sw .and. this%i_solver_sw /= ISolverCloudless) &
1014         &  .or. (this%do_lw .and. this%i_solver_lw /= ISolverCloudless)) then
1015      this%do_clouds = .true.
1016    else
1017      this%do_clouds = .false.
1018    end if
1019
1020    ! SPARTACUS only works with Exp-Ran overlap scheme
1021    if ((       this%i_solver_sw == ISolverSPARTACUS &
1022         & .or. this%i_solver_lw == ISolverSPARTACUS &
1023         & .or. this%i_solver_sw == ISolverTripleclouds &
1024         & .or. this%i_solver_lw == ISolverTripleclouds) &
1025         & .and. this%i_overlap_scheme /= IOverlapExponentialRandom) then
1026      write(nulerr,'(a)') '*** Error: SPARTACUS/Tripleclouds solvers can only do Exponential-Random overlap'
1027      call radiation_abort('Radiation configuration error')
1028
1029    end if
1030
1031    ! Set aerosol optics file name
1032    if (len_trim(this%aerosol_optics_override_file_name) > 0) then
1033      if (this%aerosol_optics_override_file_name(1:1) == '/') then
1034        this%aerosol_optics_file_name = trim(this%aerosol_optics_override_file_name)
1035      else
1036        this%aerosol_optics_file_name = trim(this%directory_name) &
1037             &  // '/' // trim(this%aerosol_optics_override_file_name)
1038      end if
1039    else
1040      ! In the IFS, the aerosol optics file should be specified in
1041      ! ifs/module/radiation_setup.F90, not here
1042      this%aerosol_optics_file_name &
1043           &   = trim(this%directory_name) // "/aerosol_ifs_rrtm_45R2.nc"
1044    end if
1045
1046    ! Set liquid optics file name
1047    if (len_trim(this%liq_optics_override_file_name) > 0) then
1048      if (this%liq_optics_override_file_name(1:1) == '/') then
1049        this%liq_optics_file_name = trim(this%liq_optics_override_file_name)
1050      else
1051        this%liq_optics_file_name = trim(this%directory_name) &
1052             &  // '/' // trim(this%liq_optics_override_file_name)
1053      end if
1054    else if (this%i_liq_model == ILiquidModelSOCRATES) then
1055      this%liq_optics_file_name &
1056           &  = trim(this%directory_name) // "/socrates_droplet_scattering_rrtm.nc"
1057    else if (this%i_liq_model == ILiquidModelSlingo) then
1058      this%liq_optics_file_name &
1059           &  = trim(this%directory_name) // "/slingo_droplet_scattering_rrtm.nc"
1060    end if
1061
1062    ! Set ice optics file name
1063    if (len_trim(this%ice_optics_override_file_name) > 0) then
1064      if (this%ice_optics_override_file_name(1:1) == '/') then
1065        this%ice_optics_file_name = trim(this%ice_optics_override_file_name)
1066      else
1067        this%ice_optics_file_name = trim(this%directory_name) &
1068             &  // '/' // trim(this%ice_optics_override_file_name)
1069      end if
1070    else if (this%i_ice_model == IIceModelFu) then
1071      this%ice_optics_file_name &
1072           &   = trim(this%directory_name) // "/fu_ice_scattering_rrtm.nc"
1073    else if (this%i_ice_model == IIceModelBaran) then
1074      this%ice_optics_file_name &
1075           &   = trim(this%directory_name) // "/baran_ice_scattering_rrtm.nc"
1076    else if (this%i_ice_model == IIceModelBaran2016) then
1077      this%ice_optics_file_name &
1078           &   = trim(this%directory_name) // "/baran2016_ice_scattering_rrtm.nc"
1079    else if (this%i_ice_model == IIceModelBaran2017) then
1080      this%ice_optics_file_name &
1081           &   = trim(this%directory_name) // "/baran2017_ice_scattering_rrtm.nc"
1082    else if (this%i_ice_model == IIceModelYi) then
1083      this%ice_optics_file_name &
1084           &   = trim(this%directory_name) // "/yi_ice_scattering_rrtm.nc"
1085    end if
1086
1087    ! Set cloud-water PDF look-up table file name
1088    if (len_trim(this%cloud_pdf_override_file_name) > 0) then
1089      if (this%cloud_pdf_override_file_name(1:1) == '/') then
1090        this%cloud_pdf_file_name = trim(this%cloud_pdf_override_file_name)
1091      else
1092        this%cloud_pdf_file_name = trim(this%directory_name) &
1093             &  // '/' // trim(this%cloud_pdf_override_file_name)
1094      end if
1095    elseif (this%i_cloud_pdf_shape == IPdfShapeLognormal) then
1096      this%cloud_pdf_file_name = trim(this%directory_name) // "/mcica_lognormal.nc"
1097    else
1098      this%cloud_pdf_file_name = trim(this%directory_name) // "/mcica_gamma.nc"
1099    end if
1100
1101    ! Aerosol data
1102    if (this%n_aerosol_types < 0 &
1103         &  .or. this%n_aerosol_types > NMaxAerosolTypes) then
1104      write(nulerr,'(a,i0)') '*** Error: number of aerosol types must be between 0 and ', &
1105           &  NMaxAerosolTypes
1106      call radiation_abort('Radiation configuration error')
1107    end if
1108
1109    if (this%use_aerosols .and. this%n_aerosol_types == 0) then
1110      if (this%iverbosesetup >= 2) then
1111        write(nulout, '(a)') 'Aerosols on but n_aerosol_types=0: optical properties to be computed outside ecRad'
1112      end if
1113    end if
1114
1115    ! In the monochromatic case we need to override the liquid, ice
1116    ! and aerosol models to ensure compatibility
1117    if (this%i_gas_model == IGasModelMonochromatic) then
1118      this%i_liq_model = ILiquidModelMonochromatic
1119      this%i_ice_model = IIceModelMonochromatic
1120      this%use_aerosols = .false.
1121    end if
1122
1123    ! McICA solver currently can't store full profiles of spectral fluxes
1124    if (this%i_solver_sw == ISolverMcICA) then
1125      this%do_save_spectral_flux = .false.
1126    end if
1127
1128    if (this%i_solver_sw == ISolverSPARTACUS .and. this%do_sw_delta_scaling_with_gases) then
1129      write(nulerr,'(a)') '*** Error: SW delta-Eddington scaling with gases not possible with SPARTACUS solver'
1130      call radiation_abort('Radiation configuration error')
1131    end if
1132
1133    if ((this%do_lw .and. this%do_sw) .and. &
1134         & (     (      this%i_solver_sw == ISolverHomogeneous  &
1135         &        .and. this%i_solver_lw /= ISolverHomogeneous) &
1136         &  .or. (      this%i_solver_sw /= ISolverHomogeneous  &
1137         &        .and. this%i_solver_lw == ISolverHomogeneous) &
1138         & ) ) then
1139      write(nulerr,'(a)') '*** Error: if one solver is "Homogeneous" then the other must be'
1140      call radiation_abort('Radiation configuration error')
1141    end if
1142
1143    ! Set is_homogeneous if the active solvers are homogeneous, since
1144    ! this affects how "in-cloud" water contents are computed
1145    if (        (this%do_sw .and. this%i_solver_sw == ISolverHomogeneous) &
1146         & .or. (this%do_lw .and. this%i_solver_lw == ISolverHomogeneous)) then
1147      this%is_homogeneous = .true.
1148    end if
1149
1150    this%is_consolidated = .true.
1151
1152    if (lhook) call dr_hook('radiation_config:consolidate',1,hook_handle)
1153
1154  end subroutine consolidate_config
1155
1156
1157  !---------------------------------------------------------------------
1158  ! This subroutine sets members of the configuration object via
1159  ! optional arguments, and any member not specified is left
1160  ! untouched. Therefore, this should be called after taking data from
1161  ! the namelist.
1162  subroutine set_config(config, directory_name, &
1163       &  do_lw, do_sw, &
1164       &  do_lw_aerosol_scattering, do_lw_cloud_scattering, &
1165       &  do_sw_direct)
1166
1167    class(config_type), intent(inout):: config
1168    character(len=*), intent(in), optional  :: directory_name
1169    logical, intent(in), optional           :: do_lw, do_sw
1170    logical, intent(in), optional           :: do_lw_aerosol_scattering
1171    logical, intent(in), optional           :: do_lw_cloud_scattering
1172    logical, intent(in), optional           :: do_sw_direct
1173
1174    if (present(do_lw)) then
1175       config%do_lw = do_lw
1176    end if
1177
1178    if(present(do_sw)) then
1179       config%do_sw = do_sw
1180    end if
1181
1182    if (present(do_sw_direct)) then
1183       config%do_sw_direct = do_sw_direct
1184    end if
1185
1186    if (present(directory_name)) then
1187       config%directory_name = trim(directory_name)
1188    end if
1189
1190    if (present(do_lw_aerosol_scattering)) then
1191       config%do_lw_aerosol_scattering = .true.
1192    end if
1193
1194    if (present(do_lw_cloud_scattering)) then
1195       config%do_lw_cloud_scattering = .true.
1196    end if
1197
1198  end subroutine set_config
1199
1200
1201  !---------------------------------------------------------------------
1202  ! Print configuration information to standard output
1203  subroutine print_config(this, iverbose)
1204
1205    use radiation_io, only : nulout
1206
1207    class(config_type), intent(in) :: this
1208
1209    integer, optional,  intent(in) :: iverbose
1210    integer                        :: i_local_verbose
1211
1212    if (present(iverbose)) then
1213      i_local_verbose = iverbose
1214    else
1215      i_local_verbose = this%iverbose
1216    end if
1217
1218    if (i_local_verbose >= 2) then
1219      !---------------------------------------------------------------------
1220      write(nulout, '(a)') 'General settings:'
1221      write(nulout, '(a,a,a)') '  Data files expected in "', &
1222           &                   trim(this%directory_name), '"'
1223      call print_logical('  Clear-sky calculations are', 'do_clear', this%do_clear)
1224      call print_logical('  Saving intermediate radiative properties', &
1225           &   'do_save_radiative_properties', this%do_save_radiative_properties)
1226      call print_logical('  Saving spectral flux profiles', &
1227           &   'do_save_spectral_flux', this%do_save_spectral_flux)
1228      call print_enum('  Gas model is', GasModelName, 'i_gas_model', &
1229           &          this%i_gas_model)
1230      call print_logical('  Aerosols are', 'use_aerosols', this%use_aerosols)
1231      call print_logical('  Clouds are', 'do_clouds', this%do_clouds)
1232
1233      !---------------------------------------------------------------------
1234      write(nulout, '(a)') 'Surface settings:'
1235      if (this%do_sw) then
1236        call print_logical('  Saving surface shortwave spectral fluxes', &
1237             &   'do_surface_sw_spectral_flux', this%do_surface_sw_spectral_flux)
1238        call print_logical('  Saving surface shortwave fluxes in abledo bands', &
1239             &   'do_canopy_fluxes_sw', this%do_canopy_fluxes_sw)
1240      end if
1241      if (this%do_lw) then
1242        call print_logical('  Saving surface longwave fluxes in emissivity bands', &
1243             &   'do_canopy_fluxes_lw', this%do_canopy_fluxes_lw)
1244        call print_logical('  Longwave derivative calculation is', &
1245             &   'do_lw_derivatives',this%do_lw_derivatives)
1246      end if
1247      if (this%do_sw) then
1248        call print_logical('  Nearest-neighbour spectral albedo mapping', &
1249             &   'do_nearest_spectral_sw_albedo', this%do_nearest_spectral_sw_albedo)
1250      end if
1251      if (this%do_lw) then
1252        call print_logical('  Nearest-neighbour spectral emissivity mapping', &
1253             &   'do_nearest_spectral_lw_emiss', this%do_nearest_spectral_lw_emiss)
1254      end if
1255      !---------------------------------------------------------------------
1256      if (this%do_clouds) then
1257        write(nulout, '(a)') 'Cloud settings:'
1258        call print_real('  Cloud fraction threshold', &
1259             &   'cloud_fraction_threshold', this%cloud_fraction_threshold)
1260        call print_real('  Cloud mixing-ratio threshold', &
1261             &   'cloud_mixing_ratio_threshold', this%cloud_mixing_ratio_threshold)
1262        call print_enum('  Liquid optics scheme is', LiquidModelName, &
1263             &          'i_liq_model',this%i_liq_model)
1264        call print_enum('  Ice optics scheme is', IceModelName, &
1265             &          'i_ice_model',this%i_ice_model)
1266        if (this%i_ice_model == IIceModelFu) then
1267          call print_logical('  Longwave ice optics bug in Fu scheme is', &
1268               &   'do_fu_lw_ice_optics_bug',this%do_fu_lw_ice_optics_bug)
1269        end if
1270        call print_enum('  Cloud overlap scheme is', OverlapName, &
1271             &          'i_overlap_scheme',this%i_overlap_scheme)
1272        call print_logical('  Use "beta" overlap parameter is', &
1273             &   'use_beta_overlap', this%use_beta_overlap)
1274        call print_enum('  Cloud PDF shape is', PdfShapeName, &
1275             &          'i_cloud_pdf_shape',this%i_cloud_pdf_shape)
1276        call print_real('  Cloud inhom decorrelation scaling', &
1277             &   'cloud_inhom_decorr_scaling', this%cloud_inhom_decorr_scaling)
1278      end if
1279
1280      !---------------------------------------------------------------------
1281      write(nulout, '(a)') 'Solver settings:'
1282      if (this%do_sw) then
1283        call print_enum('  Shortwave solver is', SolverName, &
1284             &          'i_solver_sw', this%i_solver_sw)
1285       
1286        if (this%i_gas_model == IGasModelMonochromatic) then
1287          call print_real('  Shortwave atmospheric optical depth', &
1288               &   'mono_sw_total_od', this%mono_sw_total_od)
1289          call print_real('  Shortwave particulate single-scattering albedo', &
1290               &   'mono_sw_single_scattering_albedo', &
1291               &   this%mono_sw_single_scattering_albedo)
1292          call print_real('  Shortwave particulate asymmetry factor', &
1293               &   'mono_sw_asymmetry_factor', &
1294               &   this%mono_sw_asymmetry_factor)
1295        end if
1296        call print_logical('  Shortwave delta scaling after merge with gases', &
1297             &   'do_sw_delta_scaling_with_gases', &
1298             &   this%do_sw_delta_scaling_with_gases)
1299      else
1300        call print_logical('  Shortwave calculations are','do_sw',this%do_sw)
1301      end if
1302
1303      if (this%do_lw) then
1304        call print_enum('  Longwave solver is', SolverName, 'i_solver_lw', &
1305             &          this%i_solver_lw)
1306
1307        if (this%i_gas_model == IGasModelMonochromatic) then
1308          if (this%mono_lw_wavelength > 0.0_jprb) then
1309            call print_real('  Longwave effective wavelength (m)', &
1310                 &   'mono_lw_wavelength', this%mono_lw_wavelength)
1311          else
1312            write(nulout,'(a)') '  Longwave fluxes are broadband                              (mono_lw_wavelength<=0)'
1313          end if
1314          call print_real('  Longwave atmospheric optical depth', &
1315               &   'mono_lw_total_od', this%mono_lw_total_od) 
1316          call print_real('  Longwave particulate single-scattering albedo', &
1317               &   'mono_lw_single_scattering_albedo', &
1318               &   this%mono_lw_single_scattering_albedo)
1319          call print_real('  Longwave particulate asymmetry factor', &
1320               &   'mono_lw_asymmetry_factor', &
1321               &   this%mono_lw_asymmetry_factor)
1322        end if
1323        call print_logical('  Longwave cloud scattering is', &
1324             &   'do_lw_cloud_scattering',this%do_lw_cloud_scattering)
1325        call print_logical('  Longwave aerosol scattering is', &
1326             &   'do_lw_aerosol_scattering',this%do_lw_aerosol_scattering)
1327      else
1328        call print_logical('  Longwave calculations are','do_lw', this%do_lw)
1329      end if
1330
1331      if (this%i_solver_sw == ISolverSpartacus &
1332           &  .or. this%i_solver_lw == ISolverSpartacus) then
1333        write(nulout, '(a)') '  SPARTACUS options:'
1334        call print_integer('    Number of regions', 'n_regions', this%nregions)
1335        call print_real('    Max cloud optical depth per layer', &
1336             &   'max_cloud_od', this%max_cloud_od)
1337        call print_enum('    Shortwave entrapment is', EntrapmentName, &
1338             &          'i_3d_sw_entrapment', this%i_3d_sw_entrapment)
1339        call print_logical('    Multilayer longwave horizontal transport is', &
1340             'do_3d_lw_multilayer_effects', this%do_3d_lw_multilayer_effects)
1341        call print_logical('    Use matrix exponential everywhere is', &
1342             &   'use_expm_everywhere', this%use_expm_everywhere)
1343        call print_logical('    3D effects are', 'do_3d_effects', &
1344             &             this%do_3d_effects)
1345
1346        if (this%do_3d_effects) then
1347          call print_logical('    Longwave side emissivity parameterization is', &
1348               &  'do_lw_side_emissivity', this%do_lw_side_emissivity)
1349          call print_real('    Clear-to-thick edge fraction is', &
1350               &  'clear_to_thick_fraction', this%clear_to_thick_fraction)
1351          call print_real('    Overhead sun factor is', &
1352               &  'overhead_sun_factor', this%overhead_sun_factor)
1353          call print_real('    Max gas optical depth for 3D effects', &
1354               &   'max_gas_od_3d', this%max_gas_od_3d)
1355          call print_real('    Max 3D transfer rate', &
1356               &   'max_3d_transfer_rate', this%max_3d_transfer_rate)
1357          call print_real('    Min cloud effective size (m)', &
1358               &   'min_cloud_effective_size', this%min_cloud_effective_size)
1359          call print_real('    Overhang factor', &
1360               &   'overhang_factor', this%overhang_factor)
1361        end if
1362      end if
1363           
1364    end if
1365   
1366  end subroutine print_config
1367
1368
1369
1370  !---------------------------------------------------------------------
1371  ! In order to estimate UV and photosynthetically active radiation,
1372  ! we need weighted sum of fluxes considering wavelength range
1373  ! required.  This routine returns information for how to correctly
1374  ! weight output spectral fluxes for a range of input wavelengths.
1375  ! Note that this is approximate; internally it may be assumed that
1376  ! the energy is uniformly distributed in wavenumber space, for
1377  ! example.  If the character string "weighting_name" is present, and
1378  ! iverbose>=2, then information on the weighting will be provided on
1379  ! nulout.
1380  subroutine get_sw_weights(this, wavelength1, wavelength2, &
1381       &                    nweights, iband, weight, weighting_name)
1382
1383    use parkind1, only : jprb
1384    use radiation_io, only : nulout, nulerr, radiation_abort
1385
1386    class(config_type), intent(in) :: this
1387    ! Range of wavelengths to get weights for (m)
1388    real(jprb), intent(in) :: wavelength1, wavelength2
1389    ! Output number of weights needed
1390    integer,    intent(out)   :: nweights
1391    ! Only write to the first nweights of these arrays: they contain
1392    ! the indices to the non-zero bands, and the weight in each of
1393    ! those bands
1394    integer,    intent(out)   :: iband(:)
1395    real(jprb), intent(out)   :: weight(:)
1396    character(len=*), optional, intent(in) :: weighting_name
1397
1398    ! Internally we deal with wavenumber
1399    real(jprb) :: wavenumber1, wavenumber2 ! cm-1
1400
1401    integer :: jband ! Loop index for spectral band
1402
1403    if (this%n_bands_sw <= 0) then
1404      write(nulerr,'(a)') '*** Error: get_sw_weights called before number of shortwave bands set'
1405      call radiation_abort()     
1406    end if
1407
1408    ! Convert wavelength range (m) to wavenumber (cm-1)
1409    wavenumber1 = 0.01_jprb / wavelength2
1410    wavenumber2 = 0.01_jprb / wavelength1
1411
1412    nweights = 0
1413
1414    do jband = 1,this%n_bands_sw
1415      if (wavenumber1 < this%wavenumber2_sw(jband) &
1416           &  .and. wavenumber2 > this%wavenumber1_sw(jband)) then
1417        nweights = nweights+1
1418        iband(nweights) = jband
1419        weight(nweights) = (min(wavenumber2,this%wavenumber2_sw(jband)) &
1420             &         - max(wavenumber1,this%wavenumber1_sw(jband))) &
1421             & / (this%wavenumber2_sw(jband) - this%wavenumber1_sw(jband))
1422      end if
1423    end do
1424
1425    if (nweights == 0) then
1426      write(nulerr,'(a,e8.4,a,e8.4,a)') '*** Error: wavelength range ', &
1427           &  wavelength1, ' to ', wavelength2, ' m is outside shortwave band'
1428      call radiation_abort()
1429    else if (this%iverbosesetup >= 2 .and. present(weighting_name)) then
1430      write(nulout,'(a,a,a,f6.0,a,f6.0,a)') 'Spectral weights for ', &
1431           &  weighting_name, ' (', wavenumber1, ' to ', &
1432           &  wavenumber2, ' cm-1):'
1433      do jband = 1, nweights
1434        write(nulout, '(a,i0,a,f6.0,a,f6.0,a,f8.4)') '  Shortwave band ', &
1435             &  iband(jband), ' (', this%wavenumber1_sw(iband(jband)), ' to ', &
1436             &  this%wavenumber2_sw(iband(jband)), ' cm-1): ', weight(jband)
1437      end do
1438    end if
1439
1440  end subroutine get_sw_weights
1441
1442
1443  !---------------------------------------------------------------------
1444  ! The input shortwave surface albedo coming in is likely to be in
1445  ! different spectral intervals to the gas model in the radiation
1446  ! scheme. We assume that the input albedo is defined within
1447  ! "ninterval" spectral intervals covering the wavelength range 0 to
1448  ! infinity, but allow for the possibility that two intervals may be
1449  ! indexed back to the same albedo band. 
1450  subroutine define_sw_albedo_intervals(this, ninterval, wavelength_bound, &
1451       &                                i_intervals, do_nearest)
1452
1453    use radiation_io, only : nulerr, radiation_abort
1454
1455    class(config_type),   intent(inout) :: this
1456    ! Number of spectral intervals in which albedo is defined
1457    integer,              intent(in)    :: ninterval
1458    ! Monotonically increasing wavelength bounds between intervals,
1459    ! not including the outer bounds (which are assumed to be zero and
1460    ! infinity)
1461    real(jprb),           intent(in)    :: wavelength_bound(ninterval-1)
1462    ! The albedo indices corresponding to each interval
1463    integer,              intent(in)    :: i_intervals(ninterval)
1464    logical,    optional, intent(in)    :: do_nearest
1465   
1466    if (ninterval > NMaxAlbedoIntervals) then
1467      write(nulerr,'(a,i0,a,i0)') '*** Error: ', ninterval, &
1468           &  ' albedo intervals exceeds maximum of ', NMaxAlbedoIntervals
1469      call radiation_abort();
1470    end if
1471
1472    if (present(do_nearest)) then
1473      this%do_nearest_spectral_sw_albedo = do_nearest
1474    else
1475      this%do_nearest_spectral_sw_albedo = .false.
1476    end if
1477    this%sw_albedo_wavelength_bound(1:ninterval-1) = wavelength_bound(1:ninterval-1)
1478    this%sw_albedo_wavelength_bound(ninterval:)    = -1.0_jprb
1479    this%i_sw_albedo_index(1:ninterval)            = i_intervals(1:ninterval)
1480    this%i_sw_albedo_index(ninterval+1:)           = 0
1481
1482    if (this%is_consolidated) then
1483      call this%consolidate_intervals(.true., &
1484           &  this%do_nearest_spectral_sw_albedo, &
1485           &  this%sw_albedo_wavelength_bound, this%i_sw_albedo_index, &
1486           &  this%wavenumber1_sw, this%wavenumber2_sw, &
1487           &  this%i_albedo_from_band_sw, this%sw_albedo_weights)
1488    end if
1489
1490  end subroutine define_sw_albedo_intervals
1491
1492
1493  !---------------------------------------------------------------------
1494  ! As define_sw_albedo_intervals but for longwave emissivity
1495  subroutine define_lw_emiss_intervals(this, ninterval, wavelength_bound, &
1496       &                                i_intervals, do_nearest)
1497
1498    use radiation_io, only : nulerr, radiation_abort
1499
1500    class(config_type),   intent(inout) :: this
1501    ! Number of spectral intervals in which emissivity is defined
1502    integer,              intent(in)    :: ninterval
1503    ! Monotonically increasing wavelength bounds between intervals,
1504    ! not including the outer bounds (which are assumed to be zero and
1505    ! infinity)
1506    real(jprb),           intent(in)    :: wavelength_bound(ninterval-1)
1507    ! The emissivity indices corresponding to each interval
1508    integer,              intent(in)    :: i_intervals(ninterval)
1509    logical,    optional, intent(in)    :: do_nearest
1510   
1511    if (ninterval > NMaxAlbedoIntervals) then
1512      write(nulerr,'(a,i0,a,i0)') '*** Error: ', ninterval, &
1513           &  ' emissivity intervals exceeds maximum of ', NMaxAlbedoIntervals
1514      call radiation_abort();
1515    end if
1516
1517    if (present(do_nearest)) then
1518      this%do_nearest_spectral_lw_emiss = do_nearest
1519    else
1520      this%do_nearest_spectral_lw_emiss = .false.
1521    end if
1522    this%lw_emiss_wavelength_bound(1:ninterval-1) = wavelength_bound(1:ninterval-1)
1523    this%lw_emiss_wavelength_bound(ninterval:)    = -1.0_jprb
1524    this%i_lw_emiss_index(1:ninterval)            = i_intervals(1:ninterval)
1525    this%i_lw_emiss_index(ninterval+1:)           = 0
1526
1527    if (this%is_consolidated) then
1528      call this%consolidate_intervals(.false., &
1529           &  this%do_nearest_spectral_lw_emiss, &
1530           &  this%lw_emiss_wavelength_bound, this%i_lw_emiss_index, &
1531           &  this%wavenumber1_lw, this%wavenumber2_lw, &
1532           &  this%i_emiss_from_band_lw, this%lw_emiss_weights)
1533    end if
1534
1535  end subroutine define_lw_emiss_intervals
1536
1537
1538  !---------------------------------------------------------------------
1539  ! This routine consolidates either the input shortwave albedo
1540  ! intervals with the shortwave bands, or the input longwave
1541  ! emissivity intervals with the longwave bands, depending on the
1542  ! arguments provided.
1543  subroutine consolidate_intervals(this, is_sw, do_nearest, &
1544       &  wavelength_bound, i_intervals, wavenumber1, wavenumber2, &
1545       &  i_mapping, weights)
1546
1547    use radiation_io, only : nulout, nulerr, radiation_abort
1548
1549    class(config_type),   intent(inout) :: this
1550    ! Is this the shortwave?  Otherwise longwave
1551    logical,    intent(in)    :: is_sw
1552    ! Do we find the nearest albedo interval to the centre of each
1553    ! band, or properly weight the contributions? This can be modified
1554    ! if there is only one albedo intervals.
1555    logical, intent(inout)    :: do_nearest
1556    ! Monotonically increasing wavelength bounds between intervals,
1557    ! not including the outer bounds (which are assumed to be zero and
1558    ! infinity)
1559    real(jprb), intent(in)    :: wavelength_bound(NMaxAlbedoIntervals-1)
1560    ! The albedo band indices corresponding to each interval
1561    integer,    intent(in)    :: i_intervals(NMaxAlbedoIntervals)
1562    ! Start and end wavenumber bounds for the ecRad bands (cm-1)
1563    real(jprb), intent(in)    :: wavenumber1(:), wavenumber2(:)
1564
1565    ! if do_nearest is TRUE then the result is expressed in i_mapping,
1566    ! which will be allocated to have the same length as wavenumber1,
1567    ! and contain the index of the albedo interval corresponding to
1568    ! that band
1569    integer,    allocatable, intent(inout) :: i_mapping(:)
1570    ! ...otherwise the result is expressed in "weights", of
1571    ! size(n_intervals, n_bands) containing how much of each interval
1572    ! contributes to each band.
1573    real(jprb), allocatable, intent(inout) :: weights(:,:)
1574
1575    ! Number and loop index of ecRad bands
1576    integer    :: nband, jband
1577    ! Number and index of albedo/emissivity intervals
1578    integer    :: ninterval, iinterval
1579    ! Sometimes an albedo or emissivity value will be used in more
1580    ! than one interval, so nvalue indicates how many values will
1581    ! actually be provided
1582    integer    :: nvalue
1583    ! Wavenumber bounds of the albedo/emissivity interval
1584    real(jprb) :: wavenumber1_albedo, wavenumber2_albedo
1585    ! Reciprocal of the wavenumber range of the ecRad band
1586    real(jprb) :: recip_dwavenumber ! cm
1587    ! Midpoint/bound of wavenumber band
1588    real(jprb) :: wavenumber_mid, wavenumber_bound ! cm-1
1589   
1590    nband = size(wavenumber1)
1591
1592    ! Count the number of albedo/emissivity intervals
1593    ninterval = 0
1594    do iinterval = 1,NMaxAlbedoIntervals
1595      if (i_intervals(iinterval) > 0) then
1596        ninterval = iinterval
1597      else
1598        exit
1599      end if
1600    end do
1601
1602    if (ninterval < 2) then
1603      ! Zero or one albedo/emissivity intervals found, so we index all
1604      ! bands to one interval
1605      if (allocated(i_mapping)) then
1606        deallocate(i_mapping)
1607      end if
1608      allocate(i_mapping(nband))
1609      i_mapping(:) = 1
1610      do_nearest = .true.
1611      ninterval = 1
1612      nvalue = 1
1613    else
1614      ! Check wavelength is monotonically increasing
1615      do jband = 2,ninterval-1
1616        if (wavelength_bound(jband) <= wavelength_bound(jband-1)) then
1617          if (is_sw) then
1618            write(nulerr, '(a,a)') '*** Error: wavelength bounds for shortwave albedo intervals ', &
1619                 &  'must be monotonically increasing'
1620          else
1621            write(nulerr, '(a,a)') '*** Error: wavelength bounds for longwave emissivity intervals ', &
1622                 &  'must be monotonically increasing'
1623          end if
1624          call radiation_abort()
1625        end if
1626      end do
1627
1628      ! What is the maximum index, indicating the number of
1629      ! albedo/emissivity values to expect?
1630      nvalue = maxval(i_intervals(1:ninterval))
1631     
1632      if (do_nearest) then
1633        ! Simpler nearest-neighbour mapping from band to
1634        ! albedo/emissivity interval
1635        if (allocated(i_mapping)) then
1636          deallocate(i_mapping)
1637        end if
1638        allocate(i_mapping(nband))
1639
1640        ! Loop over bands
1641        do jband = 1,nband
1642          ! Compute mid-point of band in wavenumber space (cm-1)
1643          wavenumber_mid = 0.5_jprb * (wavenumber1(jband) &
1644               &                     + wavenumber2(jband))
1645          iinterval = 1
1646          ! Convert wavelength (m) into wavenumber (cm-1) at the lower
1647          ! bound of the albedo interval
1648          wavenumber_bound = 0.01_jprb / wavelength_bound(iinterval)
1649          ! Find the albedo interval that has the largest overlap with
1650          ! the band; this approach assumes that the albedo intervals
1651          ! are larger than the spectral bands
1652          do while (wavenumber_bound >= wavenumber_mid &
1653               &    .and. iinterval < ninterval)
1654            iinterval = iinterval + 1
1655            if (iinterval < ninterval) then
1656              wavenumber_bound = 0.01_jprb / wavelength_bound(iinterval)
1657            else
1658              ! For the last interval there is no lower bound
1659              wavenumber_bound = 0.0_jprb
1660            end if
1661          end do
1662          ! Save the index of the band corresponding to the albedo
1663          ! interval and move onto the next band
1664          i_mapping(jband) = i_intervals(iinterval)
1665        end do
1666      else
1667        ! More accurate weighting
1668        if (allocated(weights)) then
1669          deallocate(weights)
1670        end if
1671        allocate(weights(nvalue,nband))
1672        weights(:,:) = 0.0_jprb
1673       
1674        ! Loop over bands
1675        do jband = 1,nband
1676          recip_dwavenumber = 1.0_jprb / (wavenumber2(jband) &
1677               &                        - wavenumber1(jband))
1678          ! Find the first overlapping albedo band
1679          iinterval = 1
1680          ! Convert wavelength (m) into wavenumber (cm-1) at the lower
1681          ! bound (in wavenumber space) of the albedo/emissivty interval
1682          wavenumber1_albedo = 0.01_jprb / wavelength_bound(iinterval)
1683          do while (wavenumber1_albedo >= wavenumber2(jband) &
1684               &    .and. iinterval < ninterval)
1685            iinterval = iinterval + 1
1686            wavenumber1_albedo = 0.01_jprb / wavelength_bound(iinterval)
1687          end do
1688         
1689          wavenumber2_albedo = wavenumber2(jband)
1690         
1691          ! Add all overlapping bands
1692          do while (wavenumber2_albedo > wavenumber1(jband) &
1693               &  .and. iinterval <= ninterval)
1694            weights(i_intervals(iinterval),jband) &
1695                 &  = weights(i_intervals(iinterval),jband) &
1696                 &  + recip_dwavenumber &
1697                 &  * (min(wavenumber2_albedo,wavenumber2(jband)) &
1698                 &   - max(wavenumber1_albedo,wavenumber1(jband)))
1699            wavenumber2_albedo = wavenumber1_albedo
1700            iinterval = iinterval + 1
1701            if (iinterval < ninterval) then
1702              wavenumber1_albedo = 0.01_jprb / wavelength_bound(iinterval)
1703            else
1704              wavenumber1_albedo = 0.0_jprb
1705            end if
1706          end do
1707        end do
1708      end if
1709    end if
1710
1711    ! Define how many bands to use for reporting surface downwelling
1712    ! fluxes for canopy radiation scheme
1713    if (is_sw) then
1714      if (this%use_canopy_full_spectrum_sw) then
1715        this%n_canopy_bands_sw = this%n_g_sw
1716      else
1717        this%n_canopy_bands_sw = nvalue
1718      end if
1719    else
1720      if (this%use_canopy_full_spectrum_lw) then
1721        this%n_canopy_bands_lw = this%n_g_lw
1722      else
1723        this%n_canopy_bands_lw = nvalue
1724      end if
1725    end if
1726
1727    if (this%iverbosesetup >= 2) then
1728      if (.not. do_nearest) then
1729        if (is_sw) then
1730          write(nulout, '(a,i0,a,i0,a)') 'Weighting of ', nvalue, ' albedo values in ', &
1731             &  nband, ' shortwave bands (wavenumber ranges in cm-1):'
1732        else
1733          write(nulout, '(a,i0,a,i0,a)') 'Weighting of ', nvalue, ' emissivity values in ', &
1734             &  nband, ' longwave bands (wavenumber ranges in cm-1):'
1735        end if
1736        do jband = 1,nband
1737          write(nulout,'(i6,a,i6,a)',advance='no') nint(wavenumber1(jband)), ' to', &
1738               &                        nint(wavenumber2(jband)), ':'
1739          do iinterval = 1,nvalue
1740            write(nulout,'(f5.2)',advance='no') weights(iinterval,jband)
1741          end do
1742          write(nulout, '()')
1743        end do
1744      else if (ninterval <= 1) then
1745        if (is_sw) then
1746          write(nulout, '(a)') 'All shortwave bands will use the same albedo'
1747        else
1748          write(nulout, '(a)') 'All longwave bands will use the same emissivty'
1749        end if
1750      else
1751        if (is_sw) then
1752          write(nulout, '(a,i0,a)',advance='no') 'Mapping from ', nband, &
1753               &  ' shortwave bands to albedo intervals:'
1754        else
1755          write(nulout, '(a,i0,a)',advance='no') 'Mapping from ', nband, &
1756               &  ' longwave bands to emissivity intervals:'
1757        end if
1758        do jband = 1,nband
1759          write(nulout,'(a,i0)',advance='no') ' ', i_mapping(jband)
1760        end do
1761        write(nulout, '()')
1762      end if
1763    end if
1764
1765  end subroutine consolidate_intervals
1766
1767
1768  !---------------------------------------------------------------------
1769  ! Return the 0-based index for str in enum_str, or abort if it is
1770  ! not found
1771  subroutine get_enum_code(str, enum_str, var_name, icode)
1772
1773    use radiation_io, only : nulerr, radiation_abort
1774
1775    character(len=*), intent(in)  :: str
1776    character(len=*), intent(in)  :: enum_str(0:)
1777    character(len=*), intent(in)  :: var_name
1778    integer,          intent(out) :: icode
1779
1780    integer :: jc
1781    logical :: is_not_found
1782
1783    ! If string is empty then we don't modify icode but assume it has
1784    ! a sensible default value
1785    if (len_trim(str) > 1) then
1786      is_not_found = .true.
1787
1788      do jc = 0,size(enum_str)-1
1789        if (trim(str) == trim(enum_str(jc))) then
1790          icode = jc
1791          is_not_found = .false.
1792          exit
1793        end if
1794      end do
1795      if (is_not_found) then
1796        write(nulerr,'(a,a,a,a,a)',advance='no') '*** Error: ', trim(var_name), &
1797             &  ' must be one of: "', enum_str(0), '"'
1798        do jc = 1,size(enum_str)-1
1799          write(nulerr,'(a,a,a)',advance='no') ', "', trim(enum_str(jc)), '"'
1800        end do
1801        write(nulerr,'(a)') ''
1802        call radiation_abort('Radiation configuration error')
1803      end if
1804    end if
1805
1806  end subroutine get_enum_code
1807
1808
1809  !---------------------------------------------------------------------
1810  ! Print one line of information: logical
1811  subroutine print_logical(message_str, name, val)
1812    use radiation_io, only : nulout
1813    character(len=*),   intent(in) :: message_str
1814    character(len=*),   intent(in) :: name
1815    logical,            intent(in) :: val
1816    character(4)                   :: on_or_off
1817    character(NPrintStringLen)     :: str
1818    if (val) then
1819      on_or_off = ' ON '
1820    else
1821      on_or_off = ' OFF'
1822    end if
1823    write(str, '(a,a4)') message_str, on_or_off
1824    write(nulout,'(a,a,a,a,l1,a)') str, ' (', name, '=', val,')'
1825  end subroutine print_logical
1826
1827
1828  !---------------------------------------------------------------------
1829  ! Print one line of information: integer
1830  subroutine print_integer(message_str, name, val)
1831    use radiation_io, only : nulout
1832    character(len=*),   intent(in) :: message_str
1833    character(len=*),   intent(in) :: name
1834    integer,            intent(in) :: val
1835    character(NPrintStringLen)     :: str
1836    write(str, '(a,a,i0)') message_str, ' = ', val
1837    write(nulout,'(a,a,a,a)') str, ' (', name, ')'
1838  end subroutine print_integer
1839
1840
1841  !---------------------------------------------------------------------
1842  ! Print one line of information: real
1843  subroutine print_real(message_str, name, val)
1844    use parkind1,     only : jprb
1845    use radiation_io, only : nulout
1846    character(len=*),   intent(in) :: message_str
1847    character(len=*),   intent(in) :: name
1848    real(jprb),         intent(in) :: val
1849    character(NPrintStringLen)     :: str
1850    write(str, '(a,a,g8.3)') message_str, ' = ', val
1851    write(nulout,'(a,a,a,a)') str, ' (', name, ')'
1852  end subroutine print_real
1853
1854
1855  !---------------------------------------------------------------------
1856  ! Print one line of information: enum
1857  subroutine print_enum(message_str, enum_str, name, val)
1858    use radiation_io, only : nulout
1859    character(len=*),   intent(in) :: message_str
1860    character(len=*),   intent(in) :: enum_str(0:)
1861    character(len=*),   intent(in) :: name
1862    integer,            intent(in) :: val
1863    character(NPrintStringLen)     :: str
1864    write(str, '(a,a,a,a)') message_str, ' "', trim(enum_str(val)), '"'
1865    write(nulout,'(a,a,a,a,i0,a)') str, ' (', name, '=', val,')'
1866  end subroutine print_enum
1867
1868
1869  !---------------------------------------------------------------------
1870  ! Return .true. if 1D allocatable array "var" is out of physical
1871  ! range specified by boundmin and boundmax, and issue a warning.
1872  ! "do_fix" determines whether erroneous values are fixed to lie
1873  ! within the physical range. To check only a subset of the array,
1874  ! specify i1 and i2 for the range.
1875  function out_of_bounds_1d(var, var_name, boundmin, boundmax, do_fix, i1, i2) result (is_bad)
1876
1877    use radiation_io,     only : nulout
1878
1879    real(jprb), allocatable, intent(inout) :: var(:)
1880    character(len=*),        intent(in) :: var_name
1881    real(jprb),              intent(in) :: boundmin, boundmax
1882    logical,                 intent(in) :: do_fix
1883    integer,       optional, intent(in) :: i1, i2
1884
1885    logical                       :: is_bad
1886
1887    real(jprb) :: varmin, varmax
1888
1889    is_bad = .false.
1890
1891    if (allocated(var)) then
1892
1893      if (present(i1) .and. present(i2)) then
1894        varmin = minval(var(i1:i2))
1895        varmax = maxval(var(i1:i2))
1896      else
1897        varmin = minval(var)
1898        varmax = maxval(var)
1899      end if
1900
1901      if (varmin < boundmin .or. varmax > boundmax) then
1902        write(nulout,'(a,a,a,g12.4,a,g12.4,a,g12.4,a,g12.4)',advance='no') &
1903             &  '*** Warning: ', var_name, ' range', varmin, ' to', varmax, &
1904             &  ' is out of physical range', boundmin, 'to', boundmax
1905        is_bad = .true.
1906        if (do_fix) then
1907          if (present(i1) .and. present(i2)) then
1908            var(i1:i2) = max(boundmin, min(boundmax, var(i1:i2)))
1909          else
1910            var = max(boundmin, min(boundmax, var))
1911          end if
1912          write(nulout,'(a)') ': corrected'
1913        else
1914          write(nulout,'(1x)')
1915        end if
1916      end if
1917
1918    end if
1919   
1920  end function out_of_bounds_1d
1921
1922
1923  !---------------------------------------------------------------------
1924  ! Return .true. if 2D allocatable array "var" is out of physical
1925  ! range specified by boundmin and boundmax, and issue a warning.  To
1926  ! check only a subset of the array, specify i1 and i2 for the range
1927  ! of the first dimension and j1 and j2 for the range of the second.
1928  function out_of_bounds_2d(var, var_name, boundmin, boundmax, do_fix, &
1929       &                    i1, i2, j1, j2) result (is_bad)
1930
1931    use radiation_io,     only : nulout
1932
1933    real(jprb), allocatable, intent(inout) :: var(:,:)
1934    character(len=*),        intent(in) :: var_name
1935    real(jprb),              intent(in) :: boundmin, boundmax
1936    logical,                 intent(in) :: do_fix
1937    integer,       optional, intent(in) :: i1, i2, j1, j2
1938
1939    ! Local copies of indices
1940    integer :: ii1, ii2, jj1, jj2
1941
1942    logical                       :: is_bad
1943
1944    real(jprb) :: varmin, varmax
1945
1946    is_bad = .false.
1947
1948    if (allocated(var)) then
1949
1950      if (present(i1) .and. present(i2)) then
1951        ii1 = i1
1952        ii2 = i2
1953      else
1954        ii1 = lbound(var,1)
1955        ii2 = ubound(var,1)
1956      end if
1957      if (present(j1) .and. present(j2)) then
1958        jj1 = j1
1959        jj2 = j2
1960      else
1961        jj1 = lbound(var,2)
1962        jj2 = ubound(var,2)
1963      end if
1964      varmin = minval(var(ii1:ii2,jj1:jj2))
1965      varmax = maxval(var(ii1:ii2,jj1:jj2))
1966
1967      if (varmin < boundmin .or. varmax > boundmax) then
1968        write(nulout,'(a,a,a,g12.4,a,g12.4,a,g12.4,a,g12.4)',advance='no') &
1969             &  '*** Warning: ', var_name, ' range', varmin, ' to', varmax,&
1970             &  ' is out of physical range', boundmin, 'to', boundmax
1971        is_bad = .true.
1972        if (do_fix) then
1973          var(ii1:ii2,jj1:jj2) = max(boundmin, min(boundmax, var(ii1:ii2,jj1:jj2)))
1974          write(nulout,'(a)') ': corrected'
1975        else
1976          write(nulout,'(1x)')
1977        end if
1978      end if
1979
1980    end if
1981   
1982  end function out_of_bounds_2d
1983
1984
1985  !---------------------------------------------------------------------
1986  ! Return .true. if 3D allocatable array "var" is out of physical
1987  ! range specified by boundmin and boundmax, and issue a warning.  To
1988  ! check only a subset of the array, specify i1 and i2 for the range
1989  ! of the first dimension, j1 and j2 for the second and k1 and k2 for
1990  ! the third.
1991  function out_of_bounds_3d(var, var_name, boundmin, boundmax, do_fix, &
1992       &                    i1, i2, j1, j2, k1, k2) result (is_bad)
1993
1994    use radiation_io,     only : nulout
1995
1996    real(jprb), allocatable, intent(inout) :: var(:,:,:)
1997    character(len=*),        intent(in) :: var_name
1998    real(jprb),              intent(in) :: boundmin, boundmax
1999    logical,                 intent(in) :: do_fix
2000    integer,       optional, intent(in) :: i1, i2, j1, j2, k1, k2
2001
2002    ! Local copies of indices
2003    integer :: ii1, ii2, jj1, jj2, kk1, kk2
2004
2005    logical                       :: is_bad
2006
2007    real(jprb) :: varmin, varmax
2008
2009    is_bad = .false.
2010
2011    if (allocated(var)) then
2012
2013      if (present(i1) .and. present(i2)) then
2014        ii1 = i1
2015        ii2 = i2
2016      else
2017        ii1 = lbound(var,1)
2018        ii2 = ubound(var,1)
2019      end if
2020      if (present(j1) .and. present(j2)) then
2021        jj1 = j1
2022        jj2 = j2
2023      else
2024        jj1 = lbound(var,2)
2025        jj2 = ubound(var,2)
2026      end if
2027      if (present(k1) .and. present(k2)) then
2028        kk1 = k1
2029        kk2 = k2
2030      else
2031        kk1 = lbound(var,3)
2032        kk2 = ubound(var,3)
2033      end if
2034      varmin = minval(var(ii1:ii2,jj1:jj2,kk1:kk2))
2035      varmax = maxval(var(ii1:ii2,jj1:jj2,kk1:kk2))
2036
2037      if (varmin < boundmin .or. varmax > boundmax) then
2038        write(nulout,'(a,a,a,g12.4,a,g12.4,a,g12.4,a,g12.4)',advance='no') &
2039             &  '*** Warning: ', var_name, ' range', varmin, ' to', varmax,&
2040             &  ' is out of physical range', boundmin, 'to', boundmax
2041        is_bad = .true.
2042        if (do_fix) then
2043          var(ii1:ii2,jj1:jj2,kk1:kk2) = max(boundmin, min(boundmax, &
2044               &                             var(ii1:ii2,jj1:jj2,kk1:kk2)))
2045          write(nulout,'(a)') ': corrected'
2046        else
2047          write(nulout,'(1x)')
2048        end if
2049      end if
2050
2051    end if
2052   
2053  end function out_of_bounds_3d
2054
2055
2056end module radiation_config
Note: See TracBrowser for help on using the repository browser.