Functions
SpecGrid
-
void load_specgrid(const char *specgrid_file_name, SpecGrid *specgrid, Stat *stat)
Create a new
SpecGridby loading data from a specgrid file.- Parameters:
specgrid_file_name – Name of the specgrid file.
specgrid – Grid object.
stat – Status code (set to NULL if not required).
-
void unload_specgrid(SpecGrid specgrid)
Unload a
specgridgrid, freeing up memory.- Parameters:
specgrid – Grid object.
-
void get_specgrid_rank(SpecGrid specgrid, int *rank)
Get the rank (dimension) of the grid.
- Parameters:
specgrid – Grid object.
rank – Rank.
-
void get_specgrid_shape(SpecGrid specgrid, int shape[])
Get the shape of the grid.
- Parameters:
specgrid – Grid object.
rank – Shape.
-
void get_specgrid_lam_min(SpecGrid specgrid, double *lam_min)
Get the minimum wavelength of the grid.
- Parameters:
specgrid – Grid object.
lam_min – Minimum wavelength (Å).
-
void get_specgrid_lam_max(SpecGrid specgrid, double *lam_max)
Get the maximum wavelength of the grid.
- Parameters:
specgrid – Grid object.
lam_max – Maximum wavelength (Å).
-
void get_specgrid_cache_lam_min(SpecGrid specgrid, double *cache_lam_min)
Get the minimum wavelength of the grid cache.
- Parameters:
specgrid – Grid object.
cache_lam_min – Minimum wavelength (Å).
-
void get_specgrid_cache_lam_max(SpecGrid specgrid, double *cache_lam_max)
Get the maximum wavelength of the grid cache.
- Parameters:
specgrid – Grid object.
cache_lam_max – Maximum wavelength (Å).
-
void get_specgrid_cache_limit(SpecGrid specgrid, int *cache_limit)
Get the maximum memory usage of the grid cache.
- Parameters:
specgrid – Grid object.
cache_limit – Maximum memory usage (MB).
-
void get_specgrid_cache_usage(SpecGrid specgrid, int *cache_usage)
Get the current memory usage of the grid cache.
- Parameters:
specgrid – Grid object.
cache_usage – Current memory usage (MB).
-
void get_specgrid_axis_x_min(SpecGrid specgrid, int i, double *x_min)
Get the minimum value of the i’th grid axis.
- Parameters:
specgrid – Grid object.
i – Axis index (beginning at 0).
x_min – Minimum value.
-
void get_specgrid_axis_x_max(SpecGrid specgrid, int i, double *x_max)
Get the maximum value of the i’th grid axis.
- Parameters:
specgrid – Grid object.
i – Axis index (beginning at 0).
x_max – Maximum value.
-
void get_specgrid_axis_label(SpecGrid specgrid, int i, char *label)
Get the label of the i’th grid axis.
- Parameters:
photgrid – Grid object.
i – Index of the label (beginning at 0).
axis_label – Buffer to store axis label (at least 17 bytes, to accomodate label plus null terminator).
-
void set_specgrid_cache_lam_min(SpecGrid specgrid, double cache_lam_min, Stat *stat)
Set the minimum wavelength of the grid cache.
- Parameters:
specgrid – Grid object.
cache_lam_min – Minimum wavelength (Å).
stat – Status code (set to NULL if not required).
-
void set_specgrid_cache_lam_max(SpecGrid specgrid, double cache_lam_max, Stat *stat)
Set the maximum wavelength of the grid cache.
- Parameters:
specgrid – Grid object.
cache_lam_max – Maximum wavelength (Å).
stat – Status code (set to NULL if not required).
-
void set_specgrid_cache_limit(SpecGrid specgrid, int cache_limit, Stat *stat)
Set the maximum memory usage of the grid cache.
- Parameters:
specgrid – Grid object.
cache_limit – Maximum memory usage (MB).
stat – Status code (set to NULL if not required).
-
void interp_specgrid_intensity(SpecGrid specgrid, int n, int r, double x_vec[], double mu, double z, double lam[], double I[], Stat *stat, bool deriv_vec[], int *order)
Interpolate the spectroscopic specific intensity for a photospheric element.
- Parameters:
specgrid – Grid object.
n – Number of wavelength points.
r – Number of photospheric parameters.
x_vec[r] – Photospheric parameter values.
mu – Cosine of angle of emergence relative to element normal.
z – Redshift of element relative to observer’s frame.
lam[n] – Wavelength abscissa (Å) in observer’s frame.
I[n-1] – Spectroscopic specific intensity (erg/cm^2/s/Å/sr) in bins delineated by lam.
stat – Status code (set to NULL if not required).
deriv_vec[r] – Derivative flags (set to NULL if not required).
order – Interpolation order (1 or 3; set to NULL if not required).
-
void interp_specgrid_E_moment(SpecGrid specgrid, int n, int r, double x_vec[], int k, double z, double lam[], double E[], Stat *stat, bool deriv_vec[], int *order)
Interpolate the spectroscopic intensity E-moment for a photospheric element.
- Parameters:
specgrid – Grid object.
n – Number of wavelength points.
r – Number of photospheric parameters.
x_vec[r] – Photospheric parameter values.
k – Degree of moment.
z – Redshift of element relative to observer’s frame.
lam[n] – Wavelength abscissa (Å) in observer’s frame.
D[n-1] – Spectroscopic intensity E-moment (erg/cm^2/s/Å) in bins delineated by lam.
stat – Status code (set to NULL if not required).
deriv_vec[r] – Derivative flags (set to NULL if not required).
order – Interpolation order (1 or 3; set to NULL if not required).
-
void interp_specgrid_P_moment(SpecGrid specgrid, int n, int r, double x_vec[], int l, double z, double lam[], double P[], Stat *stat, bool deriv_vec[], int *order)
Interpolate the spectroscopic intensity P-moment.
- Parameters:
specgrid – Grid object.
n – Number of wavelength points.
r – Number of photospheric parameters.
x_vec[r] – Photospheric parameter values.
l – Harmonic degree of moment.
z – Redshift of element relative to observer’s frame.
lam[n] – Wavelength abscissa (Å) in observer’s frame.
D[n-1] – Spectroscopic intensity P-moment (erg/cm^2/s/Å) in bins delineated by lam.
stat – Status code (set to NULL if not required).
deriv_vec[r] – Derivative flags (set to NULL if not required).
order – Interpolation order (1 or 3; set to NULL if not required).
-
void interp_specgrid_irradiance(SpecGrid specgrid, int n, int m, int r, double x_vec[], double mu[], double dOmega[], double z, double lam[], double F[], Stat *stat, bool deriv_vec[], int *order)
Interpolate the spectroscopic irradiance for an object composed of multiple photospheric elements.
- Parameters:
specgrid – Grid object.
n – Number of wavelength points.
m – Number of elements.
r – Number of photospheric parameters.
x_vec[m,r] – Photospheric parameter values.
mu[m] – Cosines of angles of emergence relative to element normals.
dOmega[m] – Solid angles (sr) of elements.
z[m] – Redshifts of elements relative to observer’s frame.
lam[n] – Wavelength abscissa (Å) in observer’s frame.
F[n-1] – Spectroscopic irradiance (erg/cm^2/s/Å) in bins delineated by lam.
stat – Status code (set to NULL if not required).
deriv_vec[r] – Derivative flags (set to NULL if not required).
order – Interpolation order (1 or 3; set to NULL if not required).
-
void interp_specgrid_flux(SpecGrid specgrid, int n, int r, double x_vec[], double z, double lam[], double F[], Stat *stat, bool deriv_vec[], int *order)
Interpolate the spectroscopic flux for a photospheric element.
- Parameters:
specgrid – Grid object.
n – Number of wavelength points.
r – Number of photospheric parameters.
x_vec[r] – Photospheric parameter values.
z – Redshift of element relative to observer’s frame.
lam[n] – Wavelength abscissa (Å) in observer’s frame.
F[n-1] – Spectroscopic flux (erg/cm^2/s/Å) in bins delineated by lam.
stat – Status code (set to NULL if not required).
deriv_vec[r] – Derivative flags (set to NULL if not required).
order – Interpolation order (1 or 3; set to NULL if not required).
-
void adjust_specgrid_x_vec(SpecGrid specgrid, int r, double x_vec[], double dx_vec[], double x_adj[], Stat *stat)
Adjust photospheric parameters in a specified direction, until they fall within a valid part of the grid.
- Parameters:
specgrid – Grid object.
r – Number of photospheric parameters.
x_vec[r] – Photospheric parameter values.
dx_vec[r] – Photospheric parameter adjustment direction. The overall scaling is unimportant, but at least one element must be non-zero.
x_adj[r] – Adjusted photospheric parameter values.
stat – Status code (set to NULL if not required).
PhotGrid
-
void load_photgrid(const char *photgrid_file_name, PhotGrid *photgrid, int *stat)
Create a new
PhotGridby loading data from a photgrid file.- Parameters:
photgrid_file_name – Name of the photgrid file.
photgrid – Grid object.
stat – Status code (set to NULL if not required).
-
void load_photgrid_from_specgrid(const char *specgrid_file_name, const char *passband_file_name, PhotGrid *photgrid, int *stat)
Create a new
PhotGridby loading data from a specgrid file, and convolving on-the-fly with a response function loaded from a passband file.
-
void unload_photgrid(PhotGrid photgrid)
Unload a photometric grid, freeing up memory.
- Parameters:
photgrid – Grid object.
-
void get_photgrid_rank(Photgrid photgrid, int *rank)
Get the rank (dimension) of the grid.
- Parameters:
photgrid – Grid object.
rank – Rank.
-
void get_photgrid_shape(Photgrid photgrid, int shaope[])
Get the shape of the grid.
- Parameters:
photgrid – Grid object.
rank – Shape.
-
void get_photgrid_cache_usage(Photgrid photgrid, int *cache_usage)
Get the current memory usage of the grid cache.
- Parameters:
photgrid – Grid object.
cache_usage – Current memory usage (MB).
-
void get_photgrid_cache_limit(Photgrid photgrid, int *cache_limit)
Get the maximum memory usage of the grid cache.
- Parameters:
photgrid – Grid object.
cache_limit – Maximum memory usage (MB).
-
void get_photgrid_axis_x_min(Photgrid photgrid, int i, double *x_min)
Get the minimum value of the i’th grid axis.
- Parameters:
photgrid – Grid object.
i – Axis index (beginning at 0).
x_min – Minimum value.
-
void get_photgrid_axis_x_max(Photgrid photgrid, int i, double *x_max)
Get the maximum value of the i’th grid axis.
- Parameters:
photgrid – Grid object.
i – Axis index (beginning at 0).
x_max – Maximum value.
-
void get_photgrid_axis_label(SpecGrid specgrid, int i, char *label)
Get the label of the i’th grid axis.
- Parameters:
photgrid – Grid object.
i – Index of the label (beginning at 0).
axis_label – Buffer to store axis label (at least 17 bytes, to accomodate label plus null terminator).
-
void set_photgrid_cache_limit(Photgrid photgrid, int cache_limit, int *stat)
Set the maximum memory usage of the grid cache.
- Parameters:
photgrid – Grid object.
cache_limit – Maximum memory usage (MB).
stat – Status code (set to NULL if not required).
-
void interp_photgrid_intensity(PhotGrid photgrid, int r, double x_vec[], double mu, double *I, int *stat, bool deriv_vec[], int *order)
Interpolate the photometric intensity for a photospheric element, normalized by the zero-point flux.
- Parameters:
photgrid – Grid object.
r – Number of photospheric parameters.
x_vec[r] – Photospheric parameter values.
mu – Cosine of angle of emergence relative to element normal.
I – Photometric specific intensity (/sr).
stat – Status code (set to NULL if not required).
deriv_vec[r] – Derivative flags (set to NULL if not required).
order – Interpolation order (1 or 3; set to NULL if not required).
-
void interp_photgrid_E_moment(PhotGrid photgrid, int r, double x_vec[], int k, double *E, int *stat, bool deriv_vec[], int *order)
Interpolate the photometric intensity E-moment for a photospheric element, normalized by the zero-point flux.
- Parameters:
photgrid – Grid object.
r – Number of photospheric parameters.
x_vec[r] – Photospheric parameter values.
k – Degree of moment.
D – Photometric intensity E-moment (dimensionless).
stat – Status code (set to NULL if not required).
deriv_vec[r] – Derivative flags (set to NULL if not required).
order – Interpolation order (1 or 3; set to NULL if not required).
-
void interp_photgrid_P_moment(PhotGrid photgrid, int r, double x_vec[], int l, double *P, int *stat, bool deriv_vec[], int *order)
Interpolate the photometric intensity P-moment for a photospheric element, normalized by the zero-point flux.
- param photgrid:
Grid object. :param r: Number of photospheric parameters. :param x_vec[r]: Photospheric parameter values. :param l: Harmonic degree of moment. :param D: Photometric intensity P-moment (dimensionless). :param stat: Status code (set to NULL if not required). :param deriv_vec[r]: Derivative flags (set to NULL if not required). :param order: Interpolation order (1 or 3; set to NULL if not required).
-
void interp_photgrid_irradiance(PhotGrid photgrid, int m, int r, double x_vec[], double mu[], double dOmega[], double *F, int *stat, bool deriv_vec[], int *order)
Interpolate the photometric irradiance for an object composed of multiple photospheric elements, normalized by the zero-point flux.
- Parameters:
photgrid – Grid object.
m – Number of elements.
r – Number of photospheric parameters.
x_vec[m,r] – Photospheric parameter values.
mu[m] – Cosines of angles of emergence relative to element normals.
dOmega[m] – Solid angles (sr) of elements.
I – Photometric irradiance (dimensionless).
stat – Status code (set to NULL if not required).
deriv_vec[r] – Derivative flags (set to NULL if not required).
order – Interpolation order (1 or 3; set to NULL if not required).
-
void interp_photgrid_flux(PhotGrid photgrid, int r, double x_vec[], double *F, int *stat, bool deriv_vec[], int *order)
Interpolate the photometric flux for a photospheric element, normalized by the zero-point flux.
- Parameters:
PhotGrid – Grid object.
r – Number of photospheric parameters.
x_vec[r] – Photospheric parameter values.
F – Photometric flux (dimensionless).
stat – Status code (set to NULL if not required).
deriv_vec[r] – Derivative flags (set to NULL if not required).
order – Interpolation order (1 or 3; set to NULL if not required).
-
void adjust_photgrid_x_vec(PhotGrid photgrid, int r, double x_vec[], double dx_vec[], double x_adj[], Stat *stat)
Adjust photospheric parameters in a specified direction, until they fall within a valid part of the grid.
- Parameters:
photgrid – Grid object.
r – Number of photospheric parameters.
x_vec[r] – Photospheric parameter values.
dx_vec[r] – Photospheric parameter adjustment direction. The overall scaling is unimportant, but at least one element must be non-zero.
x_adj[r] – Adjusted photospheric parameter values.
stat – Status code (set to NULL if not required).
Other
-
void get_msg_version(char *version)
Get the version of the MSG library.
- Parameters:
version – Buffer to store version (recommended at least 10 bytes).