-
Notifications
You must be signed in to change notification settings - Fork 68
RTE: optical properties
Class ty_optical_props
in module mo_optical_props
in the RTE library provides a representation of the optical properties of the atmosphere. The following examples assume that a variable op
of type(ty_optical_props)
has been declared.
Optical properties are defined by their spectral dependence which is described during initialization using routine
function init(band_lims_wvn, band_lims_gpt, name) result(err_message)
real(wp), dimension(:,:), intent(in ) :: band_lims_wvn
integer, dimension(:,:), &
optional, intent(in ) :: band_lims_gpt
character(len=*), optional, intent(in ) :: name
character(len = 128) :: err_message
Array band_lims_wvn
with extents (2, number-of-bands) describes the number of bands beginning and ending wavenumber of each band (in MKS units i.e. inverse meters). Optional argument band_lims_gpt
describes the beginning and ending g-point of each band; if this array isn't provided it's assumed values are available by band. An optional name
may be useful in debugging.
Variables of type ty_optical_props
can copy the spectral discretization from a previously-initialized variable e.g. err_msg = op2%init(op1)
.
Real-valued functions op%get_band_lims_wavenumber()
and op%get_band_lims_wavelength()
return arrays of extent (2, number-of-bands) containing the spectral limits of each band expressed in wavelength or wavenumber = 1/wavelength (again in MKS units).
Procedure op%is_initialized()
returns a logical value indicating if the spectral discretization has been provided. op%get_nband()
and op%get_ngpt()
return the (integer) number of band and g-points; if these numbers are the same the optical properties are defined by band.
Integer functions op%convert_band2gpt()
and op%convert_gpt2band()
provide a map between individual bands and g-points. Function op%get_gpoint_bands()
returns an integer of length op%get_ngpt()
with the band number to which each g-point belongs (i.e. it's the vector version of op%convert_gpt2band()
).
Comparison functions op1%bands_are_equal(op2)
and op1%gpoints_are_equal(op2)
return logical values indicating whether the twi sets of optical properties share the same band or g-point discretization.
Class ty_optical_props_arry
allows optical properties to be expressed as arrays: discrete values on a domain with dimensions column, layer, g-points. Class ty_optical_props_arry
extends class ty_optical_props
and inherits all the capabilities described above.
Class ty_optical_props_arry
is abstract. Three concrete sub-classes extend this type:
-
ty_optical_props_1scl
contains the data fieldtau
with dimensionsncol, nlay, ngpt
representing optical depth. This class is useful for describing problems without scattering. -
ty_optical_props_2str
contains the data fieldstau, ssa, g
with dimensionsncol, nlay, ngpt
representing optical depth, single-scattering albedo, and asymmetry parameter. This class is used to describe scattering problems for which the two-stream approximation is adequate. -
ty_optical_props_nstr
contains the data fieldstau, ssa, p
. Variablep
has leading dimensionnmom
to accommodate the moments of a phase function for use in higher-order calculations.
Variables of class ty_optical_props_arry
must be allocated as well as initialized before use. Because the argument list is slightly different among the sub-classes each allocation routine has a distinct name. Users must call one of err_msg = op%alloc_1scl(ncol, nlay)
, err_msg = op%alloc_2str(ncol, nlay)
or err_msg = op%alloc_nstr(nmom, ncol, nlay)
.
Initialization from an existing variable of type ty_optical_props
can be combined with allocation: err_msg = op%alloc_2str(ncol, nlay, op1)
, etc.
Integer functions op%get_ncol()
and op%get_nlay()
return the number of columns or layers in the domain. Class ty_optical_props_nstr
has a get_nmom()
function.
Calling err_msg = op%delta_scale(for)
will delta-scale the optical properties to remove strong forward peaks in the phase function. Optional argument for
is an real array of forward-scattering fractions with the same extents as op%tau
. If this argument is omitted the assumption for = g*g
is made.
Optical properties can be added together: err_msg = op1%increment(op2)
adds the optical properties of op1
to those of op2
making reasonable assumptions if the two variables are not of the same class.
err_msg = op%validate()
checks to be sure that all values are valid: optical depth not less than 0, single-scattering albedo between 0 and 1, asymmetry parameter between -1 and 1.
err_msg = op%subset(start, n)
extracts columns start
to start+n-1
from the class into variable subset
making reasonable assumptions if the variables are of different classes.