HEAD PREVIOUS

3  XRTM C Interface

XRTM is unlike other common RT models in that it exists as an object/instance that is created, modified, and destroyed. It maintains a valid state between simulations where only inputs that change need to be updated for subsequent simulations. This design allows input overhead to be minimized and allows redundant calculations to be saved in a transparent manor. Each XRTM instance is contained within an isolated memory scope so that multiple instances may be used in a shared memory multiprocessing environment.
The interface is made up of input configuration constants (section 3.2), functions for creating and destroying an XRTM instance (section 3.3), functions for setting and getting inputs (section 3.4), and functions that run the appropriate calculations and return outputs (section 3.5). Typical use of XRTM would be to create an instance with xrtm_create(), set inputs with the xrtm_set_*() functions, run the model and get outputs with the xrtm_calc_*() functions, loop over the last two steps until the model is no longer needed, and then finally, destroy the instance with xrtm_destroy().

3.1  (Arrays of (arrays of)) ... arrays

Several of the XRTM input and output functions take multi-dimensional arrays as arguments. In C89/90 if true multi-dimensional arrays are passed to functions all but the first dimension must be static. This is an obvious limitation for the XRTM C interface as these sizes are in fact dynamic. So instead, XRTM uses "(arrays of pointers to (arrays of pointers to)) ... 1-d arrays" or as they are referred to in this manual "(arrays of (arrays of)) ... 1-d arrays". These structures may be allocated by the user but as a convenience XRTM provides functions to do this. These functions are efficient in that internally they allocate memory with one call for all the arrays making up the entire structure and then set the pointers to the appropriate locations in memory. This method also has the advantage that the data lies contiguously in memory. So for example, with a 2 dimensional array the rows would lie one after a another where as if the rows were allocated separately they might not lie contiguously in memory.
The following is a formal description of the functions for allocating "(arrays of (arrays of)) ... 1-d arrays" and their corresponding deallocation functions:

<type name> *...*alloc_array<# of dimensions>_<type id>(int size_1,...,int size_n)

Description:
Allocate an "(array of (arrays of)) ... 1-d arrays". <# of dimensions> is either 1, 2, 3, etc. and <type> is i or d for int or double arrays, respectively. The arguments are the sizes of each dimension of the array (size_1, ... size_n) where n is the # of dimensions. The return value is a "(pointer to (a pointer to)) ... a pointer" to <type name> with a pointer depth equal to the # of dimensions.
Arguments:
size_1    size of the first dimension
size_n    size of the last (nth) dimension
Return value:
A "(pointer to (a pointer to)) ... a pointer" to <type name> represting the allocated "(array of (arrays of)) ... 1-d arrays" or NULL on error.

void free_array<# of dimensions>_<type id>(<type name> *...*array)

Description:
Deallocate (free) an "(array of (arrays of)) ... 1-d arrays". <# of dimensions> is either 1, 2, 3, etc. and <type> is i or d for int or double arrays, respectively. The argument is the "(pointer to (a pointer to)) ... a pointer" to <type name> represting the "(array of (arrays of)) ... 1-d arrays" returned by alloc_array<# of dimensions>_<type id>.
Arguments:
array    "(pointer to (a pointer to)) ... a pointer" to <type name> represting the "(array of (arrays of)) ... 1-d arrays"
Return value:
None.
So for example, if a 2 dimensional array of arrays of type double is to be allocated then the function to call would be

double **alloc_array2_d(int size_1, int size_2)

and the corresponding deallocation routine would be

void free_array2_d(double **array)

3.2  Configuration constants

3.2.1  Options

Options are turned on by setting the appropriate bit of a 32 bit wide mask which is the options argument to xrtm_create(). The appropriate bits may be set using masks (declared as enumeration constants) associated with each option. For example, Delta-M scaling and the pseudo spherical approximation may be turned on by using the bitwise inclusive OR operator with something like
options = XRTM_OPTION_DELTA_M | XRTM_OPTION_PSA.
XRTM_OPTION_CALC_DERIVS

Calculate derivatives with respect to optical property inputs. Requires n_derivs to be greater than or equal to one.
XRTM_OPTION_DELTA_M

Use Delta-M scaling [].
XRTM_OPTION_N_T_TMS

Use the Nakajima and Tanaka TMS correction [].
XRTM_OPTION_FOUR_CONV_OLD

Used for testing purposes only.
XRTM_OPTION_FOUR_CONV_NEW

Used for testing purposes only.
XRTM_OPTION_NO_AZIMUTHAL

Include only the first Fourier term of the expansion in azimuth, i.e. the azimuthal average.
XRTM_OPTION_OUTPUT_AT_LEVELS

Output at user specified levels. This is in contrast to output at user specified optical depths. Some solvers support output at TOA only, others support output at TOA and/or BOA, while some support output at any level. Check the solver descriptions for which solver supports what. Requires at least one call to xrtm_set_out_levels() once the model is created.
XRTM_OPTION_OUTPUT_AT_TAUS

Output at user specified optical depths from TOA. This is in contrast to output at user specified levels. Some solvers support output at TOA only, others support output at TOA and/or BOA, other support output at any level, while some support output at any optical depth (within layers). Check the solver descriptions for which solver supports what. Requires at least one call to xrtm_set_out_taus() once the model is created.
XRTM_OPTION_PHASE_SCALAR

Used for testing purposes only.
XRTM_OPTION_PHASE_MATRIX_GC

Used for testing purposes only.
XRTM_OPTION_PHASE_MATRIX_LC

Used for testing purposes only.
XRTM_OPTION_PSA

Use the so called pseudo-spherical approximation to model the solar beam through a spherical spherical shell atmosphere [].
XRTM_OPTION_QUAD_NORM_GAUS_LEG

Use (standard) Gauss-Legendre quadrature.
XRTM_OPTION_QUAD_DOUB_GAUS_LEG

Use double Gauss-Legendre quadrature.
XRTM_OPTION_QUAD_LOBATTO

Use Lobatto quadrature.
XRTM_OPTION_SAVE_PHASE_MATS

Save phase matrices between XRTM calls.
XRTM_OPTION_SAVE_LOCAL_R_T

Save local r and t matrices between XRTM calls.
XRTM_OPTION_SAVE_LAYER_R_T_S

XRTM_OPTION_SAVE_TOTAL_R_T_S

XRTM_OPTION_SFI

For solvers that support it use source function integration for output at arbitrary zenith angles otherwise quadrature dummy nodes are used.
XRTM_OPTION_SOURCE_SOLAR

Include solar sources.
XRTM_OPTION_SOURCE_THERMAL

Include thermal sources.
XRTM_OPTION_STACK_REUSE_ADDING

XRTM_OPTION_TOP_DOWN_ADDING

For solvers that use adding add from the top down to get output at the surface only. This is in contrast to full adding and can significantly improve run time.
XRTM_OPTION_BOTTOM_UP_ADDING

For solvers that use adding add from the bottom up to get output at TOA only. This is in contrast to full adding and can significantly improve run time.
XRTM_OPTION_UPWELLING_OUTPUT

Output upwelling values.
XRTM_OPTION_DOWNWELLING_OUTPUT

Output downwelling values.
XRTM_OPTION_VECTOR

Run the model in vector mode.

3.2.2  Solvers

XRTM may be created to use any number of solvers. Limiting the initialization to only the solvers required will in some cases lead to significant memory savings. Solvers are turned on by setting the appropriate bit of a 32 bit wide mask which is the solvers argument to xrtm_create(). The appropriate bits may be set using masks (declared as enumeration constants) associated with each solver. For example, XRTM may be create to use the Eigenmatrix/BVP solver along with the single and second order scattering solvers by using the bitwise inclusive OR operator with something like
solvers = XRTM_SOLVER_EIG_BVP | XRTM_SOLVER_SINGLE | XRTM_SOLVER_TWO_OS.
XRTM_SOLVER_DOUB_ADD

Doubling/Adding: Use doubling to get global reflection and transmission matrices for each layer. Then use adding to get global reflection and transmission matrices for the entire atmosphere [Grant and Hunt, 1969,de Haan et al., 1987,Liou, 2002].
XRTM_SOLVER_EIG_ADD

Eigenmatrix/Adding: Use the Eigenvalue problem to get global reflection and transmission matrices for each layer. Then use adding to get global reflection and transmission matrices for the entire atmosphere [Aronson, 1972,Nakajima and Tanaka, 1986,Voronovich et al., 2004,Spurr and Christi, 2007].
XRTM_SOLVER_EIG_BVP

Eigenmatrix/BVP (a.k.a. the Discrete Ordinate Method): Use the Eigenvalue problem to obtain the layer homogeneous solution. Then solve a boundary value problem for the entire atmosphere [Liou, 1973,Stamnes et al., 1988,Siewert, 2000,Spurr et al., 2001].
XRTM_SOLVER_MEM_BVP

Matrix exponential by eigenmatrix / BVP: A variant of the Discrete Ordinate Method with a matrix exponential formulation [Doicu and Trautmann, 2009a,Doicu and Trautmann, 2009b].
XRTM_SOLVER_PADE_ADD

Matrix exponential by Padé approximation / Adding (A.K.A. PARTM): Use the Padé approximation to the matrix exponential to get global reflection and transmission matrices for each layer. Then use adding to get global reflection and transmission matrices for the entire atmosphere [McGarragh and Gabriel, 2010,].
XRTM_SOLVER_SINGLE

Includes only first order scattering from the atmosphere and surface.
XRTM_SOLVER_SOS

Successive orders of scattering using an approximate integration in optical thickness [Fymat and Ueno, 1974,Min and Duan, 2004,Lenoble et al., 2007].
XRTM_SOLVER_TWO_OS

Second order scattering with the typical numerical integration over zenith and azimuth but with an analytical integration in optical thickness. [Kawabata and Ueno, 1988,Natraj and Spurr, 2007]

3.2.3  BRDF kernels

[]
XRTM_KERNEL_LAMBERTIAN

XRTM_KERNEL_ROUJEAN

[Roujean et al., 1992]
XRTM_KERNEL_LI_SPARSE

[Wanner et al., 1995]
XRTM_KERNEL_LI_DENSE

[]
XRTM_KERNEL_ROSS_THIN

[]
XRTM_KERNEL_ROSS_THICK

[]
XRTM_KERNEL_HAPKE

[Hapke, 1981,Hapke and Wells, 1981]
XRTM_KERNEL_RAHMAN

[Rahman et al., 1993]
XRTM_KERNEL_COX_MUNK

[Cox and Munk, 1954]

3.2.4  Solutions

XRTM_OUTPUT_RADIANCE

XRTM_OUTPUT_RADIANCE_MEAN

XRTM_OUTPUT_FLUX

XRTM_OUTPUT_FLUX_DIVERGENCE

3.3  Initiating and destroying XRTM

int xrtm_create(xrtm_data *d, int options, int solvers, int max_coef, int n_quad, int n_stokes, int n_derivs, int n_layers, int n_kernels, int n_kernel_quad, int *kernels, int n_out_levels, int n_out_thetas)

Description:
Create a new XRTM instance. When finished with the instance created, xrtm_destroy() must be called to free memory allocated by xrtm_create().
Arguments:
d    the xrtm_data structure which will represent the instance created
options    bit mask of XRTM configuration options
solvers    bit mask of XRTM solvers that will be used
max_coef    maximum number of phase function Legendre expansion coefficients that will be used
n_quad    number of quadrature points in one hemisphere
n_stokes    size of the stokes vector to calculate (set to one for scalar mode)
n_derivs    number of derivatives to calculate
n_layers    number of plane parallel layers in the atmosphere
n_kernels    number of BRDF kernels to use for the BRDF
n_kernel_quad    number of quadrature points to use for BRDF integration
kernels    array of BRDF kernels to use (of length n_kernels)
n_out_levels    number of user defined output levels
n_out_thetas    number of user defined output zenith angles
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

void xrtm_destroy(xrtm_data *d)

Description:
Destroy an XRTM instance which includes freeing all memory allocated by xrtm_create().
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
None.

3.4  Setting and getting inputs

int xrtm_get_options(xrtm_data *d)

Description:
Get the bit mask of XRTM options with which this instance was created.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The bit mask of XRTM options with which this instance was created or XRTM_INT_ERROR on error.

int xrtm_get_solvers(xrtm_data *d)

Description:
Get the bit mask of XRTM solvers for which this instance has been created to use.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The bit mask of XRTM solvers for which this instance has been created to use or XRTM_INT_ERROR on error.

int xrtm_get_max_coef(xrtm_data *d)

Description:
Get the maximum number of phase function Legendre expansion coefficients for which this XRTM instance has been created to handle.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The maximum number of phase function Legendre expansion coefficients for which this XRTM instance has been created to handle or XRTM_INT_ERROR on error.

int xrtm_get_n_quad(xrtm_data *d)

Description:
Get the number of quadrature points in one hemisphere used by this XRTM instance.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The number of quadrature points in one hemisphere used by this XRTM instance or XRTM_INT_ERROR on error.

int xrtm_get_n_stokes(xrtm_data *d)

Description:
Get the size of the stokes vector for which this instance has been created to calculate.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The size of the stokes vector for which this instance has been created to calculate or XRTM_INT_ERROR on error.

int xrtm_get_n_derivs(xrtm_data *d)

Description:
Get the number of derivatives for which this instance has been created to calculate.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The number of derivatives for which this instance has been created to calculate or XRTM_INT_ERROR on error.

int xrtm_get_n_layers(xrtm_data *d)

Description:
Get the number of plane parallel layers in the atmosphere modeled by this instance.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The number of plane parallel layers in the atmosphere modeled by this instance or XRTM_INT_ERROR on error.

int xrtm_get_n_kernels(xrtm_data *d)

Description:
Get the number of BRDF kernels for which this instance has been created to use.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The number of BRDF kernels for which this instance has been created to use or XRTM_INT_ERROR on error.

int xrtm_get_n_kernel_quad(xrtm_data *d)

Description:
Get the number of quadrature points for BRDF integration used by this instance.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The number of quadrature points for BRDF integration used by this instance or XRTM_INT_ERROR on error.

int xrtm_get_kernel(xrtm_data *d, int i_kernel)

Description:
Get the kernel identifier for a given kernel index. The kernel index is the index at which the kernel was given in the array kernels given as input to xrtm_create().
Arguments:
d    the xrtm_data structure which represents the instance created
i_kernel    the kernel index, where 0 ≤ i_kerneln_kernels − 1
Return value:
The kernel identifier for index i_kernel or XRTM_INT_ERROR on error.

int xrtm_get_n_out_levels(xrtm_data *d)

Description:
Get the number of levels at which this instance has been created to output.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The number of levels at which this instance has been created to output or XRTM_INT_ERROR on error.

int xrtm_get_n_out_thetas(xrtm_data *d)

Description:
Get the number of zenith angles at which this instance has been created to output.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The number of zenith angles at which this instance has been created to output or XRTM_INT_ERROR on error.

int xrtm_set_doub_d_tau(xrtm_data *d, double d_tau)

Description:
Set the initial layer optical thickness ∆τ for the doubling method. To set this value XRTM must have been created to use at least one of the following solvers: XRTM_SOLVER_DOUB_ADD, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
d_tau    the initial layer optical thickness ∆τ, where d_tau > 0.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_doub_d_tau(xrtm_data *d)

Description:
Get the initial layer optical thickness ∆τ for the doubling method. To get this value XRTM must have been created to use at least one of the following solvers: XRTM_SOLVER_DOUB_ADD, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The initial layer optical thickness ∆τ for the doubling method or XRTM_DBL_ERROR on error.

int xrtm_set_pade_params(xrtm_data *d, int pade_s, int pade_r)

Description:
Set the Padé scaling power of two (number of doublings) s and the degree of Padé approximate r. If either value is set to a value that is out of range then s and r are chosen automatically from a lookup table based on layer optical thickness τ and the maximum output zenith angle θ. To set these values, XRTM must have been created to use at least one of the following solvers: XRTM_SOLVER_PADE_ADD, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
pade_s    Padé scaling power of two s, where pade_s ≥ 0
pade_r    degree of Padé approximate r, where pade_r > 0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_get_pade_params(xrtm_data *d, int *pade_s, int *pade_r)

Description:
Get the Padé scaling power of two (number of doublings) s and the degree of Padé approximate r. To get these values, XRTM must have been created to use at least one of the following solvers: XRTM_SOLVER_PADE_ADD, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
pade_s    (output) Padé scaling power of two s
pade_r    (output) degree of Padé approximate r
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_set_sos_params(xrtm_data *d, int max_os, double max_tau, double sos_tol)

Description:
Set parameters related to successive order of scattering. They are the maximum order of scattering that will be computed, the maximum layer optical thickness used (all layers of a larger optical thickness are divided evenly into enough sub-layers so that each sub-layer has an optical thickness less than or equal to the maximum allowable value), and the successive order of scattering tolerance limit. The tolerance limit is the minimum radiance contribution from any single quadrature angle with which the succession will continue to the next order of scattering. If this limit is not met the succession will terminate. To set these values, XRTM must have been created to use at least one of the following solvers: XRTM_SOLVER_SOS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
max_os    maximum order of scattering, where max_os ≥ 0
max_tau    maximum layer optical thickness used, where max_tau > 0.0
sos_tol    successive order of scattering tolerance limit, where sos_tol ≥ 0.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_get_sos_params(xrtm_data *d, int *max_os, double *max_tau, double *sos_tol)

Description:
Get parameters related to successive order of scattering. To get these values, XRTM must have been created to use at least one of the following solvers: XRTM_SOLVER_SOS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
max_os    (output) maximum order of scattering
max_tau    (output) maximum layer optical thickness used
sos_tol    (output) successive order of scattering tolerance limit
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_set_fourier_tol(xrtm_data *d, double fourier_tol)

Description:
Set the tolerance limit for the Fourier expansion in azimuth angle. The tolerance limit is the minimum intensity contribution from any single output level and output angle with which the summation will continue to the next term. If this limit is not met for all of the output levels and output angles the summation will terminate. If a single scattering correction is to be applied (XRTM_OPTION_N_T_TMS) then the series starts with the full (untruncated) single scattering contribution while each term includes only the truncated multiple scattering contribution. If it is set to zero then the summation will include all terms (2n).
Arguments:
d    the xrtm_data structure which represents the instance created
fourier_tol    tolerance limit for the Fourier expansion in azimuth angle, where fourier_tol ≥ 0.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_fourier_tol(xrtm_data *d)

Description:
Get the tolerance limit for the Fourier expansion in azimuth angle.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The tolerance limit for the Fourier expansion in azimuth angle or XRTM_DBL_ERROR on error.

int xrtm_set_F_0(xrtm_data *d, double F_0)

Description:
Set the intensity of the incident parallel beam at TOA F0. Setting F0 to zero turns off the solar source. If a thermal source is used (XRTM_OPTION_SOURCE_THERMAL) then the units for F0 and the TOA, BOA, and level Planck radiances must be the same. Otherwise the units for F0 are arbitrary.
Arguments:
d    the xrtm_data structure which represents the instance created
F_0    intensity of the incident parallel beam at TOA F0, where F_0 ≥ 0.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_F_0(xrtm_data *d)

Description:
Get the intensity of the incident parallel beam at TOA F0.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The intensity of the incident parallel beam at TOA F0 or XRTM_DBL_ERROR on error.

int xrtm_set_theta_0(xrtm_data *d, double theta_0)

Description:
Set the zenith angle for the incident parallel beam at TOA (the solar zenith angle) θ0.
Arguments:
d    the xrtm_data structure which represents the instance created
theta_0    zenith angle for the incident parallel beam at TOA θ0, where 0.0 ≤ theta_0 < 90.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_theta_0(xrtm_data *d)

Description:
Get the zenith angle for the incident parallel beam at TOA (the solar zenith angle) θ0.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The zenith angle for incident parallel beam at TOA θ0 or XRTM_DBL_ERROR on error.

int xrtm_set_phi_0(xrtm_data *d, double phi_0)

Description:
Set the azimuth angle for the incident parallel beam at TOA (the solar azimuth angle) ϕ0.
Arguments:
d    the xrtm_data structure which represents the instance created
phi_0    azimuth angle for the incident parallel beam at TOA ϕ0, where 0.0 ≤ phi_0 < 360.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_phi_0(xrtm_data *d)

Description:
Get the azimuth angle for the incident parallel beam at TOA (the solar azimuth angle) ϕ0.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The azimuth angle for incident parallel beam at TOA ϕ0 or XRTM_DBL_ERROR on error.

int xrtm_set_out_levels(xrtm_data *d, int *out_levels)

Description:
Set the levels at which to output results given an array of length n_out_levels. Levels are defined at layer boundaries. For example a 3 layer atmosphere would have 4 levels and TOA and BOA (the surface) would be levels 0 and 3, respectively. Levels must be specified in ascending order. To set this value XRTM must have been created with the option XRTM_OPTION_OUTPUT_AT_LEVELS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
out_levels    array of output level indices in ascending order, where 0 ≤ out_levels[i]n_layers and 0 ≤ in_out_levels − 1
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_get_out_levels(xrtm_data *d, int *out_levels)

Description:
Get the levels at which results are output at as an array of length n_out_levels. To get this value XRTM must have been created with the option XRTM_OPTION_OUTPUT_AT_LEVELS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
out_levels    (output) array of output level indices in ascending order of length n_out_levels
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_set_out_taus(xrtm_data *d, double *out_taus)

Description:
Set the optical depths at which to output results given an array of length n_out_levels. Optical depths must be specified in ascending order. To set this value XRTM must have been created with the option XRTM_OPTION_OUTPUT_AT_TAUS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
out_taus    array of output optical depths in ascending order, where 0.0 ≤ out_taus[i] ≤ τs and 0 ≤ in_out_levels − 1 and τs is the optical depth to the surface
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_get_out_taus(xrtm_data *d, double *out_taus)

Description:
Get the optical depths at which results are output at as an array of length n_out_levels. To get this value XRTM must have been created with the option XRTM_OPTION_OUTPUT_AT_TAUS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
out_taus    (output) array of output optical depths in ascending order of length n_out_levels
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_set_out_thetas(xrtm_data *d, double *out_thetas)

Description:
Set the zenith angles θ at which to output results given an array of length n_out_thetas. To set this value XRTM must have been created with n_out_thetas > 0, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
out_thetas    array of output zenith angles θ, where 0.0 ≤ out_thetas[i] < 90.0 and 0 ≤ in_out_thetas − 1
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_get_out_thetas(xrtm_data *d, double *out_thetas)

Description:
Get the zenith angles θ at which results are output at as an array of length n_out_thetas. To set this value XRTM must have been created with n_out_thetas > 0, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
out_thetas    (output) array of output zenith angles θ of length n_out_thetas
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_set_top_b(xrtm_data *d, double top_b)

Description:
Set the intensity of the downward isotopic radiation at TOA. Must be in the same units as F0 and the BOA and level Planck radiances. To set this value XRTM must have been created with the option XRTM_OPTION_SOURCE_THERMAL, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
top_b    intensity of the downward isotopic radiation at TOA, where top_b ≥ 0.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_top_b(xrtm_data *d)

Description:
Get the intensity of the downward isotopic radiation at TOA. To get this value XRTM must have been created with the option XRTM_OPTION_SOURCE_THERMAL, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The intensity of the downward isotopic radiation at TOA or XRTM_DBL_ERROR on error.

int xrtm_set_planet_r(xrtm_data *d, double planet_r)

Description:
Set the planetary radius for the point located at BOA. The units for this value and the level heights must be the same. To set this value XRTM must have been created with the option XRTM_OPTION_PSA, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
planet_r    planetary radius, where planet_r ≥ 0.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_planet_r(xrtm_data *d)

Description:
Get the planetary radius for the point located at BOA. To get this value XRTM must have been created with the option XRTM_OPTION_PSA, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The planetary radius for the point located at BOA or XRTM_DBL_ERROR on error.

int xrtm_set_levels_z(xrtm_data *d, double *levels_z)

Description:
Set the level heights z as an array of length n_layers + 1. The units for this value and the planetary radius must be the same. To set this value XRTM must have been created with the option XRTM_OPTION_PSA, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
levels_z    array of level heights z, where levels_z[i] ≥ 0.0 and 0 ≤ in_layers
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_get_levels_z(xrtm_data *d, double *levels_z)

Description:
Get the level heights z as an array of length n_layers + 1. To get this value XRTM must have been created with the option XRTM_OPTION_PSA, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
levels_z    (output) array of level heights z of length n_layers + 1
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_set_levels_b(xrtm_data *d, double *levels_b)

Description:
Set the level Planck radiances B as an array of length n_layers + 1. Must be in the same units as F0 and the TOA and BOA Planck radiances. To set this value XRTM must have been created with the option XRTM_OPTION_SOURCE_THERMAL, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
levels_b    array of level Planck radiances B, where levels_b[i] ≥ 0.0 and 0 ≤ in_layers
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_get_levels_b(xrtm_data *d, double *levels_b)

Description:
Get the level Planck radiances B as an array of length n_layers + 1. To get this value XRTM must have been created with the option XRTM_OPTION_SOURCE_THERMAL, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
levels_b    (output) array of level Planck radiances B of length n_layers + 1
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_set_coef_1(xrtm_data *d, int i_layer, int n_coef, double **coef)

int xrtm_set_coef_n(xrtm_data *d, int *n_coef, double ***coef)

Description:
Set layer phase matrix expansion coefficients Bij,l, where i,j = 1 for scalar transport (n_stokes = 1) and 1 ≤ i,j ≤ 4 for vector transport (when option XRTM_OPTION_VECTOR has been specified and n_stokes > 1). In the case of scalar transport a single element B0,0,l is given to XRTM which are the Legendre expansion coefficients βl in the equation for the scattering phase function given by
P(cosΘ) = N

l = 0 
βl Pl(cosΘ),
where Θ is single scattering angle, Pl(cosΘ) are Legendre polynomials of degree l, and N is the maximum degree term in the expansion. In the case of vector transport 6 elements are given to XRTM, the so "Greek constants" [Siewert, 1982] of the matrix
Bl =





βl
γl
0
0
γl
αl
0
0
0
0
ζl
−εl
0
0
εl
δl






,
which are the coefficients for expansion in terms of generalized spherical functions of the phase matrix
F(Θ) =





a1(Θ)
b1(Θ)
0
0
b1(Θ)
a2(Θ)
0
0
0
0
a3(Θ)
b2(Θ)
0
0
−b2(Θ)
a4(Θ)






,
where
a1(Θ)
=
N

l = 0 
βl Pl0,0(cosΘ),
a2(Θ) + a3(Θ)
=
N

l = 2 
l + ζl) Pl2,2(cosΘ),
a3(Θ) − a3(Θ)
=
N

l = 2 
l − ζl) Pl2,−2(cosΘ),
a4(Θ)
=
N

l = 0 
δl Pl0,0(cosΘ),
b1(Θ)
=
N

l = 0 
γl Pl0,2(cosΘ),
b2(Θ)
=
N

l = 0 
εl Pl0,2(cosΘ),
and Plm,n(cosΘ) are generalized spherical functions. It is worth noting that a1(Θ) is equal to the scattering phase function P(cosΘ) of the scalar case and the generalized spherical function Pl0,0(cosΘ) is equal to the Legendre polynomial Pl(cosΘ). The expansion coefficients are given to XRTM as a (6 × n_coef) array of arrays according to the mapping
coef[0][l]
=
βl
coef[1][l]
=
αl
coef[2][l]
=
ζl
coef[3][l]
=
δl
coef[4][l]
=
−γl
coef[5][l]
=
−εl,
where as explanined below, n_coef is the number of coefficients, coef is the input array, and 0 ≤ ln_coef.
Any number of phase matrix coeffieicents may be supplied by the user for each layer with the restriction that the number must be less than or equal to max_coef which is set when an XRTM instance is created. (It is up to the user to set max_coef to the maximum number of coefficients to be supplied for all layers.) It is important that the user is aware of the optimal choice in different cases. With XRTM_OPTION_DELTA_M turned off 2 * n_quad − 1 coeffieicents will be used so if your phase function has less than or equal to as many coeffieicents then supply them all, otherwise you only need to supply a maximum of 2 * n_quad − 1. If XRTM_OPTION_DELTA_M is turned on then an additional coefficient must be supplied for a total of 2 * n_quad. If your phase function has less than as many coeffieicents then it does not need to be delta-M scaled and will not be affected by delta-M scaling. If XRTM_OPTION_N_T_TMS is turned on (which requires XRTM_OPTION_DELTA_M to be turned on) then the all the coefficients required to represent the full phase matrix should be supplied.
A choice of two functions are provided for this purpose of setting these elements. xrtm_set_coef_1() sets the values for a given layer index i_layer given an (n_elem × n_coef) array of arrays, where n_coef is the number of coefficients given for the layer. xrtm_set_coef_n() sets the values for all layers given an (n_layers × n_elem × n_coef[i], 0 ≤ in_layers − 1) array of arrays of arrays where n_coef is an array of length n_layers giving the number of coefficients for each layer. In both cases n_elem = 1 for scalar mode and n_elem = 6 when option XRTM_OPTION_VECTOR has been specified.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to set, where 0 ≤ i_layern_layers − 1
n_coef    number of coefficients given as a scalar (_1) or an array (_n) of length n_layers, depending on the function called, where 0 ≤ n_coef, n_coef[i]max_coef
coef    phase matrix expansion coefficients βij as an array of arrays (_1) or an array of arrays of arrays (_1), depending on the function called, where −∞ ≤ coef[i][j],coef[i][j][k] ≤ ∞, except that coef[0][0],coef[i][0][0] = 1
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_get_n_coef(xrtm_data *d, int i_layer)

Description:
Get the number of phase matrix expansion coefficients for layer index i_layer.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to get, where 0 ≤ i_layern_layers − 1
Return value:
The number of phase matrix expansion coefficients for layer index i_layer or XRTM_INT_ERROR on error.

double xrtm_get_coef(xrtm_data *d, int i_layer, int i_elem, int i_coef)

Description:
Get the phase matrix expansion coefficient for layer index i_layer, matrix element index i_elem, and coefficient index i_coef.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to get, where 0 ≤ i_layern_layers − 1
i_elem    index of phase matrix element to get, where 0 ≤ i_elemn_elem and n_elem = 0 for scalar mode and n_elem = 6 when option XRTM_OPTION_VECTOR has been specified
i_coef    index of phase matrix coefficient to get, where 0 ≤ i_coefn_coef[i] − 1 and 0 ≤ in_layers − 1
Return value:
The phase matrix expansion coefficient (i_elem, i_coef) for layer index i_layer or XRTM_DBL_ERROR on error.

int xrtm_set_coef_l_11(xrtm_data *d, int i_layer, int i_deriv, double **coef_l)

int xrtm_set_coef_l_n1(xrtm_data *d, int i_deriv, double ***coef_l)

int xrtm_set_coef_l_1n(xrtm_data *d, int i_layer, double ***coef_l)

int xrtm_set_coef_l_nn(xrtm_data *d, double ****coef_l)

Description:
Set layer linearized phase matrix expansion coefficients ∂βij / ∂x. A choice of four functions are provided for this purpose. xrtm_set_coef_l_11() sets the values for a given layer index i_layer and a given derivative index i_deriv given an (n_elem × n_coef) array of arrays. xrtm_set_coef_l_n1() sets the values for all layers for a given derivative index i_deriv given an (n_layers × n_elem × n_coef[i], 0 ≤ in_layers − 1) array of arrays of arrays. xrtm_set_coef_l_1n() sets the values for all derivatives for a given layer index i_layer given an (n_derivs × n_elem × n_coef) array of arrays of arrays. Finally, xrtm_set_coef_l_nn() sets the values for all layers and all derivatives given an (n_layers × n_derivs × n_elem × n_coef[i], 0 ≤ in_layers − 1) array of arrays of arrays of arrays. In all cases n_elem = 1 for scalar mode and n_elem = 6 when option XRTM_OPTION_VECTOR has been specified and n_coef is defined for each layer by xrtm_set_coef(). To set these values XRTM must have been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to set, where 0 ≤ i_layern_layers − 1
i_deriv    index of derivative to set, where 0 ≤ i_derivn_derivs − 1
coef_l    linearized phase matrix expansion coefficients ∂βij / ∂x as an array of arrays (_11), an array of arrays of arrays (_n1 or _1n), or an array of arrays of arrays of arrays (_nn), depending on the function called, where −∞ ≤ coef_l[i][j],coef_l[i][j][k],coef_l[i][j][k][l] ≤ ∞, except that coef_l[0][0],coef_l[i][0][0],coef_l[i][j][0][0] = 0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_coef_l(xrtm_data *d, int i_layer, int i_deriv, int i_elem, int i_coef)

Description:
Get the linearized phase matrix expansion coefficient for layer index i_layer, derivative index i_deriv, matrix element index i_elem, and coefficient index i_coef. To get these values XRTM must have been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to get, where 0 ≤ i_layern_layers − 1
i_deriv    index of derivative to get, where 0 ≤ i_derivn_derivs − 1
i_elem    index of phase matrix element to get, where 0 ≤ i_elemn_elem and n_elem = 0 for scalar mode and n_elem = 6 when option XRTM_OPTION_VECTOR has been specified
i_coef    index of phase matrix coefficient to get, where 0 ≤ i_coefn_coef[i] − 1 and 0 ≤ in_layers − 1
Return value:
The linearized phase matrix expansion coefficient (i_elem, i_coef) for layer index i_layer and derivative index i_deriv or XRTM_DBL_ERROR on error.

int xrtm_set_omega_1(xrtm_data *d, int i_layer, double omega)

int xrtm_set_omega_n(xrtm_data *d, double *omega)

Description:
Set layer single scattering albedo ω. A choice of two functions are provided for this purpose. xrtm_set_omega_1() sets the value for a given layer index i_layer given a scalar value. xrtm_set_omega_n() sets the values for all layers given an array of length n_layers.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to set, where 0 ≤ i_layern_layers − 1
omega    single scattering albedo ω as a scalar (_1) or an array (_n), depending on the function called, where 0.0 ≤ omega,omega[i] ≤ 1.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_omega(xrtm_data *d, int i_layer)

Description:
Get single scattering albedo ω for layer index i_layer.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to get, where 0 ≤ i_layern_layers − 1
Return value:
The single scattering albedo ω for layer index i_layer or XRTM_DBL_ERROR on error.

int xrtm_set_omega_l_11(xrtm_data *d, int i_layer, int i_deriv, double omega_l)

int xrtm_set_omega_l_n1(xrtm_data *d, int i_deriv, double *omega_l)

int xrtm_set_omega_l_1n(xrtm_data *d, int i_layer, double *omega_l)

int xrtm_set_omega_l_nn(xrtm_data *d, double **omega_l)

Description:
Set layer linearized single scattering albedo ∂ω/ ∂x. A choice of four functions are provided for this purpose. xrtm_set_omega_l_11() sets the value for a given layer index i_layer and a given derivative index i_deriv given a scalar value. xrtm_set_omega_l_n1() sets the values for all layers for a given derivative index i_deriv given an array of length n_layers. xrtm_set_omega_l_1n() sets the values for all derivatives for a given layer index i_layer given an array of length n_derivs. Finally, xrtm_set_omega_l_nn() sets the values for all layers and all derivatives given an (n_layers × n_derivs) array of arrays. To set these values XRTM must have been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to set, where 0 ≤ i_layern_layers − 1
i_deriv    index of derivative to set, where 0 ≤ i_derivn_derivs − 1
omega_l    linearized single scattering albedo ∂ω/ ∂x as a scalar (_1), an array (_n1 or _1n), or an array of arrays (_nn), depending on the function called, where −∞ ≤ omega_l,omega_l[i],omega_l[i][j] ≤ ∞
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_omega_l(xrtm_data *d, int i_layer, int i_deriv)

Description:
Get linearized single scattering albedo ∂ω/ ∂x for layer index i_layer and derivative index i_deriv. To get these values XRTM must have been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to get, where 0 ≤ i_layern_layers − 1
i_deriv    index of derivative to get, where 0 ≤ i_derivn_derivs − 1
Return value:
The linearized single scattering albedo ∂ω/ ∂x for layer index i_layer and derivative index i_deriv or XRTM_DBL_ERROR on error.

int xrtm_set_ltau_1(xrtm_data *d, int i_layer, double ltau)

int xrtm_set_ltau_n(xrtm_data *d, double *ltau)

Description:
Set layer optical thickness τ. A choice of two functions are provided for this purpose. xrtm_set_ltau_1() sets the value for a given layer index i_layer given a scalar value. xrtm_set_ltau_n() sets the values for all layers given an array of length n_layers.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to set, where 0 ≤ i_layern_layers − 1
ltau    optical thickness τ as a scalar (_1) or an array (_n), depending on the function called, where ltau,ltau[i] ≥ 0.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_ltau(xrtm_data *d, int i_layer)

Description:
Get optical thickness τ for layer index i_layer.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to get, where 0 ≤ i_layern_layers − 1
Return value:
The optical thickness τ for layer index i_layer or XRTM_DBL_ERROR on error.

int xrtm_set_ltau_l_11(xrtm_data *d, int i_layer, int i_deriv, double ltau_l)

int xrtm_set_ltau_l_n1(xrtm_data *d, int i_deriv, double *ltau_l)

int xrtm_set_ltau_l_1n(xrtm_data *d, int i_layer, double *ltau_l)

int xrtm_set_ltau_l_nn(xrtm_data *d, double **ltau_l)

Description:
Set linearized layer optical thickness ∂τ/ ∂x. A choice of four functions are provided for this purpose. xrtm_set_ltau_l_11() sets the value for a given layer index i_layer and a given derivative index i_deriv given a scalar value. xrtm_set_ltau_l_n1() sets the values for all layers for a given derivative index i_deriv given an array of length n_layers. xrtm_set_ltau_l_1n() sets the values for all derivatives for a given layer index i_layer given an array of length n_derivs. Finally, xrtm_set_ltau_l_nn() sets the values for all layers and all derivatives given an (n_layers × n_derivs) array of arrays. To set these values XRTM must have been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to set, where 0 ≤ i_layern_layers − 1
i_deriv    index of derivative to set, where 0 ≤ i_derivn_derivs − 1
ltau_l    linearized optical thickness ∂τ/ ∂x as a scalar (_11), an array (_n1 and _1n), or an array of arrays (_nn), depending on the function called, where −∞ ≤ ltau_l, ltau_l[i],ltau_l[i][j] ≤ ∞
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_ltau_l(xrtm_data *d, int i_layer, int i_deriv)

Description:
Get linearized optical thickness ∂τ/ ∂x for layer index i_layer and derivative index i_deriv. To get these values XRTM must have been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_layer    index of layer to get, where 0 ≤ i_layern_layers − 1
i_deriv    index of derivative to get, where 0 ≤ i_derivn_derivs − 1
Return value:
The linearized optical thickness ∂τ/ ∂x for layer index i_layer and derivative index i_deriv or XRTM_DBL_ERROR on error.

int xrtm_set_surface_b(xrtm_data *d, double surface_b)

Description:
Set the surface Planck radiance. Must be in the same units as F0 and the TOA and level Planck radiances. To set this value XRTM must have been created with the option XRTM_OPTION_SOURCE_THERMAL, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
surface_b    surface Planck radiance, where surface_b ≥ 0.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_surface_b(xrtm_data *d)

Description:
Get the surface Planck radiance. To get this value XRTM must have been created with the option XRTM_OPTION_SOURCE_THERMAL, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
The surface Planck radiance or XRTM_DBL_ERROR on error.

int xrtm_set_kernel_ampfac(xrtm_data *d, int i_kernel, double ampfac)

Description:
Set the amplitude factor for kernel index i_kernel. The kernel index is the index at which the kernel was given in the array kernels given as input to xrtm_create(). To set this value XRTM must have been created with n_kernels > 0, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_kernel    index of kernel to set, where 0 ≤ i_kerneln_kernels − 1
ampfac    the amplitude factor for kernel i_kernel, where 0.0 ≤ ampfac ≤ 1.0
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_kernel_ampfac(xrtm_data *d, int i_kernel)

Description:
Get the amplitude factor for kernel index i_kernel. The kernel index is the index at which the kernel was given in the array kernels given as input to xrtm_create().
Arguments:
d    the xrtm_data structure which represents the instance created
i_kernel    index of kernel to get, where 0 ≤ i_kerneln_kernels − 1
Return value:
The amplitude factor for kernel index i_kernel or XRTM_DBL_ERROR on error.

int xrtm_set_kernel_params_1(xrtm_data *d, int i_kernel, int i_param, double param)

int xrtm_set_kernel_params_n(xrtm_data *d, int i_kernel, double *params)

Description:
Set the kernel parameters for kernel index i_kernel. The kernel index is the index at which the kernel was given in the array kernels given as input to xrtm_create(). A choice of two functions are provided for this purpose. xrtm_set_kernel_params_1() sets the value for a given parameter index i_param given a scalar value. xrtm_set_kernel_params_n() sets the values for all parameters given an array of length n_params, where n_params is the number of parameters required for kernel i_kernel. The number of parameters and a description of each parameter is given for each kernel in section 3.2.3.
Arguments:
d    the xrtm_data structure which represents the instance created
i_kernel    index of kernel to set, where 0 ≤ i_kerneln_kernels − 1
i_param    index of parameter to set, where 0 ≤ i_paramn_params − 1
param    kernel parameter as a scalar (_1) or an array (_n), depending on the function called, where −∞ ≤ param,param[i] ≤ ∞
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_kernel_params(xrtm_data *d, int i_kernel, int i_param)

Description:
Get the kernel parameter for kernel index i_kernel and parameter index i_param. The kernel index is the index at which the kernel was given in the array kernels given as input to xrtm_create(). The parameter indices are described for each kernel in section 3.2.3.
Arguments:
d    the xrtm_data structure which represents the instance created
i_kernel    index of kernel to get, where 0 ≤ i_kerneln_kernels − 1
i_param    index of parameter to set, where 0 ≤ i_paramn_params − 1
Return value:
The kernel parameter for kernel index i_kernel and parameter index i_param or XRTM_DBL_ERROR on error.

int xrtm_set_kernel_ampfac_l_1(xrtm_data *d, int i_kernel, int i_deriv, double ampfac_l)

int xrtm_set_kernel_ampfac_l_n(xrtm_data *d, int i_kernel, double *ampfac_l)

Description:
Set the linearized amplitude factor for kernel index i_kernel. The kernel index is the index at which the kernel was given in the array kernels given as input to xrtm_create(). A choice of two functions are provided for this purpose. xrtm_set_kernel_ampfac_l_1() sets the value for a given derivative index i_deriv given a scalar value. xrtm_set_kernel_ampfac_l_n() sets the values for all derivatives given an array of length n_derivs. To set this value XRTM must have been created with n_kernels > 0, otherwise it is an error. To set these values XRTM must have been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_kernel    index of kernel to set, where 0 ≤ i_kerneln_kernels − 1
i_deriv    index of derivative to set, where 0 ≤ i_derivn_derivs − 1
ampfac_l    the linearized amplitude factor for kernel i_kernel, where −∞ ≤ ampfac,ampfac[i] ≤ ∞
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_kernel_ampfac_l(xrtm_data *d, int i_kernel, int i_deriv)

Description:
Get the linearized amplitude factor for kernel index i_kernel and derivative index i_deriv . The kernel index is the index at which the kernel was given in the array kernels given as input to xrtm_create(). To get these values XRTM must have been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_kernel    index of kernel to get, where 0 ≤ i_kerneln_kernels − 1
i_deriv    index of derivative to get, where 0 ≤ i_derivn_derivs − 1
Return value:
The linearized amplitude factor for kernel index i_kernel or XRTM_DBL_ERROR on error.

int xrtm_set_kernel_params_l_11(xrtm_data *d, int i_kernel, int i_deriv, int i_param, double param_l)

int xrtm_set_kernel_params_l_1n(xrtm_data *d, int i_kernel, int i_deriv, double *params_l)

int xrtm_set_kernel_params_l_n1(xrtm_data *d, int i_kernel, int i_param, double *params_l)

int xrtm_set_kernel_params_l_nn(xrtm_data *d, int i_kernel, double **params_l)

Description:
Set the linearized kernel parameters for kernel index i_kernel. The kernel index is the index at which the kernel was given in the array kernels given as input to xrtm_create(). A choice of four functions are provided for this purpose. xrtm_set_kernel_params_l_11() sets the value for a given derivative index i_deriv and a given parameter index i_param given a scalar value. xrtm_set_kernel_params_l_1n() sets the value for all parameters given a derivative index i_deriv given an array of length n_params. xrtm_set_kernel_params_l_n1() sets the value for all derivatives given a parameter index i_param given an array of length n_derivs. xrtm_set_kernel_params_l_nn() sets the value for all derivatives given and all parameters given a (n_derivs × n_params) array of arrays. n_params is the number of parameters required for kernel i_kernel. The number of parameters and a description of each parameter is given for each kernel in section 3.2.3. To set these values XRTM must have been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_kernel    index of kernel to set, where 0 ≤ i_kerneln_kernels − 1
i_deriv    index of derivative to set, where 0 ≤ i_derivn_derivs − 1
i_param    index of parameter to set, where 0 ≤ i_paramn_params − 1
param_l    the linearized kernel parameter for kernel i_kernel, where −∞ ≤ param_l,param_l[i],param_l[i][j] ≤ ∞
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

double xrtm_get_kernel_params_l(xrtm_data *d, int i_kernel, int i_deriv, int i_param)

Description:
Get the linearized kernel parameter for kernel index i_kernel, derivative index i_deriv, and parameter index i_param. The kernel index is the index at which the kernel was given in the array kernels given as input to xrtm_create(). The parameter indices are described for each kernel in section 3.2.3. To get these values XRTM must have been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
i_kernel    index of kernel to get, where 0 ≤ i_kerneln_kernels − 1
i_deriv    index of derivative to get, where 0 ≤ i_derivn_derivs − 1
i_param    index of parameter to set, where 0 ≤ i_paramn_params − 1
Return value:
The linearized kernel parameter for kernel index i_kernel and parameter index i_param or XRTM_DBL_ERROR on error.

3.5  Running the model and getting output

int xrtm_update_varied_layers(xrtm_data *d)

Description:
This function must be called every time the set of layers that vary, including the surface, changes. A layer "varies" if at least one of the linearized values associated with that layer is nonzero. For example, if at least one of the linearized values for a particular layer has been set to a nonzero value where before all the values for that layer were zero then that layer has been added to the set of layers that varies and xrtm_update_varied_layers() must be called. Equivalently, If all the linearized values for a layer are set to zero where before at least one of the values was nonzero then that layer has been removed from the set of layers that vary so xrtm_update_varied_layers() must be called. On the other hand, if the values of nonzero linearized values are only changed to other nonzero values then xrtm_update_varied_layers() does not need to be called. For efficiency, it is a good idea to call the function once, after all changes to all layers and all values have been changed just before running the model. Calling this function is only valid when XRTM has been created with the option XRTM_CALC_DERIVS, otherwise it is an error.
Arguments:
d    the xrtm_data structure which represents the instance created
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_solution(xrtm_data *d, enum xrtm_solver_mask solver, int solutions, int n_out_phis, double **out_phis, double ****I_p, double ****I_m, double *****I_p_l, double *****I_m_l, double *mean_p, double *mean_m, double **mean_p_l, double **mean_m_l, double *flux_p, double *flux_m, double **flux_p_l, double **flux_m_l, double *flux_div, double **flux_div_l)

Description:
Run XRTM using a specified solver and return various result types.
Arguments:
d    the xrtm_data structure which represents the instance created
solver    bit mask for the solver to be used
solutions    bit mask for the solutions to return results for
n_out_phis    number azimuth angles ϕ to return results for
out_phis    array of arrays of output azimuth angles with dimensions (n_out_thetas, n_out_phis), where 0.0 ≤ out_phis[i][j] ≤ 360.0
I_p    (output) array of arrays of arrays of arrays of upward radiances with dimensions (n_out_levels, n_out_thetas, n_out_phis, n_stokes)
I_m    (output) array of arrays of arrays of arrays of downward radiances with dimen-sions (n_out_levels, n_out_thetas, n_out_phis, n_stokes)
I_p_l    (output) array of arrays of arrays of arrays of arrays of upward radi-ance derivatives with dimensions (n_out_levels, n_derivs, n_out_thetas, n_out_phis, n_stokes)
I_m_l    (output) array of arrays of arrays of arrays of arrays of downward radiance derivatives with dimensions (n_out_levels, n_derivs, n_out_thetas, n_out_phis, n_stokes)
mean_p    (output) array of arrays of upward mean radiances with dimensions (n_out_levels, n_stokes)
mean_m    (output) array of arrays of downward mean radiances with dimensions (n_out_levels, n_stokes)
mean_p_l    (output) array of arrays of arrays of upward mean radiance derivatives with dimensions (n_out_levels, n_derivs, n_stokes)
mean_m_l    (output) array of arrays of arrays of downward mean radiance derivatives with dimensions (n_out_levels, n_derivs, n_stokes)
flux_p    (output) array of arrays of upward fluxes with dimensions (n_out_levels, n_stokes)
flux_m    (output) array of arrays of downward fluxes with dimensions (n_out_levels, n_stokes)
flux_p_l    (output) array of arrays of arrays of upward flux derivatives with dimen- sions (n_out_levels, n_derivs, n_stokes)
flux_m_l    (output) array of arrays of arrays of downward flux derivatives with dimensions (n_out_levels, n_derivs, n_stokes)
flux_div    (output) array of arrays of flux divergence with dimensions (n_out_levels, n_stokes)
flux_div_l    (output) array of arrays of arrays of flux divergence derivatives with dimensions (n_out_levels, n_derivs, n_stokes)
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_radiance(xrtm_data *d, enum xrtm_solver_mask solver, int n_out_phis, double **out_phis, double ****I_p, double ****I_m, double *****I_p_l, double *****I_m_l)

Description:
Run XRTM using a specified solver and return radiance results.
Arguments:
d    the xrtm_data structure which represents the instance created
solver    bit mask for the solver to be used
n_out_phis    number azimuth angles ϕ to return results for
out_phis    array of arrays of output azimuth angles with dimensions (n_out_thetas, n_out_phis), where 0.0 ≤ out_phis[i][j] ≤ 360.0
I_p    (output) array of arrays of arrays of arrays of upward radiances with dimensions (n_out_levels, n_out_thetas, n_out_phis, n_stokes)
I_m    (output) array of arrays of arrays of arrays of downward radiances with dimen-sions (n_out_levels, n_out_thetas, n_out_phis, n_stokes)
I_p_l    (output) array of arrays of arrays of arrays of arrays of upward radi-ance derivatives with dimensions (n_out_levels, n_derivs, n_out_thetas, n_out_phis, n_stokes)
I_m_l    (output) array of arrays of arrays of arrays of arrays of downward radiance derivatives with dimensions (n_out_levels, n_derivs, n_out_thetas, n_out_phis, n_stokes)
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_mean_radiance(xrtm_data *d, enum xrtm_solver_mask solver, double *mean_p, double *mean_m, double **mean_p_l, double **mean_m_l)

Description:
Run XRTM using a specified solver and return mean radiance results.
Arguments:
d    the xrtm_data structure which represents the instance created
solver    bit mask for the solver to be used
mean_p    (output) array of arrays of upward mean radiances with dimensions (n_out_levels, n_stokes)
mean_m    (output) array of arrays of downward mean radiances with dimensions (n_out_levels, n_stokes)
mean_p_l    (output) array of arrays of arrays of upward mean radiance derivatives with dimensions (n_out_levels, n_derivs, n_stokes)
mean_m_l    (output) array of arrays of arrays of downward mean radiance derivatives with dimensions (n_out_levels, n_derivs, n_stokes)
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_flux(xrtm_data *d, enum xrtm_solver_mask solver, double *flux_p, double *flux_m, double **flux_p_l, double **flux_m_l)

Description:
Run XRTM using a specified solver and return flux results.
Arguments:
d    the xrtm_data structure which represents the instance created
solver    bit mask for the solver to be used
flux_p    (output) array of arrays of upward fluxes with dimensions (n_out_levels, n_stokes)
flux_m    (output) array of arrays of downward fluxes with dimensions (n_out_levels, n_stokes)
flux_p_l    (output) array of arrays of arrays of upward flux derivatives with dimen- sions (n_out_levels, n_derivs, n_stokes)
flux_m_l    (output) array of arrays of arrays of downward flux derivatives with dimensions (n_out_levels, n_derivs, n_stokes)
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

int xrtm_flux_divergence(xrtm_data *d, enum xrtm_solver_mask solver, double *flux_div, double **flux_div_l)

Description:
Run XRTM using a specified solver and return flux divergence results.
Arguments:
d    the xrtm_data structure which represents the instance created
solver    bit mask for the solver to be used
flux_div    (output) array of arrays of flux divergence with dimensions (n_out_levels, n_stokes)
flux_div_l    (output) array of arrays of arrays of flux divergence derivatives with dimensions (n_out_levels, n_derivs, n_stokes)
Return value:
Zero with successful completion or XRTM_INT_ERROR on error.

3.6  Miscellaneous functions

const char *xrtm_get_version()

Description:
Get the XRTM version sstring.
Arguments:
None.
Return value:
Pointer to the XRTM version string.

3.7  Error Handling

All functions in the XRTM C interface return either an int or a double. If an error occurs when calling these functions then the return value will be set to one of the follwing error constants depending on the function's return type:
XRTM_INT_ERROR

Returned from functions returning an int when an error has occurred while calling the function.
XRTM_DBL_ERROR

Returned from functions returning a double when an error has occurred while calling the function.
When an error occurs XRTM will print an error message followed by a function call stack to the standard error (stderr) stream.

3.8  Example C program using XRTM

An example program using the C interface is at
     examples/example_c.c

and when the XRTM code is compiled properly the C example program will be compiled as
     examples/example_c


HEAD NEXT