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

Last change on this file since 5260 was 4853, checked in by idelkadi, 8 months ago

Ecrad update in LMDZ / Implementation of Ecrad double call in LMDZ

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