![]() |
JURASSIC
|
JURASSIC library definitions. More...
#include "jurassic.h"Go to the source code of this file.
Functions | |
| void | analyze_avk (const ret_t *ret, const ctl_t *ctl, const atm_t *atm, const int *iqa, const int *ipa, const gsl_matrix *avk) |
| Analyze averaging kernel (AVK) matrix for retrieval diagnostics. More... | |
| void | analyze_avk_quantity (const gsl_matrix *avk, const int iq, const int *ipa, const size_t *n0, const size_t *n1, double *cont, double *res) |
| Analyze averaging kernel submatrix for a specific retrieved quantity. More... | |
| size_t | atm2x (const ctl_t *ctl, const atm_t *atm, gsl_vector *x, int *iqa, int *ipa) |
| Convert atmospheric data to state vector elements. More... | |
| void | atm2x_help (const double value, const int value_iqa, const int value_ip, gsl_vector *x, int *iqa, int *ipa, size_t *n) |
| Append a single atmospheric value to the state vector. More... | |
| void | cart2geo (const double *x, double *z, double *lon, double *lat) |
| Converts Cartesian coordinates to geographic coordinates. More... | |
| void | climatology (const ctl_t *ctl, atm_t *atm) |
| Initializes atmospheric climatology profiles. More... | |
| double | cos_sza (const double sec, const double lon, const double lat) |
| Calculates the cosine of the solar zenith angle. More... | |
| double | cost_function (const gsl_vector *dx, const gsl_vector *dy, const gsl_matrix *s_a_inv, const gsl_vector *sig_eps_inv) |
| Compute the normalized quadratic cost function for optimal estimation. More... | |
| double | ctmco2 (const double nu, const double p, const double t, const double u) |
| Compute carbon dioxide continuum (optical depth). More... | |
| double | ctmh2o (const double nu, const double p, const double t, const double q, const double u) |
| Compute water vapor continuum (optical depth). More... | |
| double | ctmn2 (const double nu, const double p, const double t) |
| Compute N₂ collision-induced absorption coefficient. More... | |
| double | ctmo2 (const double nu, const double p, const double t) |
| Compute O₂ collision-induced absorption coefficient. More... | |
| void | copy_atm (const ctl_t *ctl, atm_t *atm_dest, const atm_t *atm_src, const int init) |
| Copy or initialize atmospheric profile data. More... | |
| void | copy_obs (const ctl_t *ctl, obs_t *obs_dest, const obs_t *obs_src, const int init) |
| Copy or initialize observation geometry and radiance data. More... | |
| int | find_emitter (const ctl_t *ctl, const char *emitter) |
| Find gas species index by name. More... | |
| void | formod (const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs) |
| Execute the selected forward model. More... | |
| void | formod_continua (const ctl_t *ctl, const los_t *los, const int ip, double *beta) |
| Compute total extinction including gaseous continua. More... | |
| void | formod_fov (const ctl_t *ctl, obs_t *obs) |
| Apply field-of-view (FOV) convolution to modeled radiances. More... | |
| void | formod_pencil (const ctl_t *ctl, const tbl_t *tbl, const atm_t *atm, obs_t *obs, const int ir) |
| Compute line-of-sight radiances using the pencil-beam forward model. More... | |
| void | formod_rfm (const ctl_t *ctl, const atm_t *atm, obs_t *obs) |
| Interface routine for the Reference Forward Model (RFM). More... | |
| void | formod_srcfunc (const ctl_t *ctl, const tbl_t *tbl, const double t, double *src) |
| Interpolate the source function (Planck radiance) at a given temperature. More... | |
| void | geo2cart (const double z, const double lon, const double lat, double *x) |
| Converts geographic coordinates (longitude, latitude, altitude) to Cartesian coordinates. More... | |
| void | hydrostatic (const ctl_t *ctl, atm_t *atm) |
| Adjust pressure profile using the hydrostatic equation. More... | |
| void | idx2name (const ctl_t *ctl, const int idx, char *quantity) |
| Convert a quantity index to a descriptive name string. More... | |
| void | init_srcfunc (const ctl_t *ctl, tbl_t *tbl) |
| Initialize the source-function (Planck radiance) lookup table. More... | |
| void | intpol_atm (const ctl_t *ctl, const atm_t *atm, const double z, double *p, double *t, double *q, double *k) |
| Interpolate atmospheric state variables at a given altitude. More... | |
| void | intpol_tbl_cga (const ctl_t *ctl, const tbl_t *tbl, const los_t *los, const int ip, double tau_path[ND][NG], double tau_seg[ND]) |
| Interpolate emissivities and transmittances using the Curtis–Godson approximation (CGA). More... | |
| void | intpol_tbl_ega (const ctl_t *ctl, const tbl_t *tbl, const los_t *los, const int ip, double tau_path[ND][NG], double tau_seg[ND]) |
| Interpolate emissivities and transmittances using the Emissivity Growth Approximation (EGA). More... | |
| double | intpol_tbl_eps (const tbl_t *tbl, const int ig, const int id, const int ip, const int it, const double u) |
| Interpolate emissivity from lookup tables as a function of column density. More... | |
| double | intpol_tbl_u (const tbl_t *tbl, const int ig, const int id, const int ip, const int it, const double eps) |
| Interpolate column density from lookup tables as a function of emissivity. More... | |
| void | jsec2time (const double jsec, int *year, int *mon, int *day, int *hour, int *min, int *sec, double *remain) |
| Converts Julian seconds to calendar date and time components. More... | |
| void | kernel (const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs, gsl_matrix *k) |
| Compute the Jacobian (kernel) matrix by finite differences. More... | |
| int | locate_irr (const double *xx, const int n, const double x) |
| Locate index for interpolation on an irregular grid. More... | |
| int | locate_reg (const double *xx, const int n, const double x) |
| Locate index for interpolation on a regular (uniform) grid. More... | |
| int | locate_tbl (const float *xx, const int n, const double x) |
| Locate index for interpolation within emissivity table grids. More... | |
| void | matrix_invert (gsl_matrix *a) |
| Invert a square matrix, optimized for diagonal or symmetric positive-definite matrices. More... | |
| void | matrix_product (const gsl_matrix *a, const gsl_vector *b, const int transpose, gsl_matrix *c) |
| Compute structured matrix products of the form \(A^T B A\) or \(A B A^T\). More... | |
| size_t | obs2y (const ctl_t *ctl, const obs_t *obs, gsl_vector *y, int *ida, int *ira) |
| Convert observation radiances into a measurement vector. More... | |
| void | optimal_estimation (ret_t *ret, ctl_t *ctl, tbl_t *tbl, obs_t *obs_meas, obs_t *obs_i, atm_t *atm_apr, atm_t *atm_i, double *chisq) |
| Perform optimal estimation retrieval using Levenberg–Marquardt minimization. More... | |
| void | raytrace (const ctl_t *ctl, const atm_t *atm, obs_t *obs, los_t *los, const int ir) |
| Perform line-of-sight (LOS) ray tracing through the atmosphere. More... | |
| void | read_atm (const char *dirname, const char *filename, const ctl_t *ctl, atm_t *atm) |
| Read atmospheric input data from a file. More... | |
| void | read_atm_asc (FILE *in, const ctl_t *ctl, atm_t *atm) |
| Read atmospheric data in ASCII format. More... | |
| void | read_atm_bin (FILE *in, const ctl_t *ctl, atm_t *atm) |
| Read atmospheric data in binary format. More... | |
| void | read_ctl (int argc, char *argv[], ctl_t *ctl) |
| Read model control parameters from command-line and configuration input. More... | |
| void | read_matrix (const char *dirname, const char *filename, gsl_matrix *matrix) |
| Read a numerical matrix from an ASCII file. More... | |
| void | read_obs (const char *dirname, const char *filename, const ctl_t *ctl, obs_t *obs) |
| Read observation data from an input file. More... | |
| void | read_obs_asc (FILE *in, const ctl_t *ctl, obs_t *obs) |
| Read ASCII-formatted observation data from an open file stream. More... | |
| void | read_obs_bin (FILE *in, const ctl_t *ctl, obs_t *obs) |
| Read binary-formatted observation data from an open file stream. More... | |
| double | read_obs_rfm (const char *basename, const double z, const double *nu, const double *f, const int n) |
| Read and spectrally convolve an RFM output spectrum. More... | |
| void | read_ret (int argc, char *argv[], const ctl_t *ctl, ret_t *ret) |
| Read retrieval configuration and error parameters. More... | |
| void | read_rfm_spec (const char *filename, double *nu, double *rad, int *npts) |
| Read a Reference Forward Model (RFM) ASCII spectrum. More... | |
| void | read_shape (const char *filename, double *x, double *y, int *n) |
| Read a two-column shape function from an ASCII file. More... | |
| tbl_t * | read_tbl (const ctl_t *ctl) |
| Read all emissivity lookup tables for all gases and frequencies. More... | |
| void | read_tbl_asc (const ctl_t *ctl, tbl_t *tbl, const int id, const int ig) |
| Read a single ASCII emissivity lookup table. More... | |
| void | read_tbl_bin (const ctl_t *ctl, tbl_t *tbl, const int id, const int ig) |
| Read a single binary emissivity lookup table. More... | |
| void | read_tbl_gas (const ctl_t *ctl, tbl_t *tbl, const int id, const int ig) |
| Read one frequency block from a per-gas binary table file. More... | |
| int | read_tbl_gas_close (tbl_gas_t *g) |
| Close a per-gas binary table file and optionally rewrite metadata. More... | |
| int | read_tbl_gas_open (const char *path, tbl_gas_t *g) |
| Open a per-gas binary table file for reading and writing. More... | |
| int | read_tbl_gas_single (const tbl_gas_t *g, const double freq, tbl_t *tbl, const int id, const int ig) |
| Read one emissivity table block from a per-gas table file. More... | |
| double | scan_ctl (int argc, char *argv[], const char *varname, const int arridx, const char *defvalue, char *value) |
| Scan control file or command-line arguments for a configuration variable. More... | |
| void | set_cov_apr (const ret_t *ret, const ctl_t *ctl, const atm_t *atm, const int *iqa, const int *ipa, gsl_matrix *s_a) |
| Construct the a priori covariance matrix \(\mathbf{S_a}\) for retrieval parameters. More... | |
| void | set_cov_meas (const ret_t *ret, const ctl_t *ctl, const obs_t *obs, gsl_vector *sig_noise, gsl_vector *sig_formod, gsl_vector *sig_eps_inv) |
| Construct measurement error standard deviations and their inverse. More... | |
| void | tangent_point (const los_t *los, double *tpz, double *tplon, double *tplat) |
| Determine the tangent point along a line of sight (LOS). More... | |
| void | time2jsec (const int year, const int mon, const int day, const int hour, const int min, const int sec, const double remain, double *jsec) |
| Converts time components to seconds since January 1, 2000, 12:00:00 UTC. More... | |
| void | timer (const char *name, const char *file, const char *func, int line, int mode) |
| Simple wall-clock timer for runtime diagnostics. More... | |
| void | write_atm (const char *dirname, const char *filename, const ctl_t *ctl, const atm_t *atm) |
| Write atmospheric data to a file. More... | |
| void | write_atm_asc (FILE *out, const ctl_t *ctl, const atm_t *atm) |
| Write atmospheric data to an ASCII file. More... | |
| void | write_atm_bin (FILE *out, const ctl_t *ctl, const atm_t *atm) |
| Write atmospheric data to a binary file. More... | |
| void | write_atm_rfm (const char *filename, const ctl_t *ctl, const atm_t *atm) |
| Write atmospheric profile in RFM-compatible format. More... | |
| void | write_matrix (const char *dirname, const char *filename, const ctl_t *ctl, const gsl_matrix *matrix, const atm_t *atm, const obs_t *obs, const char *rowspace, const char *colspace, const char *sort) |
| Write a fully annotated matrix (e.g., Jacobian or gain matrix) to file. More... | |
| void | write_obs (const char *dirname, const char *filename, const ctl_t *ctl, const obs_t *obs) |
| Write observation data to an output file in ASCII or binary format. More... | |
| void | write_obs_asc (FILE *out, const ctl_t *ctl, const obs_t *obs) |
| Write observation data to an ASCII text file. More... | |
| void | write_obs_bin (FILE *out, const ctl_t *ctl, const obs_t *obs) |
| Write observation data in binary format to an output file stream. More... | |
| void | write_shape (const char *filename, const double *x, const double *y, const int n) |
| Write tabulated shape function data to a text file. More... | |
| void | write_stddev (const char *quantity, const ret_t *ret, const ctl_t *ctl, const atm_t *atm, const gsl_matrix *s) |
| Write retrieval standard deviation profiles to disk. More... | |
| void | write_tbl (const ctl_t *ctl, const tbl_t *tbl) |
| Write all emissivity lookup tables in the format specified by the control structure. More... | |
| void | write_tbl_asc (const ctl_t *ctl, const tbl_t *tbl) |
| Write all lookup tables in human-readable ASCII format. More... | |
| void | write_tbl_bin (const ctl_t *ctl, const tbl_t *tbl) |
| Write all lookup tables in compact binary format. More... | |
| void | write_tbl_gas (const ctl_t *ctl, const tbl_t *tbl) |
| Write lookup tables into per-gas binary table files with indexed blocks. More... | |
| int | write_tbl_gas_create (const char *path) |
| Create a new per-gas table file with an empty index. More... | |
| int | write_tbl_gas_single (tbl_gas_t *g, const double freq, const tbl_t *tbl, const int id, const int ig) |
| Append or overwrite a single frequency-table block in a per-gas file. More... | |
| void | x2atm (const ctl_t *ctl, const gsl_vector *x, atm_t *atm) |
| Map retrieval state vector back to atmospheric structure. More... | |
| void | x2atm_help (double *value, const gsl_vector *x, size_t *n) |
| Helper function to extract a single value from the retrieval state vector. More... | |
| void | y2obs (const ctl_t *ctl, const gsl_vector *y, obs_t *obs) |
Copy elements from the measurement vector y into the observation structure. More... | |
JURASSIC library definitions.
Definition in file jurassic.c.
| void analyze_avk | ( | const ret_t * | ret, |
| const ctl_t * | ctl, | ||
| const atm_t * | atm, | ||
| const int * | iqa, | ||
| const int * | ipa, | ||
| const gsl_matrix * | avk | ||
| ) |
Analyze averaging kernel (AVK) matrix for retrieval diagnostics.
Decomposes and evaluates the averaging kernel matrix \(\mathbf{A}\) to quantify the retrieval sensitivity and vertical resolution for each retrieved quantity (pressure, temperature, trace gases, extinction, etc.). The results are written to diagnostic atmospheric files for visualization.
| [in] | ret | Retrieval control and configuration structure (ret_t). |
| [in] | ctl | Global control structure (ctl_t) defining retrieval setup and quantities. |
| [in] | atm | Retrieved atmospheric state structure (atm_t). |
| [in] | iqa | Array mapping state vector indices to physical quantities (e.g., IDXP, IDXT, IDXQ(...)). |
| [in] | ipa | Array mapping state vector indices to atmospheric grid points. |
| [in] | avk | Averaging kernel matrix \(\mathbf{A}\) (gsl_matrix, n×n). |
The averaging kernel matrix \(\mathbf{A}\) describes the sensitivity of the retrieved state \(\hat{x}\) to the true state \(x\):
\[ \hat{x} - x_a = \mathbf{A} (x - x_a) \]
where \(x_a\) is the a priori state.
This function separates and evaluates:
Internally, this function:
analyze_avk_quantity() for each retrieved parameter type to compute quantitative measures of contribution and resolution.atm_cont.tab — contribution functions (sensitivity).atm_res.tab — resolution functions (vertical response width).iqa[i] to the quantity constants (e.g., IDXT, IDXQ(ig), etc.).ret->dir).ret->dir.Definition at line 29 of file jurassic.c.

| void analyze_avk_quantity | ( | const gsl_matrix * | avk, |
| const int | iq, | ||
| const int * | ipa, | ||
| const size_t * | n0, | ||
| const size_t * | n1, | ||
| double * | cont, | ||
| double * | res | ||
| ) |
Analyze averaging kernel submatrix for a specific retrieved quantity.
Computes the contribution (sensitivity) and resolution (information density) profiles for a given physical quantity from its corresponding submatrix of the full averaging kernel matrix \(\mathbf{A}\).
| [in] | avk | Averaging kernel matrix \(\mathbf{A}\) (gsl_matrix, n×n), describing the sensitivity of the retrieved state to the true state. |
| [in] | iq | Quantity index identifier (e.g., IDXP, IDXT, IDXQ(ig), IDXK(iw), IDXCLZ, etc.). |
| [in] | ipa | Array mapping state vector indices to atmospheric grid indices. |
| [in] | n0 | Array of starting indices for each quantity sub-block in the state vector. |
| [in] | n1 | Array of lengths (number of elements) for each quantity sub-block. |
| [out] | cont | Array of contribution values, representing the total sensitivity or response of the retrieval at each grid point. |
| [out] | res | Array of resolution measures, representing the vertical information density or averaging kernel width at each grid point. |
For a given quantity \(q\), the function extracts its corresponding block \(\mathbf{A}_q\) from the full averaging kernel matrix and computes:
\[ \text{cont}_i = \sum_j A_{ij}^{(q)} \]
giving the total response of retrieved element \(i\) to perturbations in all true elements \(j\) of the same quantity.\[ \text{res}_i = \frac{1}{A_{ii}^{(q)}} \]
approximating the degree of vertical resolution or smoothing at level \(i\).The results are stored in the provided arrays cont[] and res[], indexed according to the atmospheric profile index via ipa[].
iq) are processed.atm->np).A_ii is zero or near-zero, the corresponding resolution value may be undefined or numerically unstable.n0[iq] and n1[iq] must have been computed by analyze_avk().Definition at line 86 of file jurassic.c.
Convert atmospheric data to state vector elements.
Extracts selected quantities from an atmospheric profile (atm_t) according to retrieval settings in ctl_t, and appends them to the state vector x. For each included quantity, the function also stores its quantity index (iqa) and profile index (ipa).
The function respects retrieval altitude limits defined in ctl (e.g., retp_zmin/zmax, rett_zmin/zmax, etc.) and includes only variables flagged for retrieval (e.g., ret_clz, ret_sft, etc.).
| [in] | ctl | Control settings defining retrieval configuration and limits. |
| [in] | atm | Atmospheric profile data to extract from. |
| [out] | x | GSL vector to store state-vector elements. |
| [out] | iqa | Quantity index array corresponding to elements in x. |
| [out] | ipa | Profile index array corresponding to elements in x. |
Definition at line 110 of file jurassic.c.

| void atm2x_help | ( | const double | value, |
| const int | value_iqa, | ||
| const int | value_ip, | ||
| gsl_vector * | x, | ||
| int * | iqa, | ||
| int * | ipa, | ||
| size_t * | n | ||
| ) |
Append a single atmospheric value to the state vector.
Helper routine for atm2x(). Inserts one scalar value and its corresponding quantity and profile indices into the state vector and tracking arrays, then increments the element counter.
| [in] | value | Value to add to the state vector. |
| [in] | value_iqa | Quantity index (e.g., IDXP, IDXT, etc.). |
| [in] | value_ip | Profile index within the atmospheric profile. |
| [out] | x | GSL vector containing state-vector elements (may be NULL). |
| [out] | iqa | Quantity index array corresponding to x (may be NULL). |
| [out] | ipa | Profile index array corresponding to x (may be NULL). |
| [in,out] | n | Current number of elements in the state vector; incremented on return. |
Definition at line 164 of file jurassic.c.
| void cart2geo | ( | const double * | x, |
| double * | z, | ||
| double * | lon, | ||
| double * | lat | ||
| ) |
Converts Cartesian coordinates to geographic coordinates.
This function converts a point from Cartesian coordinates (x, y, z) to geographic coordinates (longitude, latitude, and altitude). It uses the spherical Earth approximation for the conversion.
| x | Pointer to an array containing the Cartesian coordinates (x, y, z) in kilometers. |
| z | Pointer to a double where the computed altitude (above the reference ellipsoid) will be stored, in kilometers. |
| lon | Pointer to a double where the computed longitude (in degrees) will be stored. |
| lat | Pointer to a double where the computed latitude (in degrees) will be stored. |
Definition at line 185 of file jurassic.c.
Initializes atmospheric climatology profiles.
This function populates the atmospheric state (atm) with standard climatological profiles of pressure, temperature, and trace gas concentrations (e.g., H2O, CH4, CO, O3, etc.) as a function of altitude. The profiles are based on reference climatological datasets and are used for atmospheric modeling and radiative transfer calculations.
| [in] | ctl | Control parameters structure. |
| [out] | atm | Atmospheric state structure to be populated with climatological data. |
Definition at line 200 of file jurassic.c.

| double cos_sza | ( | const double | sec, |
| const double | lon, | ||
| const double | lat | ||
| ) |
Calculates the cosine of the solar zenith angle.
This function computes the cosine of the solar zenith angle (SZA), which describes the angle between the local zenith (straight up) and the line connecting the observer to the center of the Sun. The cosine of the SZA is often used directly in radiative transfer and photochemical calculations to avoid unnecessary use of trigonometric inverse functions.
| sec | Seconds elapsed since 2000-01-01T12:00Z. |
| lon | Observer's longitude in degrees. |
| lat | Observer's latitude in degrees. |
The cosine of the solar zenith angle is computed based on the observer's position (longitude and latitude) and the specified time in seconds elapsed since 2000-01-01T12:00Z.
Definition at line 1138 of file jurassic.c.
| double cost_function | ( | const gsl_vector * | dx, |
| const gsl_vector * | dy, | ||
| const gsl_matrix * | s_a_inv, | ||
| const gsl_vector * | sig_eps_inv | ||
| ) |
Compute the normalized quadratic cost function for optimal estimation.
Evaluates the cost function
\[ J = \frac{1}{m} \left[ (\mathbf{y} - \mathbf{y_a})^T \mathbf{S_\epsilon^{-1}} (\mathbf{y} - \mathbf{y_a}) + (\mathbf{x} - \mathbf{x_a})^T \mathbf{S_a^{-1}} (\mathbf{x} - \mathbf{x_a}) \right] \]
where \(\mathbf{dx} = \mathbf{x} - \mathbf{x_a}\) and \(\mathbf{dy} = \mathbf{y} - \mathbf{y_a}\) represent the deviations from a priori state and measurement vectors, respectively.
| [in] | dx | State deviation vector (x - x_a). |
| [in] | dy | Measurement deviation vector (y - y_a). |
| [in] | s_a_inv | Inverse of the a priori covariance matrix ( \(\mathbf{S_a^{-1}}\)). |
| [in] | sig_eps_inv | Vector of inverse measurement uncertainties ( \(\mathbf{S_\epsilon^{-1/2}}\) diagonal elements). |
Definition at line 1179 of file jurassic.c.
| double ctmco2 | ( | const double | nu, |
| const double | p, | ||
| const double | t, | ||
| const double | u | ||
| ) |
Compute carbon dioxide continuum (optical depth).
Definition at line 1212 of file jurassic.c.
| double ctmh2o | ( | const double | nu, |
| const double | p, | ||
| const double | t, | ||
| const double | q, | ||
| const double | u | ||
| ) |
Compute water vapor continuum (optical depth).
Definition at line 2075 of file jurassic.c.
| double ctmn2 | ( | const double | nu, |
| const double | p, | ||
| const double | t | ||
| ) |
Compute N₂ collision-induced absorption coefficient.
Calculates the nitrogen (N₂) absorption coefficient due to collision-induced absorption (CIA) near the 4.3 µm CO₂ band using tabulated laboratory data for the absorption strength (B) and temperature exponent (β).
The function linearly interpolates B and β as a function of wavenumber, then applies a temperature- and pressure-dependent scaling relation to compute the absorption coefficient.
| [in] | nu | Wavenumber [cm⁻¹]. |
| [in] | p | Pressure [hPa]. |
| [in] | t | Temperature [K]. |
Definition at line 3127 of file jurassic.c.

| double ctmo2 | ( | const double | nu, |
| const double | p, | ||
| const double | t | ||
| ) |
Compute O₂ collision-induced absorption coefficient.
Calculates the molecular oxygen (O₂) absorption coefficient due to collision-induced absorption (CIA) using tabulated laboratory data for the absorption strength (B) and temperature exponent (β).
The function linearly interpolates B and β as a function of wavenumber and applies a pressure- and temperature-dependent scaling relation to compute the absorption coefficient.
| [in] | nu | Wavenumber [cm⁻¹]. |
| [in] | p | Pressure [hPa]. |
| [in] | t | Temperature [K]. |
Definition at line 3196 of file jurassic.c.

Copy or initialize atmospheric profile data.
Copies all fields from one atmospheric structure (atm_src) to another (atm_dest), including geolocation, thermodynamic, gas, extinction, cloud, and surface parameters. If init is nonzero, the destination fields are instead initialized to default values (zeros for most fields, unity for surface emissivity).
| [in] | ctl | Control structure defining array dimensions. |
| [out] | atm_dest | Destination atmospheric structure. |
| [in] | atm_src | Source atmospheric structure. |
| [in] | init | Initialization flag:
|
atm_src->np) determines the amount of data copied. The function performs shallow copies using memcpy for efficiency.Definition at line 3258 of file jurassic.c.
Copy or initialize observation geometry and radiance data.
Copies all observation fields from a source structure (obs_src) to a destination structure (obs_dest), including observer, view, and tangent point geometry, as well as radiance and transmittance data. If init is nonzero, radiance and transmittance values are reset to zero wherever finite values are found.
| [in] | ctl | Control structure defining the number of channels (ctl_t::nd) and other dimensions. |
| [out] | obs_dest | Destination observation structure. |
| [in] | obs_src | Source observation structure. |
| [in] | init | Initialization flag:
|
obs_src->nr) defines the copied data size. Shallow copies are performed via memcpy for efficiency.Definition at line 3308 of file jurassic.c.
| int find_emitter | ( | const ctl_t * | ctl, |
| const char * | emitter | ||
| ) |
Find gas species index by name.
Searches the list of emitter (gas) names defined in the control structure for a case-insensitive match to the given string. Returns the corresponding gas index if found, or -1 otherwise.
| [in] | ctl | Control structure containing the list of gas emitters. |
| [in] | emitter | Name of the gas species to search for (e.g. "H2O", "CO2"). |
strcasecmp().Definition at line 3346 of file jurassic.c.
Execute the selected forward model.
Computes synthetic radiances or brightness temperatures for the given atmospheric state and observation geometry using the selected forward-model method (CGA, EGA, or RFM). The function also applies hydrostatic adjustment, field-of-view convolution, and optional brightness-temperature conversion.
| [in] | ctl | Control structure defining model settings and options. |
| [in] | tbl | Emissivity and source-function lookup tables. |
| [in,out] | atm | Atmospheric profile; may be adjusted for hydrostatic balance. |
| [in,out] | obs | Observation geometry and radiance data; populated with model output. |
obs->rad elements marked as invalid (NaN) by applying an internal observation mask.Definition at line 3359 of file jurassic.c.

Compute total extinction including gaseous continua.
Calculates the total extinction coefficient for a given line-of-sight point by summing spectrally dependent extinction and optional gaseous continuum contributions (CO₂, H₂O, N₂, O₂).
The function updates the extinction array beta for each spectral channel based on line-by-line extinction and enabled continua flags specified in ctl.
| [in] | ctl | Control structure defining model setup and continuum options. |
| [in] | los | Line-of-sight data containing pressure, temperature, gas concentrations, and extinction coefficients. |
| [in] | ip | Index of the line-of-sight point to process. |
| [out] | beta | Array of total extinction coefficients [km⁻¹] per channel. |
Definition at line 3409 of file jurassic.c.

Apply field-of-view (FOV) convolution to modeled radiances.
Convolves pencil-beam radiances and transmittances along each ray path with the instrument field-of-view weighting function defined in the control structure. This simulates finite FOV effects on measured radiances and transmittances.
| [in] | ctl | Control structure containing FOV parameters (offsets ctl_t::fov_dz and weights ctl_t::fov_w). |
| [in,out] | obs | Observation structure; input pencil-beam data are replaced with FOV-convolved radiances and transmittances. |
| ERRMSG | if insufficient data are available for convolution. |
Definition at line 3445 of file jurassic.c.

| void formod_pencil | ( | const ctl_t * | ctl, |
| const tbl_t * | tbl, | ||
| const atm_t * | atm, | ||
| obs_t * | obs, | ||
| const int | ir | ||
| ) |
Compute line-of-sight radiances using the pencil-beam forward model.
Simulates monochromatic radiances and transmittances along a single line of sight using a layer-by-layer pencil-beam approximation. The model includes gaseous absorption, continuum extinction, surface emission, reflection, and optional solar illumination.
| [in] | ctl | Control structure defining model configuration, gas setup, surface type, and spectral parameters. |
| [in] | tbl | Emissivity and source-function lookup tables. |
| [in] | atm | Atmospheric state containing pressure, temperature, and gas profiles. |
| [in,out] | obs | Observation data; updated with modeled radiances and transmittances for the specified ray path. |
| [in] | ir | Index of the current ray path in obs. |
Definition at line 3510 of file jurassic.c.

Interface routine for the Reference Forward Model (RFM).
Prepares input data, executes the RFM executable, and imports simulated radiances and transmittances into the observation structure. The routine converts the atmospheric and geometric configuration from internal JURASSIC data structures into RFM-compatible driver and atmosphere files, then reads the RFM output spectra.
| [in] | ctl | Control structure defining model configuration, RFM executable path, HITRAN database, and spectral setup. |
| [in] | atm | Atmospheric state structure containing pressure, temperature, and gas profiles. |
| [in,out] | obs | Observation geometry and radiance data to be populated with RFM-computed radiances and transmittances. |
atm. It automatically determines whether the geometry is limb, nadir, or observer-based.RAD, TRA, MIX, CTM) based on the control settings. Temporary files such as rfm.drv, rfm.atm, and rad_*.asc are created and removed automatically.| ERRMSG | on inconsistent geometry, failed I/O, or system call errors. |
Definition at line 3645 of file jurassic.c.

Interpolate the source function (Planck radiance) at a given temperature.
Computes source function values by linearly interpolating between precomputed Planck radiances in the lookup table. The resulting radiance spectrum corresponds to the input temperature and is stored in src for all spectral channels.
| [in] | ctl | Control structure defining the number of spectral channels. |
| [in] | tbl | Emissivity and source-function lookup table containing Planck radiances (tbl_t::sr) and corresponding temperatures (tbl_t::st). |
| [in] | t | Temperature [K] for which the source function is evaluated. |
| [out] | src | Output array of interpolated source-function values [W·m⁻²·sr⁻¹·cm⁻¹] per spectral channel. |
Definition at line 3786 of file jurassic.c.

| void geo2cart | ( | const double | z, |
| const double | lon, | ||
| const double | lat, | ||
| double * | x | ||
| ) |
Converts geographic coordinates (longitude, latitude, altitude) to Cartesian coordinates.
This function converts geographic coordinates specified by longitude, latitude, and altitude into Cartesian coordinates. The Earth is approximated as a sphere with radius defined by the constant RE.
| z | The altitude above the Earth's surface in kilometers. |
| lon | The longitude in degrees. |
| lat | The latitude in degrees. |
| x | Pointer to an array of three doubles where the computed Cartesian coordinates (x, y, z) will be stored. |
The function computes the Cartesian coordinates using the given altitude, longitude, and latitude. It assumes the Earth is a perfect sphere and uses the following formulas:
RE is defined as the Earth's radius in kilometers. Definition at line 3803 of file jurassic.c.
Adjust pressure profile using the hydrostatic equation.
Recomputes the atmospheric pressure field to ensure hydrostatic equilibrium with respect to altitude and temperature. Starting from a reference altitude (ctl_t::hydz), the routine integrates the hydrostatic balance equation both upward and downward through the profile using small interpolation steps.
The air density is corrected for humidity using the mean molecular mass of dry air and water vapor.
| [in] | ctl | Control structure providing model constants and reference altitude (ctl_t::hydz) and H₂O index (ctl_t::ig_h2o). |
| [in,out] | atm | Atmospheric state; input temperatures, heights, and humidities are used to update pressure [hPa]. |
Definition at line 3823 of file jurassic.c.
| void idx2name | ( | const ctl_t * | ctl, |
| const int | idx, | ||
| char * | quantity | ||
| ) |
Convert a quantity index to a descriptive name string.
Translates a model quantity index (e.g., pressure, temperature, gas mixing ratio, extinction, or surface parameter) into a human-readable name. The function uses the index mapping defined in the control structure to assign appropriate labels.
| [in] | ctl | Control structure containing gas, window, cloud, and surface setup information. |
| [in] | idx | Quantity index (see IDXP, IDXT, IDXQ, IDXK, etc.). |
| [out] | quantity | Character buffer to receive the descriptive quantity name (e.g., "PRESSURE", "H2O", "CLOUD_HEIGHT", "SURFACE_EMISSIVITY_1000.0"). |
quantity using sprintf(). The caller must ensure sufficient buffer size (≥ LEN).Definition at line 3883 of file jurassic.c.
Initialize the source-function (Planck radiance) lookup table.
Computes channel-averaged Planck radiances for a range of temperatures and stores them in the source-function table. For each spectral channel, the Planck function is integrated over the instrument filter function defined in the corresponding filter file (*.filt).
| [in] | ctl | Control structure defining spectral channels and table base name. |
| [out] | tbl | Emissivity and source-function lookup table to populate. |
<tblbase>_<wavenumber>.filt.Definition at line 3922 of file jurassic.c.

| void intpol_atm | ( | const ctl_t * | ctl, |
| const atm_t * | atm, | ||
| const double | z, | ||
| double * | p, | ||
| double * | t, | ||
| double * | q, | ||
| double * | k | ||
| ) |
Interpolate atmospheric state variables at a given altitude.
Computes pressure, temperature, volume mixing ratios, and extinction coefficients at the specified altitude by interpolating between adjacent model levels in the atmospheric profile.
| [in] | ctl | Control structure defining the number of gases (ctl_t::ng) and spectral windows (ctl_t::nw). |
| [in] | atm | Atmospheric profile providing altitude, pressure, temperature, gas mixing ratios, and extinction data. |
| [in] | z | Target altitude [km]. |
| [out] | p | Interpolated pressure [hPa]. |
| [out] | t | Interpolated temperature [K]. |
| [out] | q | Interpolated gas volume mixing ratios [ppv], length ctl_t::ng. |
| [out] | k | Interpolated extinction coefficients [km⁻¹], length ctl_t::nw. |
Definition at line 3976 of file jurassic.c.

| void intpol_tbl_cga | ( | const ctl_t * | ctl, |
| const tbl_t * | tbl, | ||
| const los_t * | los, | ||
| const int | ip, | ||
| double | tau_path[ND][NG], | ||
| double | tau_seg[ND] | ||
| ) |
Interpolate emissivities and transmittances using the Curtis–Godson approximation (CGA).
Computes gas transmittance along a line-of-sight segment by interpolating precomputed emissivity values from lookup tables. The interpolation is performed in pressure, temperature, and column density space using bilinear (in p, T) and logarithmic (in p) interpolation.
| [in] | ctl | Control structure defining number of gases (ctl_t::ng) and channels (ctl_t::nd). |
| [in] | tbl | Emissivity lookup tables (tbl_t) containing tabulated pressure, temperature, and column density grids. |
| [in] | los | Line-of-sight structure providing Curtis–Godson mean parameters and column densities. |
| [in] | ip | Index of the current LOS point. |
| [in,out] | tau_path | Path transmittance array [nd][ng]; updated cumulatively for each gas. |
| [out] | tau_seg | Total segment transmittance per channel [nd]. |
tbl->eps) for each gas and channel.Definition at line 4001 of file jurassic.c.

| void intpol_tbl_ega | ( | const ctl_t * | ctl, |
| const tbl_t * | tbl, | ||
| const los_t * | los, | ||
| const int | ip, | ||
| double | tau_path[ND][NG], | ||
| double | tau_seg[ND] | ||
| ) |
Interpolate emissivities and transmittances using the Emissivity Growth Approximation (EGA).
Computes gas transmittance along a line-of-sight segment by interpolating emissivity values from lookup tables based on the Emissivity Growth Approximation (EGA). The interpolation is performed in pressure, temperature, and effective column density space derived from the local LOS properties.
| [in] | ctl | Control structure defining number of gases (ctl_t::ng) and channels (ctl_t::nd). |
| [in] | tbl | Emissivity lookup tables (tbl_t) containing tabulated pressure, temperature, and column density grids. |
| [in] | los | Line-of-sight structure providing local pressure, temperature, and column density data. |
| [in] | ip | Index of the current LOS point. |
| [in,out] | tau_path | Path transmittance array [nd][ng]; updated cumulatively for each gas. |
| [out] | tau_seg | Total segment transmittance per channel [nd]. |
tbl->eps) and performs bilinear interpolation in pressure and temperature.Definition at line 4090 of file jurassic.c.

|
inline |
Interpolate emissivity from lookup tables as a function of column density.
Retrieves the emissivity corresponding to a given column density u for a specific gas, channel, pressure level, and temperature index from the precomputed emissivity tables.
| [in] | tbl | Emissivity lookup tables (tbl_t). |
| [in] | ig | Gas index. |
| [in] | id | Channel index. |
| [in] | ip | Pressure level index. |
| [in] | it | Temperature level index. |
| [in] | u | Column density [molecules/cm²]. |
u for u < u_min.eps → 1 as u → ∞).tbl->u and tbl->eps.Definition at line 4186 of file jurassic.c.

|
inline |
Interpolate column density from lookup tables as a function of emissivity.
Returns the column density corresponding to a given emissivity eps for a specific gas, channel, pressure level, and temperature index from the precomputed emissivity tables.
| [in] | tbl | Emissivity lookup tables (tbl_t). |
| [in] | ig | Gas index. |
| [in] | id | Channel index. |
| [in] | ip | Pressure level index. |
| [in] | it | Temperature level index. |
| [in] | eps | Emissivity value (0–1). |
eps < eps_min, applies linear extrapolation proportional to emissivity.eps > eps_max, applies exponential extrapolation following the emissivity growth law.tbl->eps and tbl->u.Definition at line 4218 of file jurassic.c.

| void jsec2time | ( | const double | jsec, |
| int * | year, | ||
| int * | mon, | ||
| int * | day, | ||
| int * | hour, | ||
| int * | min, | ||
| int * | sec, | ||
| double * | remain | ||
| ) |
Converts Julian seconds to calendar date and time components.
This function converts Julian seconds to calendar date and time components, including year, month, day, hour, minute, and second. It also calculates the fractional part of the seconds.
| jsec | Julian seconds to convert. |
| year | Pointer to store the year. |
| mon | Pointer to store the month. |
| day | Pointer to store the day. |
| hour | Pointer to store the hour. |
| min | Pointer to store the minute. |
| sec | Pointer to store the second. |
| remain | Pointer to store the fractional part of seconds. |
The function initializes a time structure t0 with a fixed starting date and time. It then converts the Julian seconds to a time_t type by adding the seconds to the epoch time. Next, it converts the time_t value to a UTC time structure t1. Finally, it extracts the year, month, day, hour, minute, and second components from t1 and calculates the fractional part of seconds, which is stored in remain.
Definition at line 4250 of file jurassic.c.
Compute the Jacobian (kernel) matrix by finite differences.
Evaluates the sensitivity of the simulated radiances to each element of the atmospheric state vector by perturbing one parameter at a time and re-running the forward model. The result is the Jacobian matrix \( K = \partial y / \partial x \), where y is the measurement vector and x is the state vector.
| [in] | ctl | Control structure defining retrieval configuration and model setup. |
| [in] | tbl | Emissivity lookup tables used by the forward model. |
| [in] | atm | Atmospheric state vector and profile data. |
| [in] | obs | Observation geometry and radiance data. |
| [out] | k | Jacobian matrix [m×n], where m is the number of measurements and n the number of state variables. |
h depending on its physical type (pressure, temperature, VMR, etc.).Definition at line 4283 of file jurassic.c.

| int locate_irr | ( | const double * | xx, |
| const int | n, | ||
| const double | x | ||
| ) |
Locate index for interpolation on an irregular grid.
Finds the lower index ilo such that \( xx[ilo] \le x < xx[ilo+1] \) for monotonically increasing or decreasing grids.
Used in interpolation routines for altitude, pressure, temperature, or wavenumber profiles that are not evenly spaced.
| [in] | xx | Array of monotonic grid values (increasing or decreasing). |
| [in] | n | Number of grid points. |
| [in] | x | Target value to locate within the grid range. |
ilo of the lower grid point surrounding x.x lies within the range of xx; no bounds checking beyond the first and last grid points is performed.Definition at line 4378 of file jurassic.c.
| int locate_reg | ( | const double * | xx, |
| const int | n, | ||
| const double | x | ||
| ) |
Locate index for interpolation on a regular (uniform) grid.
Computes the lower index i such that \( xx[i] \le x < xx[i+1] \) for evenly spaced grid points. Used for quick index lookup when the grid spacing is constant.
| [in] | xx | Array of regularly spaced grid values. |
| [in] | n | Number of grid points. |
| [in] | x | Target value to locate within the grid range. |
i of the lower grid point surrounding x.[0, n - 2] to avoid out-of-bounds indices.Definition at line 4408 of file jurassic.c.
|
inline |
Locate index for interpolation within emissivity table grids.
Finds the lower index ilo such that \( xx[ilo] \le x < xx[ilo+1] \) in a monotonically increasing single-precision grid.
Used for emissivity and column density interpolation in table-based routines such as intpol_tbl_eps and intpol_tbl_u.
| [in] | xx | Monotonic (increasing) single-precision grid array. |
| [in] | n | Number of grid points. |
| [in] | x | Target value to locate within the grid range. |
ilo of the lower grid point surrounding x.float to minimize memory use.x lies within the table range.Definition at line 4427 of file jurassic.c.
| void matrix_invert | ( | gsl_matrix * | a | ) |
Invert a square matrix, optimized for diagonal or symmetric positive-definite matrices.
Performs in-place inversion of the matrix \(\mathbf{A}\) using either:
| [in,out] | a | Square matrix (gsl_matrix) to be inverted in place. |
The function first checks whether the input matrix is diagonal by testing all off-diagonal elements. If diagonal, each diagonal element \(a_{ii}\) is replaced by its reciprocal \(1/a_{ii}\).
For non-diagonal matrices, a Cholesky decomposition is performed:
\[ \mathbf{A} = \mathbf{L}\mathbf{L}^T \]
followed by inversion using the Cholesky factors, yielding \(\mathbf{A}^{-1}\).
This approach assumes \(\mathbf{A}\) is symmetric and positive-definite.
Definition at line 4449 of file jurassic.c.
| void matrix_product | ( | const gsl_matrix * | a, |
| const gsl_vector * | b, | ||
| const int | transpose, | ||
| gsl_matrix * | c | ||
| ) |
Compute structured matrix products of the form \(A^T B A\) or \(A B A^T\).
Evaluates matrix products commonly used in covariance propagation and optimal estimation, depending on the specified transpose mode:
The vector \(\mathbf{b}\) represents the diagonal elements of \(\mathbf{B}\), i.e. a diagonal weighting or covariance matrix.
| [in] | a | Input matrix \(\mathbf{A}\) (size m×n). |
| [in] | b | Vector representing the diagonal of \(\mathbf{B}\) (length m or n). |
| [in] | transpose | Operation selector:
|
| [out] | c | Output matrix to store the resulting product. |
dgemm routines for efficiency: \[ A^T B A = (B^{1/2}A)^T (B^{1/2}A), \quad A B A^T = (A B^{1/2}) (A B^{1/2})^T \]
transpose is not 1 or 2, the function performs no operation.Definition at line 4479 of file jurassic.c.
Convert observation radiances into a measurement vector.
Extracts all finite radiance values from the observation structure and stores them sequentially in a GSL vector.
Optionally records detector (id) and ray path (ir) indices for each measurement element.
| [in] | ctl | Control structure containing observation setup (e.g. number of detectors). |
| [in] | obs | Observation data structure containing radiances. |
| [out] | y | Measurement vector to store radiances (may be NULL). |
| [out] | ida | Optional array to store detector indices (may be NULL). |
| [out] | ira | Optional array to store ray path indices (may be NULL). |
nd) and ray paths (nr).NaN or Inf) radiance values.ida and ira must be preallocated with sufficient size to hold all finite radiances if provided.Definition at line 4524 of file jurassic.c.
| void optimal_estimation | ( | ret_t * | ret, |
| ctl_t * | ctl, | ||
| tbl_t * | tbl, | ||
| obs_t * | obs_meas, | ||
| obs_t * | obs_i, | ||
| atm_t * | atm_apr, | ||
| atm_t * | atm_i, | ||
| double * | chisq | ||
| ) |
Perform optimal estimation retrieval using Levenberg–Marquardt minimization.
This function performs an optimal estimation of atmospheric state variables based on measured observations, a priori information, and forward radiative transfer modeling. The estimation follows the Rodgers (2000) formalism and uses a Levenberg–Marquardt algorithm to iteratively minimize the cost function:
χ² = (x - x_a)^T S_a⁻¹ (x - x_a) + (y - F(x))^T S_ε⁻¹ (y - F(x)),
where x is the atmospheric state vector, x_a is its a priori estimate, S_a is the a priori covariance matrix, y is the measurement vector, F(x) is the forward model, and S_ε is the measurement error covariance.
The routine updates the atmospheric state until convergence criteria are met or the maximum number of iterations is reached. Optionally, the full retrieval error budget and averaging kernel analysis are computed.
| [out] | ret | Retrieval configuration and output container. Determines convergence, kernel recomputation frequency, and error analysis options. |
| [in] | ctl | Control parameters describing problem setup (grids, species, etc.). |
| [in] | tbl | Lookup tables required by the forward model. |
| [in] | obs_meas | Measured observations used as input. |
| [out] | obs_i | Intermediate and final modeled observations corresponding to the retrieved state. |
| [in] | atm_apr | A priori atmospheric state used as reference. |
| [out] | atm_i | Atmospheric state vector to be iteratively retrieved and updated. |
| [out] | chisq | Final value of the cost function (χ²) upon convergence. |
ret->err_ana), the function produces:Definition at line 4551 of file jurassic.c.

Perform line-of-sight (LOS) ray tracing through the atmosphere.
Computes the geometric path of a viewing ray from the observer to the atmosphere (and possibly the surface), accounting for spherical geometry, optional refraction, and cloud or surface interactions.
Fills the LOS structure with pressure, temperature, gas concentrations, extinction, and path length at each step.
| [in] | ctl | Control structure containing model and numerical settings. |
| [in] | atm | Atmospheric state structure (profiles of p, T, q, k, etc.). |
| [in,out] | obs | Observation geometry and radiance data; updated tangent point. |
| [out] | los | Line-of-sight structure to be populated with sampled quantities. |
| [in] | ir | Index of the current ray path in the observation set. |
ds determined by altitude and user-specified controls (rayds, raydz).n(p, T).NLOS.Definition at line 4823 of file jurassic.c.

Read atmospheric input data from a file.
This function reads atmospheric data from the file specified by filename, optionally prefixed by the directory dirname. The file format (ASCII or binary) is determined by the control structure ctl. The data are stored in the atmospheric structure atm.
Supported file formats are:
ctl->atmfmt == 1)ctl->atmfmt == 2)The function initializes the atmospheric data container, opens the specified file, reads its contents according to the selected format, and performs sanity checks on the number of data points. It also logs basic statistical information about the loaded atmospheric fields (time, altitude, longitude, latitude, pressure, temperature, and species/emitter mixing ratios).
| dirname | Optional directory path where the atmospheric file resides. If NULL, only filename is used. |
| filename | Name of the atmospheric data file to read. |
| ctl | Pointer to a control structure specifying input parameters, file format, number of emitters, spectral windows, and additional atmospheric configuration. |
| atm | Pointer to an atmospheric data structure that will be filled with the values read from the file. The structure must be allocated before calling this function. |
atm structure has been properly allocated, and that the control parameters in ctl are valid before calling this function.Definition at line 5059 of file jurassic.c.

Read atmospheric data in ASCII format.
This function parses atmospheric input data from an opened ASCII file stream (in) and stores the values in the atmospheric structure atm. The number and type of fields to read are determined by the control structure ctl. Each line of the ASCII file corresponds to a single atmospheric data point.
The expected order of fields in each line is:
ctl->ng values)ctl->nw values)Additionally, if cloud or surface layer parameters are enabled in ctl, they are read once from the first line only:
ctl->ncl values)ctl->nsf values)| in | Pointer to an already opened input file stream containing ASCII atmospheric data. |
| ctl | Pointer to a control structure defining the number of gases, spectral windows, and whether cloud or surface layer information should be read. |
| atm | Pointer to an initialized atmospheric structure where the parsed data points will be stored. The function updates atm->np to reflect the number of successfully read data records. |
ERRMSG if more than NP data points are encountered, or if the input format deviates from expectations.Definition at line 5142 of file jurassic.c.
Read atmospheric data in binary format.
This function reads atmospheric input data from a binary file stream (in) and stores the decoded values in the atmospheric structure atm. The expected binary format is predefined and must match the configuration provided in the control structure ctl. The function reads a header containing metadata describing the dataset, followed by the atmospheric fields and optional cloud/surface layer properties.
The binary file layout is expected to follow this structure:
ng)nw)ncl)nsf) These must match the corresponding values in ctl.np)np for time, altitude, longitude, latitude, pressure, temperaturenpnpctl->ncl values)ctl->nsf values)| in | Pointer to an opened binary input file stream. |
| ctl | Pointer to a control structure specifying expected dimensions of atmospheric fields and optional layers. Used for header validation. |
| atm | Pointer to an allocated atmospheric data structure that will be filled with the contents of the binary file. All arrays must be allocated prior to calling this function. |
ERRMSG if:Definition at line 5186 of file jurassic.c.
| void read_ctl | ( | int | argc, |
| char * | argv[], | ||
| ctl_t * | ctl | ||
| ) |
Read model control parameters from command-line and configuration input.
Parses all numerical and string parameters required to initialize a JURASSIC simulation, including atmospheric composition, radiative channels, cloud and surface options, continua, ray-tracing setup, retrieval parameters, and output settings.
Populates the ctl_t structure with all configuration values.
| [in] | argc | Argument count from the command line. |
| [in] | argv | Argument vector containing user-specified options. |
| [out] | ctl | Control structure to be filled with parsed settings. |
NG, EMITTER),ND, NU, WINDOW),NCL, CLNU),NSF, SFNU, SFTYPE, SFSZA),HYDZ),CTM_CO2, CTM_H2O, CTM_N2, CTM_O2),REFRAC, RAYDS, RAYDZ),FOV),RETP_ZMIN, RETQ_ZMAX, etc.),WRITE_BBT, WRITE_MATRIX),RFMBIN, RFMHIT, RFMXSC).NG > NG_MAX).NCL > 1, NSF > 1).ERRMSG() aborts.Definition at line 5272 of file jurassic.c.

| void read_matrix | ( | const char * | dirname, |
| const char * | filename, | ||
| gsl_matrix * | matrix | ||
| ) |
Read a numerical matrix from an ASCII file.
Loads values into a GSL matrix from a text file containing sparse or indexed entries in tabular format.
Each valid line is parsed for row and column indices and the corresponding value.
| [in] | dirname | Directory path containing the matrix file (may be NULL). |
| [in] | filename | Name of the matrix file to read. |
| [out] | matrix | Pointer to the GSL matrix to be filled with values. |
gsl_matrix_set().Definition at line 5390 of file jurassic.c.
Read observation data from an input file.
This function reads atmospheric observation data from the specified file and stores the results in the provided obs_t structure. The file may be in ASCII or binary format, depending on the control settings passed via ctl_t. After reading the data, the routine performs basic validation (e.g., verifies that at least one observation entry was loaded) and logs diagnostic statistics such as ranges of times, observer coordinates, view point coordinates, tangent point coordinates, radiance or brightness temperature values, and transmittances.
The input file path is constructed from the provided directory and filename. If a directory name is given, the file is assumed to reside within it; otherwise, the filename is used as-is. Depending on the value of ctl->obsfmt, the function dispatches either to read_obs_asc() for ASCII files or read_obs_bin() for binary files.
| [in] | dirname | Directory containing the input file, or NULL to use only filename. |
| [in] | filename | Name of the observation file to read. |
| [in] | ctl | Pointer to a control structure specifying file format and other options. |
| [out] | obs | Pointer to an observation structure where the data will be stored. |
obs structure must be properly allocated before calling this function. The function assumes its arrays are large enough to store all values contained in the input file.Definition at line 5430 of file jurassic.c.

Read ASCII-formatted observation data from an open file stream.
This function parses atmospheric observation data from an ASCII text file and stores it in the provided obs_t structure. Each line in the input file is expected to contain numerical values representing a single observation record, including time, observer coordinates, view point coordinates, tangent point coordinates, radiance or brightness temperature values, and transmittances. The number of radiance and transmittance values per record is determined by ctl->nd.
The function reads the file line by line, tokenizes the data fields, and fills the corresponding observation arrays. The number of successfully read observation entries is stored in obs->nr.
| [in] | in | Open file pointer from which the ASCII observation data will be read. The file must already be opened in read mode. |
| [in] | ctl | Control structure containing metadata such as the number of spectral channels (nd). |
| [out] | obs | Observation structure where parsed data will be stored. |
obs structure has been preallocated with sufficient space for all records and spectral channels. No memory allocation is performed inside this routine.NR.Definition at line 5516 of file jurassic.c.
Read binary-formatted observation data from an open file stream.
This C function reads observation data stored in a compact binary format and initializes the provided obs_t structure with the values retrieved from the input file. The binary format begins with a header that contains a magic identifier and the expected number of spectral channels. The number of channels in the file must match ctl->nd, otherwise the routine aborts with an error.
After verifying the header, the function reads the number of ray paths and then sequentially loads arrays corresponding to observation time, observer location, view point location, tangent point location, radiance (or brightness temperature), and transmittance data. The number of ray paths is assigned to obs->nr. All arrays must have been allocated prior to calling this function.
| [in] | in | Open file stream positioned at the beginning of the binary observation data. The file must be opened in binary mode. |
| [in] | ctl | Pointer to a control structure specifying the number of spectral channels (nd) and other configuration settings. |
| [out] | obs | Pointer to an observation structure where the decoded binary data will be stored. |
obs have sufficient capacity for the data being read.NR is encountered, or if any read operation fails.Definition at line 5553 of file jurassic.c.
| double read_obs_rfm | ( | const char * | basename, |
| const double | z, | ||
| const double * | nu, | ||
| const double * | f, | ||
| const int | n | ||
| ) |
Read and spectrally convolve an RFM output spectrum.
Opens the appropriate RFM ASCII spectrum file for a given tangent or observation altitude and convolves the high-resolution spectrum with a provided instrument filter function.
Returns the integrated (filtered) radiance value.
| [in] | basename | Base filename of the RFM output (e.g. "rad" or "tra"). |
| [in] | z | Altitude [km] used to select the corresponding RFM file. |
| [in] | nu | Wavenumber grid [cm⁻¹] of the filter function. |
| [in] | f | Filter transmission values corresponding to nu. |
| [in] | n | Number of points in nu and f arrays. |
basename_<altitude_in_meters>.asc (e.g. rad_04500.asc).nurfm) and radiances (rad).\[ R = \frac{\int f(\nu) \, I(\nu) \, d\nu}{\int f(\nu) \, d\nu} \]
nu and f arrays are monotonic and have at least two points.Definition at line 5623 of file jurassic.c.

Read retrieval configuration and error parameters.
Initializes the retrieval control structure (ret_t) by reading all iteration and uncertainty parameters from the command line or an input control file using the scan_ctl() interface.
| [in] | argc | Number of command-line arguments. |
| [in] | argv | Command-line argument vector. |
| [in] | ctl | Pointer to global control structure (ctl_t) defining the number of emitters, detectors, windows, clouds, etc. |
| [out] | ret | Pointer to retrieval configuration structure (ret_t) to be populated with iteration and error settings. |
The function performs the following initialization steps:
KERNEL_RECOMP — number of iterations between kernel recomputations. CONV_ITMAX — maximum number of retrieval iterations. CONV_DMIN — minimum normalized step size for convergence.ERR_ANA — enables or disables retrieval error analysis (0 = off, 1 = on).ERR_FORMOD[id] — relative (%) forward model uncertainty per detector channel. ERR_NOISE[id] — absolute instrument noise per detector channel write_bbt.ERR_PRESS, ERR_PRESS_CZ, ERR_PRESS_CH — pressure error [%] and correlation lengths [km]. ERR_TEMP, ERR_TEMP_CZ, ERR_TEMP_CH — temperature error [K] and correlation lengths [km].ERR_Q[ig], ERR_Q_CZ[ig], ERR_Q_CH[ig] — per gas [%] and correlation lengths [km].ERR_K[iw], ERR_K_CZ[iw], ERR_K_CH[iw] — per spectral window [km⁻¹] and correlation lengths [km].ERR_CLZ — cloud top height error [km]. ERR_CLDZ — cloud depth error [km]. ERR_CLK[icl] — cloud extinction error per frequency [km⁻¹].ERR_SFT — surface temperature error [K]. ERR_SFEPS[isf] — surface emissivity errors (dimensionless).-999 indicate uncorrelated (diagonal) treatment.ctl and ret dimensions.Definition at line 5681 of file jurassic.c.

| void read_rfm_spec | ( | const char * | filename, |
| double * | nu, | ||
| double * | rad, | ||
| int * | npts | ||
| ) |
Read a Reference Forward Model (RFM) ASCII spectrum.
Parses an RFM output file containing high-resolution spectral radiances and fills the provided arrays with wavenumber and radiance values.
| [in] | filename | Name of the RFM ASCII spectrum file (e.g. "rad_04500.asc"). |
| [out] | nu | Array to receive the spectral wavenumber grid [cm⁻¹]. |
| [out] | rad | Array to receive the corresponding radiances. |
| [out] | npts | Pointer to integer receiving the number of spectral points. |
npts),nu0 [cm⁻¹],dnu [cm⁻¹],nu1 [cm⁻¹].nu0 and nu1: \[ \nu_i = \nu_0 + i \, \frac{\nu_1 - \nu_0}{N - 1}, \quad i = 0,\dots,N-1 \]
RFMNPTS..asc output.Definition at line 5734 of file jurassic.c.
| void read_shape | ( | const char * | filename, |
| double * | x, | ||
| double * | y, | ||
| int * | n | ||
| ) |
Read a two-column shape function from an ASCII file.
Loads tabulated x–y data pairs (e.g., filter transmission or field-of-view weighting function) into the provided arrays.
| [in] | filename | Name of the ASCII file containing the shape function. |
| [out] | x | Array to receive the abscissa values. |
| [out] | y | Array to receive the ordinate values. |
| [out] | n | Pointer to integer receiving the number of data points read. |
x), typically wavenumber [cm⁻¹] or angular offset [deg].y), typically transmission or weighting value.NSHAPE.x values.Definition at line 5795 of file jurassic.c.
Read all emissivity lookup tables for all gases and frequencies.
This function allocates a new tbl_t structure and fills it by reading emissivity lookup tables for each trace gas (ig) and each frequency index (id) specified in the control structure. The lookup tables may be read from ASCII, binary, or per-gas table files depending on ctl->tblfmt.
After loading all tables, the source function is initialized with init_srcfunc().
| ctl | Pointer to control structure containing table metadata, number of gases, number of frequencies, filenames, etc. |
tbl_t structure containing all loaded lookup-table data. The caller owns the returned pointer and must free it when done.ERRMSG() if unexpected table formats or dimension overflows occur.Definition at line 5837 of file jurassic.c.

Read a single ASCII emissivity lookup table.
This reads one ASCII table corresponding to frequency index id and gas index ig. The table format is:
pressure temperature column_density emissivity
The function automatically determines the pressure, temperature, and column-density indices based on new values appearing in the file.
Out-of-range values for u or eps are skipped and counted.
| ctl | Pointer to control structure specifying filenames and grids. |
| tbl | Pointer to the table structure to be filled. |
| id | Frequency index. |
| ig | Gas index. |
ERRMSG() if table dimensions exceed TBLNP/TBLNT/TBLNU.Definition at line 5889 of file jurassic.c.
Read a single binary emissivity lookup table.
Reads the binary table stored as:
The function fills the corresponding entries of the tbl_t structure.
| ctl | Pointer to control structure specifying filenames and grids. |
| tbl | Pointer to the table structure to be filled. |
| id | Frequency index. |
| ig | Gas index. |
ERRMSG() if table dimensions exceed TBLNP/TBLNT/TBLNU.Definition at line 5987 of file jurassic.c.
Read one frequency block from a per-gas binary table file.
Opens the gas-specific table file (e.g., base_emitter.tbl) and reads the table block corresponding to frequency ctl->nu[id]. The block is appended to the in-memory tbl_t.
| ctl | Pointer to control structure containing table metadata. |
| tbl | Pointer to table structure to populate. |
| id | Frequency index. |
| ig | Gas index. |
Definition at line 6047 of file jurassic.c.

| int read_tbl_gas_close | ( | tbl_gas_t * | g | ) |
Close a per-gas binary table file and optionally rewrite metadata.
If the table was modified (g->dirty != 0), the header and index are rewritten before closing the file. After closing, memory associated with the table index is freed.
| g | Pointer to an open gas-table handle. |
Definition at line 6077 of file jurassic.c.
| int read_tbl_gas_open | ( | const char * | path, |
| tbl_gas_t * | g | ||
| ) |
Open a per-gas binary table file for reading and writing.
Reads and validates the file header, then loads the entire index of table blocks. The resulting tbl_gas_t structure tracks the file pointer, index, and table count.
| path | Path to the .tbl file. |
| g | Output parameter: populated table-file handle. |
ERRMSG() on invalid magic or format.Definition at line 6110 of file jurassic.c.
| int read_tbl_gas_single | ( | const tbl_gas_t * | g, |
| const double | freq, | ||
| tbl_t * | tbl, | ||
| const int | id, | ||
| const int | ig | ||
| ) |
Read one emissivity table block from a per-gas table file.
Locates the index entry corresponding to the requested frequency freq. If found, seeks to the stored offset and reads:
The data are stored into tbl[id][ig].
| g | Pointer to an open gas-table handle. |
| freq | Frequency to be read. |
| tbl | Pointer to output table structure. |
| id | Frequency index. |
| ig | Gas index. |
Definition at line 6139 of file jurassic.c.
| double scan_ctl | ( | int | argc, |
| char * | argv[], | ||
| const char * | varname, | ||
| const int | arridx, | ||
| const char * | defvalue, | ||
| char * | value | ||
| ) |
Scan control file or command-line arguments for a configuration variable.
Searches for a named variable in the JURASSIC control file or command-line arguments, returning its value as a double. Optionally stores the value as a string and supports array-style parameters (e.g., EMITTER[0], NU[5]).
| [in] | argc | Number of command-line arguments. |
| [in] | argv | Command-line argument vector. |
| [in] | varname | Name of the control variable to read. |
| [in] | arridx | Array index (use -1 for scalar variables). |
| [in] | defvalue | Default value if variable is not found (can be empty). |
| [out] | value | Optional pointer to a string buffer receiving the value (may be NULL if only numeric output is required). |
double.argv[1]), unless it starts with '-'.VAR (scalar)VAR[index] (explicit array index)VAR[*] (wildcard entry applying to all indices)defvalue is used (if non-empty).ctl_t checks to ensure consistency.Definition at line 6217 of file jurassic.c.
| void set_cov_apr | ( | const ret_t * | ret, |
| const ctl_t * | ctl, | ||
| const atm_t * | atm, | ||
| const int * | iqa, | ||
| const int * | ipa, | ||
| gsl_matrix * | s_a | ||
| ) |
Construct the a priori covariance matrix \(\mathbf{S_a}\) for retrieval parameters.
Builds the full a priori covariance matrix based on specified retrieval error assumptions and correlation lengths defined in the ret_t structure. Each diagonal element represents the variance of an individual state vector element, while off-diagonal terms encode spatial correlations between parameters of the same type.
| [in] | ret | Retrieval configuration and error parameters (ret_t). |
| [in] | ctl | Control structure defining retrieval setup and state vector mapping (ctl_t). |
| [in] | atm | Atmospheric profile structure containing geolocation and altitude information (atm_t). |
| [in] | iqa | Index array linking state vector elements to physical quantities (e.g. pressure, temperature, gas, extinction, etc.). |
| [in] | ipa | Index array linking state vector elements to atmospheric grid points. |
| [out] | s_a | Output a priori covariance matrix \(\mathbf{S_a}\) (gsl_matrix), dimension n×n. |
atm2x) and scales them according to their a priori uncertainties:\[ \rho_{ij} = \exp\left(-\frac{d_{ij}}{L_h} - \frac{|z_i - z_j|}{L_v}\right) \]
where:iqa[i] == iqa[j]) are assumed to share correlation properties.ret_t, differing by parameter type:err_press_cz, err_press_cherr_temp_cz, err_temp_cherr_q_cz[], err_q_ch[]err_k_cz[], err_k_ch[]Definition at line 6286 of file jurassic.c.

| void set_cov_meas | ( | const ret_t * | ret, |
| const ctl_t * | ctl, | ||
| const obs_t * | obs, | ||
| gsl_vector * | sig_noise, | ||
| gsl_vector * | sig_formod, | ||
| gsl_vector * | sig_eps_inv | ||
| ) |
Construct measurement error standard deviations and their inverse.
Builds the total measurement uncertainty vector used in the optimal estimation retrieval, accounting for both instrument noise and forward model (systematic) errors.
| [in] | ret | Retrieval configuration and error parameters (ret_t). |
| [in] | ctl | Control structure defining spectral channels and setup (ctl_t). |
| [in] | obs | Observation dataset (obs_t), containing measured radiances or brightness temperatures. |
| [out] | sig_noise | Vector of instrument noise standard deviations (gsl_vector), length m. |
| [out] | sig_formod | Vector of forward model error standard deviations (gsl_vector), length m. |
| [out] | sig_eps_inv | Vector of inverse total standard deviations, \(\sigma_\epsilon^{-1}\), used for normalization. |
\[ \sigma_{\epsilon,i}^2 = \sigma_{\text{noise},i}^2 + \sigma_{\text{formod},i}^2 \]
and stores its reciprocal square root:\[ (\sigma_{\epsilon,i}^{-1}) = \frac{1}{\sqrt{\sigma_{\epsilon,i}^2}} \]
sig_noise)** ret->err_noise[id] for each spectral channel. The noise term is always included in the fit.sig_formod)** ret->err_formod[id]) of the measured radiance (or brightness temperature) per channel. This represents uncertainty due to imperfect forward modeling.sig_eps_inv) is used to normalize the measurement residuals \((y - F(x))\) in the cost function.NAN.obs and ctl are consistent in dimension and indexing.Definition at line 6397 of file jurassic.c.

| void tangent_point | ( | const los_t * | los, |
| double * | tpz, | ||
| double * | tplon, | ||
| double * | tplat | ||
| ) |
Determine the tangent point along a line of sight (LOS).
Computes the location of the tangent point — the point of minimum altitude — along the current line of sight, based on the LOS geometry stored in los_t.
| [in] | los | Pointer to the line-of-sight (LOS) structure containing altitude, longitude, latitude, and segment length data. |
| [out] | tpz | Pointer to variable receiving tangent point altitude [km]. |
| [out] | tplon | Pointer to variable receiving tangent point longitude [deg]. |
| [out] | tplat | Pointer to variable receiving tangent point latitude [deg]. |
ds) must be consistent with the geometric spacing between altitude points for the interpolation to be accurate.Definition at line 6441 of file jurassic.c.

| void time2jsec | ( | const int | year, |
| const int | mon, | ||
| const int | day, | ||
| const int | hour, | ||
| const int | min, | ||
| const int | sec, | ||
| const double | remain, | ||
| double * | jsec | ||
| ) |
Converts time components to seconds since January 1, 2000, 12:00:00 UTC.
This function calculates the number of seconds elapsed since January 1, 2000, 12:00:00 UTC, based on the provided year, month, day, hour, minute, and second. It also includes a fractional part to represent the remaining seconds.
| year | The year. |
| mon | The month (1-12). |
| day | The day of the month (1-31). |
| hour | The hour of the day (0-23). |
| min | The minute (0-59). |
| sec | The second (0-59). |
| remain | The fractional part of seconds. |
| jsec | Pointer to store the calculated number of seconds since January 1, 2000, 12:00:00 UTC. |
The function calculates the time elapsed since January 1, 2000, 12:00:00 UTC, up to the specified time and includes any fractional seconds indicated by the "remain" parameter.
Definition at line 6485 of file jurassic.c.
| void timer | ( | const char * | name, |
| const char * | file, | ||
| const char * | func, | ||
| int | line, | ||
| int | mode | ||
| ) |
Simple wall-clock timer for runtime diagnostics.
Provides a lightweight timing utility based on omp_get_wtime() to measure wall-clock durations between marked code regions. The function supports up to ten concurrent nested timers.
| [in] | name | Name or label of the timed code section. |
| [in] | file | Source file name (usually __FILE__ macro). |
| [in] | func | Function name (usually __func__ macro). |
| [in] | line | Source line number (usually __LINE__ macro). |
| [in] | mode | Timer operation mode:
|
mode == 1 starts a new timer instance and stores its start time and corresponding source line.mode == 2 or mode == 3 is called, the elapsed wall-clock time (in seconds) is computed using: \[ \Delta t = t_{\text{now}} - t_{\text{start}} \]
and written to the log via the LOG macro.mode == 2 or 3 without a prior start causes an internal error.Definition at line 6516 of file jurassic.c.
| void write_atm | ( | const char * | dirname, |
| const char * | filename, | ||
| const ctl_t * | ctl, | ||
| const atm_t * | atm | ||
| ) |
Write atmospheric data to a file.
This function writes the atmospheric dataset stored in atm to the file specified by filename, optionally prefixed by dirname. The output format (ASCII or binary) is selected based on the atmospheric format flag ctl->atmfmt. The function creates the output file, delegates the writing process to the appropriate format-specific routine, and logs summary statistics of the written atmospheric data.
Supported output formats:
ctl->atmfmt == 1), written by write_atm_asc()ctl->atmfmt == 2), written by write_atm_bin()The function writes:
ctl->ng)ctl->nw)ctlAfter writing, the function logs minimum and maximum values of the written fields for verification and diagnostic purposes.
| dirname | Optional directory in which the output file will be created. If NULL, only filename is used. |
| filename | Name of the output file to be created and populated with atmospheric data. |
| ctl | Pointer to a control structure defining the output format and the sizes of gas, spectral, cloud, and surface parameter arrays. |
| atm | Pointer to the atmospheric data structure whose contents will be written. All required fields must be initialized and contain atm->np valid data points. |
ERRMSG if the file cannot be created or if an unsupported output format is requested.Definition at line 6554 of file jurassic.c.

Write atmospheric data to an ASCII file.
This function writes the contents of an atmospheric structure atm to an ASCII-formatted output stream out. A descriptive column header is generated first, documenting the meaning, units, and ordering of each data field. Atmospheric data points are then written line by line, with optional cloud and surface layer parameters appended if they are enabled in the control structure ctl.
The output columns include, in order:
ctl->ng) [ppv]ctl->nw) [km^-1]If cloud layer properties are enabled (ctl->ncl > 0), the following are added:
ctl->ncl) [km^-1]If surface layer properties are enabled (ctl->nsf > 0), the following are added:
ctl->nsf)| out | Pointer to an open output file stream where the ASCII data is written. |
| ctl | Pointer to a control structure defining the number of gases, spectral windows, and whether cloud or surface layer information should be included. |
| atm | Pointer to the atmospheric structure containing the data to be written. The function writes all atm->np data points. |
atm are properly allocated and populated. No validation of data ranges is performed here.Definition at line 6630 of file jurassic.c.
Write atmospheric data to a binary file.
This function writes the atmospheric dataset contained in atm to a binary file stream out. The output format is compact and includes a file header followed by the serialized atmospheric fields. The format is compatible with read_atm_bin(), ensuring that files written by this function can be read back without loss of information.
The binary file structure written is as follows:
"ATM1" (4 bytes)ctl->ng)ctl->nw)ctl->ncl)ctl->nsf)npnp containing:ctl->ng × np)ctl->nw × np)ctl:ctl->ncl)ctl->nsf)| out | Pointer to an already opened binary output file stream where the atmospheric data will be written. |
| ctl | Pointer to a control structure specifying the number of gases, spectral windows, and whether cloud or surface layer parameters must be included. |
| atm | Pointer to the atmospheric data structure containing values to be written. All arrays must be populated and atm->np must contain the number of valid atmospheric records. |
atm contents. It assumes that the memory layout matches expectations.read_atm_bin(); modifying either implementation requires updating the other accordingly.Definition at line 6691 of file jurassic.c.
Write atmospheric profile in RFM-compatible format.
Exports the current atmospheric state to a file formatted for use with the Reference Forward Model (RFM). The file includes altitude, pressure, temperature, and volume mixing ratio profiles for each active emitter.
| [in] | filename | Output file name for the RFM atmosphere file. |
| [in] | ctl | Pointer to the control structure defining active emitters. |
| [in] | atm | Pointer to the atmospheric profile to export. |
atm->np).*HGT, *PRE, etc.) and lists one value per line.Definition at line 6767 of file jurassic.c.
| void write_matrix | ( | const char * | dirname, |
| const char * | filename, | ||
| const ctl_t * | ctl, | ||
| const gsl_matrix * | matrix, | ||
| const atm_t * | atm, | ||
| const obs_t * | obs, | ||
| const char * | rowspace, | ||
| const char * | colspace, | ||
| const char * | sort | ||
| ) |
Write a fully annotated matrix (e.g., Jacobian or gain matrix) to file.
Outputs a numerical matrix along with detailed metadata describing the row and column spaces. Depending on configuration, the rows and columns may correspond to measurement or state variables.
| [in] | dirname | Output directory path (may be NULL). |
| [in] | filename | Output file name. |
| [in] | ctl | Pointer to control structure defining model setup and metadata. |
| [in] | matrix | Pointer to GSL matrix to write (e.g., Jacobian, kernel, or covariance). |
| [in] | atm | Pointer to atmospheric data structure (used when state-space indexing applies). |
| [in] | obs | Pointer to observation data structure (used when measurement-space indexing applies). |
| [in] | rowspace | Selects row labeling: "y" = measurement space, otherwise state space. |
| [in] | colspace | Selects column labeling: "y" = measurement space, otherwise state space. |
| [in] | sort | Determines writing order: "r" = row-major, otherwise column-major. |
ctl->write_matrix — output is skipped if disabled.sort to support row-major or column-major output.Definition at line 6805 of file jurassic.c.

| void write_obs | ( | const char * | dirname, |
| const char * | filename, | ||
| const ctl_t * | ctl, | ||
| const obs_t * | obs | ||
| ) |
Write observation data to an output file in ASCII or binary format.
This C function constructs the full output file path from the provided directory and filename, opens the file for writing, and exports the contents of the obs_t structure in either ASCII or binary format, depending on the observation format specified by ctl->obsfmt. The actual writing of formatted data is delegated to write_obs_asc() or write_obs_bin().
After writing, the function prints diagnostic information showing ranges of times, observer coordinates, view point coordinates, tangent point coordinates, radiance or brightness temperature values (depending on ctl->write_bbt), and transmittances. These diagnostics provide useful verification that the output data is valid and consistent.
| [in] | dirname | Optional directory path. If NULL, only filename is used. |
| [in] | filename | Name of the output observation file. |
| [in] | ctl | Control structure specifying output format, spectral channel configuration, and brightness-temperature mode. |
| [in] | obs | Observation structure containing the data to be written. |
ctl->obsfmt specifies an unsupported format.Definition at line 6983 of file jurassic.c.

Write observation data to an ASCII text file.
This C function writes the contents of the obs_t observation structure as human-readable ASCII text to the given output stream. It first prints a descriptive header that documents each column of the output format, including observation time, observer and view geometry, tangent point information, and spectral values. The number and meaning of spectral fields depend on ctl->nd and whether brightness temperature output is enabled via ctl->write_bbt.
The function then writes one line of data per ray path, including the base geometric information followed by radiance or brightness temperature values and transmittances for each spectral channel. Blank lines are inserted whenever the time stamp changes, providing visual separation of distinct observation groups.
| [in] | out | Output file stream opened in text mode. |
| [in] | ctl | Control structure specifying the number of spectral channels (nd), wavenumbers (nu), and output mode (write_bbt). |
| [in] | obs | Observation structure containing the data to be written. |
out is valid and writable. No attempt is made to reopen or validate the file stream.Definition at line 7065 of file jurassic.c.
Write observation data in binary format to an output file stream.
This C function serializes the contents of the obs_t structure into a compact binary format and writes it to the file stream provided via out. The binary format begins with a header consisting of a magic identifier ("OBS1") and the number of spectral channels (ctl->nd). This header is used by read_obs_bin() to validate compatibility when reading.
Following the header, the function writes the number of ray paths and then sequentially outputs arrays of observation metadata, geometric parameters, radiance or brightness temperature values, and transmittances. All values are written in native binary representation using the FWRITE() macro, which performs buffered writes and error checking.
| [out] | out | Output file stream opened in binary mode. |
| [in] | ctl | Control structure specifying the number of spectral channels (nd) and corresponding configuration parameters. |
| [in] | obs | Observation structure containing the data to be written. |
out is writable and already opened in binary mode. The function does not validate stream state.Definition at line 7113 of file jurassic.c.
| void write_shape | ( | const char * | filename, |
| const double * | x, | ||
| const double * | y, | ||
| const int | n | ||
| ) |
Write tabulated shape function data to a text file.
Exports a shape function (typically a weighting or field-of-view profile) defined by paired arrays of x and y values to an ASCII file.
| [in] | filename | Output file name. |
| [in] | x | Pointer to array of x-values (independent variable). |
| [in] | y | Pointer to array of y-values (dependent variable). |
| [in] | n | Number of data points to write. |
read_shape() for re-import.Definition at line 7173 of file jurassic.c.
| void write_stddev | ( | const char * | quantity, |
| const ret_t * | ret, | ||
| const ctl_t * | ctl, | ||
| const atm_t * | atm, | ||
| const gsl_matrix * | s | ||
| ) |
Write retrieval standard deviation profiles to disk.
Extracts the diagonal elements of a covariance matrix (a priori, posterior, or error covariance) to obtain the standard deviations of retrieved quantities and writes them as an atmospheric profile file.
| [in] | quantity | Name of the retrieved quantity (e.g., "apr", "pos", "err"), used to label the output file. |
| [in] | ret | Retrieval configuration structure (ret_t), providing the working directory for output files. |
| [in] | ctl | Global control structure (ctl_t) defining retrieval setup and quantities. |
| [in] | atm | Reference atmospheric state (atm_t) for spatial/geometric metadata. |
| [in] | s | Covariance matrix (gsl_matrix, n×n) from which standard deviations are derived (typically posterior covariance \(\mathbf{S}\)). |
This function performs the following operations:
atm) into an auxiliary structure (atm_aux) to preserve coordinate and geometric metadata.x2atm(), thereby mapping elements of the state vector to the corresponding atmospheric quantities.\[ \texttt{<ret->dir>/atm\_err\_<quantity>.tab} \]
using the standard JURASSIC atmospheric file format.atm_err_<quantity>.tab (e.g., atm_err_apr.tab, atm_err_pos.tab).ctl.s must be symmetric and positive-definite.x2atm) must correspond to the matrix ordering.Definition at line 7203 of file jurassic.c.

Write all emissivity lookup tables in the format specified by the control structure.
This function dispatches to one of three table writers depending on ctl->tblfmt:
1: ASCII tables written by write_tbl_asc()2: Binary tables written by write_tbl_bin()3: Per-gas binary tables written by write_tbl_gas()If an unknown format is given, the function aborts via ERRMSG().
| ctl | Control structure specifying table format, filenames, number of gases, number of frequencies, etc. |
| tbl | Fully populated lookup-table structure to be written. |
Definition at line 7236 of file jurassic.c.

Write all lookup tables in human-readable ASCII format.
For every gas index (ig) and frequency index (id), the function generates a file of the form:
<base>_<nu[id]>_<emitter[ig]>.tab
The ASCII file contains four columns:
1. pressure [hPa] 2. temperature [K] 3. column density [molecules/cm²] 4. emissivity [-]
Table dimensions are taken from the tbl_t structure.
Missing files cause the program to abort via ERRMSG().
| ctl | Control structure providing grid metadata and filename base. |
| tbl | Table data to be written. |
Definition at line 7259 of file jurassic.c.
Write all lookup tables in compact binary format.
For each gas index (ig) and frequency index (id), a binary file named
<base>_<nu[id]>_<emitter[ig]>.bin
is created. The format is:
| ctl | Control structure containing filename base and spectral grid. |
| tbl | Table data to be serialized. |
Definition at line 7304 of file jurassic.c.
Write lookup tables into per-gas binary table files with indexed blocks.
This function creates (if necessary) and updates gas-specific files of the form:
<base>_<emitter>.tbl
Each file contains:
For each frequency index (id), a block is appended (or overwritten) using write_tbl_gas_single(), which stores both the serialized table and its offset/size in the on-disk index.
| ctl | Control structure containing spectral grid, emitters, and filenames. |
| tbl | Table data from which individual frequency blocks are extracted. |
Definition at line 7359 of file jurassic.c.

| int write_tbl_gas_create | ( | const char * | path | ) |
Create a new per-gas table file with an empty index.
Writes the “GTL1” magic header, initializes the table count to zero, and creates a MAX_TABLES-sized index whose entries are zeroed.
The resulting file layout is:
magic[4] = "GTL1" ntables = 0 index[MAX_TABLES] (all zero)
| path | Path to the table file to create. |
Definition at line 7399 of file jurassic.c.
| int write_tbl_gas_single | ( | tbl_gas_t * | g, |
| const double | freq, | ||
| const tbl_t * | tbl, | ||
| const int | id, | ||
| const int | ig | ||
| ) |
Append or overwrite a single frequency-table block in a per-gas file.
Searches the in-memory index for an entry matching freq. If found, the corresponding block is updated. Otherwise a new entry is created (subject to the MAX_TABLES limit) and the block is appended to the end of the file.
The block format written is identical to the binary format used in write_tbl_bin():
The index entry is then updated with:
| g | Open gas-table handle obtained from read_tbl_gas_open(). |
| freq | Frequency associated with the table block. |
| tbl | Full lookup table from which one block is extracted. |
| id | Frequency index into tbl. |
| ig | Gas index into tbl. |
Definition at line 7430 of file jurassic.c.
Map retrieval state vector back to atmospheric structure.
Updates the atmospheric data structure (atm_t) from the contents of a retrieval state vector (x). This function performs the inverse transformation of atm2x(), assigning retrieved quantities such as pressure, temperature, gas volume mixing ratios, extinction, and cloud/surface parameters to the corresponding atmospheric fields.
| [in] | ctl | Pointer to control structure defining retrieval settings, vertical range limits, and active retrieval flags. |
| [in] | x | Pointer to retrieval state vector containing the updated values. |
| [out] | atm | Pointer to atmospheric data structure to be updated. |
ctl->ret*_zmin/ctl->ret*_zmax).x2atm_help() is called to sequentially read the next element from the state vector.atm2x() to ensure one-to-one correspondence between state vector indices and atmospheric fields.Quantities mapped:
p[zmin:zmax])t[zmin:zmax])q[ig][zmin:zmax])k[iw][zmin:zmax])clz, cldz, clk)sft, sfeps)ctl (e.g. ret_sft, ret_clk) are modified.n) advances automatically as each element is read.x2atm_help() abstracts sequential access to vector elements.Definition at line 7527 of file jurassic.c.

| void x2atm_help | ( | double * | value, |
| const gsl_vector * | x, | ||
| size_t * | n | ||
| ) |
Helper function to extract a single value from the retrieval state vector.
Retrieves the next element from the state vector x and assigns it to the provided scalar variable. This function is used by x2atm() to sequentially map the contents of the retrieval vector into the corresponding fields of the atmospheric structure.
| [out] | value | Pointer to the scalar variable to be updated. |
| [in] | x | Pointer to the retrieval state vector (gsl_vector). |
| [in,out] | n | Pointer to the current index in the state vector. The index is incremented after each extraction. |
atm2x() and x2atm().*n after reading a value, maintaining the correct position in the vector for subsequent calls.Definition at line 7577 of file jurassic.c.
Copy elements from the measurement vector y into the observation structure.
Decomposes the 1-D measurement vector y into its radiance components and writes them into the 2-D observation array obs->rad[id][ir], using the same ordering as produced by the corresponding forward model.
Only entries for which obs->rad[id][ir] is finite are updated. This allows missing or masked radiances to remain untouched in the observation structure.
| [in] | ctl | Control settings defining the number of detector channels (ctl->nd) and other retrieval configuration parameters. |
| [in] | y | Measurement vector containing radiances in forward-model order. |
| [out] | obs | Observation structure whose radiance array (obs->rad) is to be filled with values from y. |
The function loops over all ray paths (obs->nr) and detector channels (ctl->nd). For each pair (detector id, ray ir) where an existing value in obs->rad[id][ir] is finite, the next element from the measurement vector y is inserted. The counter m tracks progression through y.
This function is the inverse operation of the packing performed when constructing the measurement vector from an obs_t structure.
y must contain as many finite elements as the number of finite entries in obs->rad, in the same scanning order.Definition at line 7589 of file jurassic.c.