source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_config.F90 @ 4773

Last change on this file since 4773 was 4773, checked in by idelkadi, 5 months ago
  • Update of Ecrad in LMDZ The same organization of the Ecrad offline version is retained in order to facilitate the updating of Ecrad in LMDZ and the comparison between online and offline results. version 1.6.1 of Ecrad (https://github.com/lguez/ecrad.git)
  • Implementation of the double call of Ecrad in LMDZ


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