MSG
Multidimensional Spectral Grids (MSG) is an open-source software package that synthesizes stellar spectra and photometric colors via interpolation in pre-calculated grids. MSG is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, version 3.
Preliminaries
Why use MSG?
Synthesizing stellar spectra from first principles is a complicated endeavor, requiring a detailed understanding of radiative transfer and atomic physics, together with significant computational resources. Therefore, in most circumstances its better to use one of the many pre-calculated grids of spectra published in the astrophysical literature (see, e.g., Lanz & Hubeny, 2003; Lanz & Hubeny, 2007; Kirby, 2011; de Laverny et al., 2012; Husser et al., 2013; Allende Prieto et al., 2018; Chiavassa et al., 2018; Zsargó et al., 2020). However, even with these grids a significant obstacle remains: when photospheric parameters fall between the grid nodes, some kind of interpolation is necessary in order to evaluate a spectrum.
MSG is designed to solve this problem. It’s not the first software package that offers stellar spectral interpolation (see, e.g., FERRE, Starfish and stsynphot); however, with spectral interpolation as its sole focus, it offers a combination of features unmatched by other packages:
scalability — MSG handles grids that are much larger (on disk) than available computer memory.
extensibility — MSG handles grids with an arbitrary number of dimensions.
portability — MSG is platform-agnostic and provide APIs for the programming languages (Fortran, C, Python) most commonly used in Astronomy.
performance — MSG provides smooth and accurate interpolates with minimal computational cost.
robustness — MSG gracefully handles missing data caused by holes and/or ragged boundaries in the grid.
Together, these features mean that MSG is flexible and powerful while remaining straightforward to use: it’s the perfect condiment to add flavor to your science!.
Obtaining MSG
The source code for MSG is hosted in the rhdtownsend/msg git repository on GitHub. Instructions for downloading and installing the software can be found in the Quick Start chapter.
Citing MSG
If you use MSG in your research, please cite the following papers:
Townsend R. H. D., Lopez A., 2022, Journal of Open-Source Software, in preparation
Be sure to also cite the source of the grid data that you’re using with MSG. For instance, if you’re working with one of the CAP18 grids, you should cite Allende Prieto et al. (2018).
Development Team
MSG remains under active development by the following team:
Rich Townsend (University of Wisconsin-Madison); project leader
Acknowledgments
MSG has been developed with financial support from the following grants:
NSF awards ACI-1663696 and AST-1716436;
NASA award 80NSSC20K0515.
Quick Start
To get started with MSG, follow these simple steps:
Install the MESA Software Development Kit.
Download the MSG source code.
Unpack the source code using the tar utility.
Set the
MSG_DIR
environment variable to point to the newly created source directory.Compile MSG using the command make -C $MSG_DIR.
Test MSG using the command make -C $MSG_DIR test.
Install the
pymsg
Python module using the command pip install $MSG_DIR/python.
(the last step can be skipped if you don’t plan to use the Python interface).
For a more in-depth installation guide that covers alternative use-cases, refer to the Installation chapter. If the code doesn’t compile properly, consult the Troubleshooting chapter. Otherwise, proceed to the next chapter where you’ll undertake your first MSG calculations.
Walkthroughs
The best way to familiarize yourself with MSG is use it. This chapter provides example walkthroughs in Python, Fortran, and C, each focused on calculating the spectrum and colors of Sirius (\(\alpha\) Canis Majoris A).
Python Walkthrough
This section walks through using the MSG Python interface
to evaluate spectra and photometric colors for a model of Sirius. The code fragments are presented as a
sequence of Jupyter notebook cells, but pasting all of them
into a Python script should work also. Alternatively, the entire notebook (including this text) is available in the source distribution at $MSG_DIR/docs/source/user-guide/walkthroughs/python-walkthrough.ipynb
.
Preparation
Before starting Jupyter, make sure that you have installed the pymsg
package as described in the Quick Start and/or Installation chapters. Also, ensure that the NumPy, Matplotlib and
Astropy packages are installed on your
system. Finally, confirm that the MSG_DIR
environment
variable is set.
Initialization
Launch Jupyter, create a new Python3 notebook, and initialize the notebook with the following statements:
[1]:
# Import standard modules
import os
import numpy as np
import astropy.constants as con
import astropy.units as unt
import matplotlib.pyplot as plt
# Import pymsg
import pymsg
# Set plot parameters
%matplotlib inline
plt.rcParams.update({'font.size': 16})
Loading & Inspecting the Demo Grid
For demonstration purposes MSG includes a low-resolution specgrid (spectroscopic grid) file, $MSG_DIR/data/grids/sg-demo.h5
1. Load this demo grid into memory2 by creating a
new pymsg.SpecGrid
object:
[2]:
# Load the SpecGrid
MSG_DIR = os.environ['MSG_DIR']
GRID_DIR = os.path.join(MSG_DIR, 'data', 'grids')
specgrid_file_name = os.path.join(GRID_DIR, 'sg-demo.h5')
specgrid = pymsg.SpecGrid(specgrid_file_name)
This object has a number of (read-only) properties that inform us about the parameter space spanned by the grid:
[3]:
# Inspect grid parameters
print('Grid parameters:')
for label in specgrid.axis_labels:
print(f' {label} [{specgrid.axis_x_min[label]} -> {specgrid.axis_x_max[label]}]')
print(f' lam [{specgrid.lam_min} -> {specgrid.lam_max}]')
Grid parameters:
Teff [3500.0 -> 50000.0]
log(g) [0.0 -> 5.0]
lam [2999.9999999999977 -> 9003.490078514782]
Here, Teff
and log(g)
correspond (respectively) to the
\(\Teff/\kelvin\) and \(\log_{10}(g/\cm\,\second^{-2})\) photospheric parameters, while
lam
is wavelength \(\lambda/\angstrom\).
Plotting the Flux
With the grid loaded, let’s evaluate and plot a flux spectrum for Sirius A. First, store the photospheric parameters for the star in a dict:
[4]:
# Set Sirius A's fundamental parameters
M = 2.02 *con.M_sun
R = 1.71 * con.R_sun
L = 25.4 * con.L_sun
p = 0.379
d = 1/p * con.pc
# Set photospheric parameters dict
Teff = (L/(4*np.pi*R**2*con.sigma_sb))**0.25/unt.K
logg = np.log10(con.G*M/R**2/(unt.cm/unt.s**2))
x = {'Teff': Teff, 'log(g)': logg}
print(Teff, logg)
9909.16258994498 4.277426775965449
(the mass, radius and luminosity data are taken from Liebert et al., 2005, while the parallax comes from Simbad). Then set up a uniformly-spaced wavelength abcissa for a spectrum spanning the visible range, \(3\,000\,\angstrom\) to \(7\,000\,\angstrom\).
[5]:
# Set up the wavelength abscissa
lam_min = 3000.
lam_max = 7000.
lam = np.linspace(lam_min, lam_max, 501)
lam_c = 0.5*(lam[1:] + lam[:-1])
Here, the array lam
defines the boundaries of 500 wavelength bins
\(\{\lambda_{i},\lambda_{i+1}\}\) (\(i=1,\ldots,500\)) and the
array lam_c
stores the central wavelength of each bin.
With all our parameters defined, evaluate the observed flux using a
call to the pymsg.SpecGrid.flux()
method, and then plot it:
[6]:
# Evaluate the observed flux
F_obs = (R/d)**2*specgrid.flux(x, lam)
# Plot
plt.figure(figsize=[8,8])
plt.plot(lam_c, F_obs)
plt.xlabel(r'$\lambda ({\AA})$')
plt.ylabel(r'$F^{\mathrm{o}}_{\lambda}\ ({\rm erg\,cm^{-2}\,s^{-1}}\,\AA^{-1})$')
[6]:
Text(0, 0.5, '$F^{\\mathrm{o}}_{\\lambda}\\ ({\\rm erg\\,cm^{-2}\\,s^{-1}}\\,\\AA^{-1})$')

(See equation (7) of the MSG Fundamentals chapter to understand the origin of the expression for F_obs
). The plot looks about right — we can clearly see the Balmer series, starting with H\(\alpha\) at \(6563\,\angstrom\).
Plotting the Intensity
Sometimes we want to know the specific intensity of the radiation
emerging from each element of a star’s surface; an example might be when we’re
modeling a rotationally distorted star, whose surface geometry and properties vary with latitude. For this, we can use the pymsg.SpecGrid.intensity()
method.
Here’s a demonstration of this method in action, plotting the specific intensity across the H\(\alpha\) line profile for ten different values of the emergence direction cosine \(\mu=0.1,0.2,\ldots,1.0\).
[7]:
# Set up the wavelength abscissa
lam_min = 6300.
lam_max = 6800.
lam = np.linspace(lam_min, lam_max, 100)
lam_c = 0.5*(lam[1:] + lam[:-1])
# Loop over mu
plt.figure(figsize=[8,8])
for mu in np.linspace(1.0, 0.1, 10):
# Evaluate the intensity
I = specgrid.intensity(x, mu, lam)
# Plot
if mu==0.1 or mu==1.0:
label=r'$\mu={:3.1f}$'.format(mu)
else:
label=None
plt.plot(lam_c, I, label=label)
plt.xlabel(r'$\lambda ({\AA})$')
plt.ylabel(r'$I_{\lambda}\ ({\rm erg\,cm^{-2}\,s^{-1}}\,\AA^{-1}\,sr^{-1})$')
plt.legend()
[7]:
<matplotlib.legend.Legend at 0x7ff643465670>

We can clearly see that limb-darkening in the line core is much weaker than in the continuum — exactly what we expect from such a strong line.
Evaluating Magnitudes & Colors
As a final step in this walkthrough, let’s evaluate the magnitude and
colors of Sirius A in the Johnson photometric system. We can do this by
creating a new pymsg.PhotGrid
object for each passband:
[8]:
# Load the PhotGrids
PASS_DIR = os.path.join(MSG_DIR, 'data', 'passbands')
filters = ['U', 'B', 'V']
photgrids = {}
for filter in filters:
passband_file_name = os.path.join(PASS_DIR, f'pb-Generic-Johnson.{filter}-Vega.h5')
photgrids[filter] = pymsg.PhotGrid(specgrid_file_name, passband_file_name)
(for convenience, we store the pymsg.PhotGrid
objects in a
dict, keyed by filter name). In the calls to the object constructor
pymsg.PhotGrid()
, the first argument is the name of a
specgrid file (here, the demo grid), and the second argument is the
name of a passband file; a limited set of these files is
provided in the $MSG_DIR/data/passbands
subdirectory 3.
The observed fluxes of Sirius A are then found using the pymsg.PhotGrid.flux()
method:
[9]:
# Evaluate the observed fluxes (each normalized to the passband
# zero-point flux)
F_obs = {}
for filter in filters:
F_obs[filter] = (R/d)**2*photgrids[filter].flux(x)
To convert these into apparent magnitudes, we apply Pogson’s logarithmic formula:
[10]:
# Evaluate apparent magnitudes and print out magnitude & colors
mags = {}
for filter in filters:
mags[filter] = -2.5*np.log10(F_obs[filter])
print(f"V={mags['V']}, U-B={mags['U']-mags['B']}, B-V={mags['B']-mags['V']}")
V=-1.4466779855804992, U-B=-0.05964681835037444, B-V=0.001121934559611848
Reassuringly, the resulting values are within 10 millimag of Sirius A’s apparent V-band magnitude and colors, as given by Simbad.
Footnotes
- 1
Larger grids can be downloaded separately; see the Grid Files appendix.
- 2
Behind the scenes, the data in the specgrid file are loaded on demand; see the Data Caching chapter for further details.
- 3
Passband definition files for other instruments/photometric systems can be downloaded separately; see the Passband Files appendix.
Fortran Walkthrough
This section reprises the steps of the Python Walkthrough — evaluating spectra and photometric colors for Sirius A — but now using the MSG Fortran interface.
Preparation
In your working directory, create a new file
fortran-walkthrough.f90
with the following source code:
program fortran_walkthrough
! Uses
use forum_m
use fmsg_m
! No implicit typing
implicit none
! Parameters
real(RD), parameter :: lam_min = 3000._RD
real(RD), parameter :: lam_max = 7000._RD
integer, parameter :: n_lam = 501
! Variables
type(specgrid_t) :: specgrid
type(axis_t) :: axis
character(LABEL_LEN) :: label
integer :: i
real(RD) :: R2_d2
real(RD) :: lam(n_lam)
real(RD) :: lam_c(n_lam-1)
real(RD) :: x_vec(2)
real(RD) :: F(n_lam-1)
real(RD) :: F_obs(n_lam-1)
integer :: unit
type(photgrid_t) :: photgrid_U
type(photgrid_t) :: photgrid_B
type(photgrid_t) :: photgrid_V
real(RD) :: F_U
real(RD) :: F_B
real(RD) :: F_V
real(RD) :: F_U_obs
real(RD) :: F_B_obs
real(RD) :: F_V_obs
real(RD) :: U
real(RD) :: B
real(RD) :: V
! Load the specgrid
call load_specgrid('sg-demo.h5', specgrid)
! Set photospheric parameters to correspond to Sirius A
do i = 1, 2
call specgrid%get_axis(i, axis)
call axis%get_label(label)
select case(label)
case('log(g)')
x_vec(i) = 4.2774_RD
case('Teff')
x_vec(i) = 9909.2_RD
case default
stop 'unrecognized axis label'
end select
end do
! Set the dilution factor R2_d2 = R**2/d**2, where R is Sirius A's
! radius and d its distance
R2_d2 = 2.1351E-16_RD
! Set up the wavelength abscissa
lam = [((lam_min*(n_lam-i) + lam_max*(i-1))/(n_lam-1), i=1,n_lam)]
lam_c = 0.5_RD*(lam(:n_lam-1) + lam(2:))
! Evaluate the observed flux
call specgrid%interp_flux(x_vec, lam, F)
F_obs = R2_d2*F
! Write it to a file
open(NEWUNIT=unit, FILE='spectrum.dat', STATUS='REPLACE')
do i = 1, n_lam-1
write(unit, *) lam_c(i), F_obs(i)
end do
close(unit)
! Load the photgrids
call load_photgrid_from_specgrid('sg-demo.h5', 'pb-Generic-Johnson.U-Vega.h5', photgrid_U)
call load_photgrid_from_specgrid('sg-demo.h5', 'pb-Generic-Johnson.B-Vega.h5', photgrid_B)
call load_photgrid_from_specgrid('sg-demo.h5', 'pb-Generic-Johnson.V-Vega.h5', photgrid_V)
! Evaluate fluxes
call photgrid_U%interp_flux(x_vec, F_U)
call photgrid_B%interp_flux(x_vec, F_B)
call photgrid_V%interp_flux(x_vec, F_V)
F_U_obs = R2_d2*F_U
F_B_obs = R2_d2*F_B
F_V_obs = R2_d2*F_V
! Evaluate apparent magnitudes
U = -2.5_RD*LOG10(F_U_obs)
B = -2.5_RD*LOG10(F_B_obs)
V = -2.5_RD*LOG10(F_V_obs)
print *, ' V=', V
print *, 'U-B=', U-B
print *, 'B-V=', B-V
! Finish
end program fortran_walkthrough
A few brief comments on the code:
The
use forum_m
statement provides access to the Fortran Utility Module (ForUM). For the purposes of the demo program, this module defines the RD kind type parameter for double precision real variables.The
use fmsg_m
statement provides access to the MSG Fortran interface. Primarily, this interface serves to define thespecgrid_t
andphotgrid_t
datatypes.Because Fortran doesn’t have
dict
datatypes, the photospheric parameters must be passed to MSG as a plain array (here, stored in the variablex_vec
). Aselect case
construct is used to make sure the correct values are stored in each array element.
Compiling
The next step is to compile the demo program. Make sure the
MSG_DIR
environment variable is set, as described in the
Quick Start chapter. Then, enter the following on the command line:
gfortran -o fortran-walkthrough fortran-walkthrough.f90 -I$MSG_DIR/include `$MSG_DIR/scripts/fmsg_link`
The -I$MSG_DIR/include
option tells the compiler where to find
the module definition (.mod
) files, while the
$MSG_DIR/scripts/fmsg_link
clause (note the enclosing
backticks) runs a link script that returns the compiler flags
necessary to link the program against the appropriate libraries.
Running
To run the code, first create a symbolic link to the demo grid:
ln -s $MSG_DIR/data/grids/sg-demo.h5
Then, execute the command
./fortran-walkthrough
The code will create a file spectrum.dat
containing the flux
spectrum for Sirius A (as an ASCII table), and print out the
apparent V-band magnitude and colors of the star.
C Walkthrough
This section reprises the steps of the Python Walkthrough — evaluating spectra and photometric colors for Sirius A — but now using the MSG C interface.
Preparation
In your working directory, create a new file
c-walkthrough.c
with the following source code:
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <string.h>
#include "cmsg.h"
#define N_LAM 501
#define LAM_MIN 3000.
#define LAM_MAX 7000.
int main(int argc, char *argv[]) {
SpecGrid specgrid;
PhotGrid photgrid_U;
PhotGrid photgrid_B;
PhotGrid photgrid_V;
char label[17];
double x_vec[2];
double R2_d2;
double lam[N_LAM];
double lam_c[N_LAM-1];
double F[N_LAM-1];
double F_obs[N_LAM-1];
double F_U;
double F_B;
double F_V;
double F_U_obs;
double F_B_obs;
double F_V_obs;
double U;
double B;
double V;
FILE *fptr;
// Load the specgrid
load_specgrid("sg-demo.h5", &specgrid, NULL);
// Set photospheric parameters to correspond to Sirius A
for(int i=0; i < 2; i++) {
get_specgrid_axis_label(specgrid, i, label);
if (strcmp(label, "log(g)") == 0) {
x_vec[i] = 4.2774;
}
else if (strcmp(label, "Teff") == 0) {
x_vec[i] = 9909.2;
}
else {
printf("invalid label\n");
exit(1);
}
}
// Set the dilution factor R2_d2 = R**2/d**2, where R is Sirius A's
// radius and d its distance
R2_d2 = 2.1351E-16;
// Set up the wavelength abscissa
for(int i=0; i < N_LAM; i++) {
lam[i] = (LAM_MIN*(N_LAM-1-i) + LAM_MAX*i)/(N_LAM-1);
}
for(int i=0; i < N_LAM-1; i++) {
lam_c[i] = 0.5*(lam[i] + lam[i+1]);
}
// Evaluate the flux
interp_specgrid_flux(specgrid, x_vec, N_LAM, lam, F, NULL, NULL);
for(int i=0; i < N_LAM-1; i++) {
F_obs[i] = R2_d2*F[i];
}
// Write it to a file
fptr = fopen("spectrum.dat", "w");
for(int i=0; i < N_LAM-1; i++) {
fprintf(fptr, "%.17e %.17e\n", lam_c[i], F_obs[i]);
}
fclose(fptr);
// Load the photgrids
load_photgrid_from_specgrid("sg-demo.h5", "pb-Generic-Johnson.U-Vega.h5", &photgrid_U, NULL);
load_photgrid_from_specgrid("sg-demo.h5", "pb-Generic-Johnson.B-Vega.h5", &photgrid_B, NULL);
load_photgrid_from_specgrid("sg-demo.h5", "pb-Generic-Johnson.V-Vega.h5", &photgrid_V, NULL);
// Evaluate fluxes
interp_photgrid_flux(photgrid_U, x_vec, &F_U, NULL, NULL);
interp_photgrid_flux(photgrid_B, x_vec, &F_B, NULL, NULL);
interp_photgrid_flux(photgrid_V, x_vec, &F_V, NULL, NULL);
F_U_obs = R2_d2*F_U;
F_B_obs = R2_d2*F_B;
F_V_obs = R2_d2*F_V;
// Evaluate apparent magnitudes
U = -2.5*log10(F_U_obs);
B = -2.5*log10(F_B_obs);
V = -2.5*log10(F_V_obs);
printf(" V= %.17e\n", V);
printf("U-B= %.17e\n", U-B);
printf("B-V= %.17e\n", B-V);
// Clean up
unload_specgrid(specgrid);
unload_photgrid(photgrid_U);
unload_photgrid(photgrid_B);
unload_photgrid(photgrid_V);
// Finish
exit(0);
}
A few brief comments on the code:
The
#include "cmsg.h"
statement includes the header definitions for the MSG C interface.Because C doesn’t have
dict
datatypes, the photospheric parameters must be passed to MSG as a plain array (here, stored in the variablex_vec
). A loop withstrcmp()
calls is used to make sure the correct values are stored in each array element.Many of the calls to MSG routines (e.g.,
load_specgrid()
,interp_specgrid_flux()
) containNULL
trailing arguments; these correspond to omitted optional arguments.
Compiling
The next step is to compile the demo program. Make sure the
MSG_DIR
environment variable is set, as described in the
Quick Start chapter. Then, enter the following on the command line:
gcc -o c-walkthrough c-walkthrough.c -I$MSG_DIR/include `$MSG_DIR/scripts/cmsg_link`
The -I$MSG_DIR/include
option tells the compiler where to find
the header file, while the $MSG_DIR/scripts/cmsg_link
clause
(note the enclosing backticks) runs a link script that returns the
compiler flags necessary to link the program against the appropriate
libraries.
Running
To run the code, first create a symbolic link to the demo grid:
ln -s $MSG_DIR/data/grids/sg-demo.h5
Then, execute the command
./c-walkthrough
The code will create a file spectrum.dat
containing the flux
spectrum for Sirius A (as an ASCII table), and print out the
apparent V-band magnitude and colors of the star.
MSG Fundamentals
This chapter expands on the Walkthroughs, by describing in detail how MSG evaluates stellar spectra and photometric colors. It also aims to clarify the variety of different concepts that are often lumped together under the name ‘spectrum’.
Elemental Spectra
The radiation emitted by a small element of a star’s surface is most completely characterized by the specific intensity \(\intsy(\lambda; \vshat; \vx)\). This quantity is defined such that energy passing through the element into a solid angle \(\diff{\Omega}\) oriented along the unit direction vector \(\vshat\), and within wavelength interval \([\lambda, \lambda+\diff{\lambda}]\) and time interval \(\diff{t}\), is
Here, \(\diff{\area}\) is the area of the surface element, \(\vnhat\) is the unit surface normal vector, and \(\vx\) is a vector specifying the photospheric parameters of the element — for instance, its effective temperature \(T_{\rm eff}\) and gravity \(g\). Integrating this equation over all solid angles yields the net energy passing through the element in the wavelength and time intervals:
where the flux is introduced as
Typically the radiation field is axisymmetric around \(\vnhat\), and so \(\intsy\) depends on direction solely via the angle parameter \(\mu \equiv \vshat \cdot \vnhat\). Then, the flux simplifies to
(the lower bound on \(\mu\) is set to 0 rather than -1 under the assumption that there is no external radiation at the stellar surface).
Both \(\intsy\) and \(\flux\) can reasonably be dubbed a ‘spectrum’, as they each represent the distribution of electromagnetic energy with respect to wavelength. However, one should be careful to distinguish these elemental spectra from the observed spectrum of a star; this distinction is further clarified below.
Evaluating an elemental spectrum requires solution of the radiative transfer equation throughout the atmospheric layers composing the surface element. This is often far too computationally expensive to do on-the-fly. An alternative approach is to pre-calculate spectra across a multi-dimensional grid spanning a range of photospheric parameters, and then interpolate within this grid when an elemental spectrum is required for a specific \(\vx\). This is the fundamental purpose of MSG.
Evaluating Elemental Spectra
To evaluate an elemental specific intensity spectrum, MSG models the dependence of \(\intsy\) on each of its three arguments as follows:
Wavelength dependence is represented using piecewise-constant functions.
Angle dependence is represented using limb-darkening laws.
Photospheric parameter dependence is represented using tensor surface interpolation.
The following subsections discuss these in greater detail.
Wavelength Dependence
The \(\lambda\) dependence of the specific intensity is represented as a piecewise-constant function on a wavelength abscissa \(\lambda = \{\lambda_{1},\lambda_{2},\ldots,\lambda_{M}\}\):
(for brevity, the dependence of \(\intsy\) on \(\mu\) and \(\vx\) has been suppressed). Mapping intensity data onto a new abscissa \(\lambda' = \{\lambda'_{1},\lambda'_{2},\ldots\,\lambda'_{M'}\}\) is performed conservatively, according to the expression
Beyond its simplicity, the advantage of this approach (as compared to higher-order interpolations) is that the equivalent width of line profiles is preserved.
Angle Dependence
The \(\mu\) dependence of the specific intensity is represented using limb-darkening laws. Most familiar is the linear law
where \(\intsy(1)\) represents the normally emergent (\(\mu=1\)) intensity and \(c\) is the linear limb-darkening coefficient (as before, the dependence of the intensity on other parameters has been suppressed). An improved characterization involves additional \(\mu\)-dependent terms on the right-hand side; for instance, the four-coefficient law devised by Claret (2000) is
where there are now four limb-darkening coefficients \(c_{k}\).
The advantage of using limb-darkening laws is the ease with which other useful quantities can be calculated. For instance, the flux (3) can be evaluated analytically, as can any of the Eddington (1926) intensity moments (or E-moments, as MSG terms them):
MSG supports the following limb-darkening laws:
- CONST
Constant law, where \(I_{\lambda}\) has no dependence on \(\mu\) whatsoever. This is discussed further below.
- LINEAR
Linear law given in equation (4) above.
- SQRT
Square-root law introduced by Diaz-Cordoves & Gimenez (1992).
- QUAD
Quadratic law introduced by Wade & Rucinski (1985).
- CLARET
Four-coefficient law introduced by Claret (2000) and given in equation (5) above.
The choice of law is made during grid construction (see the Grid Tools appendix for more details). The coefficients appearing in the limb-darkening laws (e.g., \(c\) and \(c_{k}\)) are typically determined from least-squares fits to tabulations of the specific intensity. In cases where these tabulations include flux but not specific intensity data, the CONST law is used; the angle-independent specific intensity is determined so that it produces the correct flux when evaluated using equation (3).
Photospheric Parameter Dependence
The photospheric parameter dependence of the specific intensity is represented using cubic Hermite tensor product interpolation. The appendices provide a (relatively) gentle introduction to tensor product interpolation. The short version is that the intensity is modeled via piecewise-cubic functions of each component of \(\vx\), constructed to be continuous and smooth at the join between each piecewise region. The derivatives at these joins are estimated using second-order finite difference approximations involving neighboring points (or first-order at grid boundaries).
Grids often contain holes and/or ragged boundaries (the latter typically arising near the edge of the region of the \(\Teff-g\) plane corresponding to super-Eddington luminosity). When an interpolation tries to access such missing data, MSG either switches to a lower-order scheme, or (if there simply aren’t sufficient data to interpolate) signals an exception (see the Exception Handling chapter for further details).
Observed Spectra
Suppose we observe a star from Earth, at a distance \(d\) along unit direction vector \(\vdhat\). The energy measured by a detector of area \(\areao\), within the usual wavelength and time intervals, is
(here and subsequently the superscript \(^{\obs}\) should be read as ‘observed’), where the observed flux is introduced as
The integral here is similar to that in equation (2), but \(\vshat\) has been replaced by \(-\vdhat\), the solid angle element \(\diff{\Omega}\) has been replaced by \(\diff{\area}/d^{2}\), and the bounds of the integral are limited to the parts of the stellar surface that are visible from Earth.
However, for stars that are spherical and have photospheric parameters that don’t vary across their surface, further simplifications can be made. Let \(\theta\) and \(\phi\) be the colatitude and azimuth angles in a spherical coordinate system centered on the star and with polar axis antiparallel to \(\vdhat\). Then, the observed flux becomes
where \(R\) is the stellar radius. Assuming an axisymmetric radiation field, this further reduces to
With the substitution \(\mu = \cos\theta\) (replacing a spatial coordinate \(\theta\) with a directional one \(\mu\)) the result pops out that
Don’t be fooled by the apparent triviality of this result: it means that we need only the elemental flux spectrum, and not the specific intensity, to calculate the observed spectrum of a star. This is why many spectral grids in the literature include flux spectra instead of specific intensity spectra.
However, remember that equation (7) applies only for spherically symmetric and fixed-\(\vx\) stars. In more complex situations, for instance when a star is rotationally oblate, spotted, pulsating or even eclipsed, evaluation of \(\fluxo\) must proceed via the visible-disk integration appearing in equation (6), which requires the specific intensity.
Photometric Colors
To evaluate a photometric color, MSG convolves stellar spectra with an appropriate passband response function \(S'(\lambda)\). This function represents the combined sensitivity of the optical pathway, filter and detector. The passband-averaged specific intensity is defined as
meaning that \(S'(\lambda)\) is interpreted as an energy response function (see appendix A of Bessell & Murphy, 2012 for a discussion of the relationship between \(S'\) and the corresponding photon response function \(S\)). The passband-averaged observed flux follows from equation (6) as
and the apparent magnitude of the star is
where the normalizing flux \(\fluxz\) is set by the zero-point of the photometric system.
The convolution in (8) can be performed before or after the interpolations discussed above:
the ‘before’ option performs the convolution as a pre-processing step using the specgrid_to_photgrid tool to create a photgrid file from a specgrid file (as discussed in the Importing Data section). This is computationally more efficient, but requires a separate photgrid file to be created for each passband.
the ‘after’ option loads data from a specgrid file, but performs the convolution on-the-fly after each spectrum is interpolated. This is computationally less efficient, but incurs no storage requirements beyond the specgrid file.
Data Files
This chapter discusses the files in which MSG stores its data. These files adopt the HDF5 data format, a platform-neutral binary storage format with advanced features such as transparent decompression.
File Types
There are five types of MSG data files, distinguished by their differing content:
specint files store spectroscopic intensity data for a single combination of photospheric parameters.
photint files store photometric intensity data for a single combination of photospheric parameters.
specgrid files store spectroscopic intensity data over a grid of photospheric parameters.
photgrid files store photometric intensity data over a grid of photospheric parameters.
passband files store passband response functions, used to convert spectroscopic intensities into corresponding photometric intensities.
For a detailed description of each file type, refer to the Data Schema chapter.
Obtaining Data
MSG ships with a limited set of data files, sufficient to enable the walkthroughs. Additional files can be downloaded separately from the Grid Files and Passband Files appendices.
Importing Data
To import an existing spectroscopic grid into MSG, first convert the individual spectra into corresponding specint files. MSG provides a number of tools to assist with this conversion; see the Grid Tools appendix for further details.
The next step is to create a manifest (named, say,
manifest.txt
) listing all the specint files
composing the grid. This is a simple text file with each line naming
one file; for instance:
specint-0001.h5
specint-0002.h5
specint-0003.h5
Then, run the specint_to_specgrid tool to create a specgrid file:
$MSG_DIR/bin/specint_to_specgrid manifest.txt specgrid.h5
To build a photgrid file from the data in a specgrid file, run the specgrid_to_photgrid tool:
$MSG_DIR/bin/specgrid_to_photgrid specgrid.h5 passband.h5 photgrid.h5
…where passband.h5
is the name of the passband
file to use. Note that it’s not always necessary to create a
photgrid file, as MSG can convolve with passbands on the
fly (as discussed in the Photometric Colors section).
Data Caching
This chapter discusses MSG’s data caching, a key component of its performance and scaleability. Each instance of a grid data type (e.g., specgrid_t
in Fortran, or pymsg.PhotGrid
in Python) has an associated cache that temporarily stores spectroscopic and/or photometric data for use in interpolations. This grid cache is initially empty when the instance is created, but grows in size as data are loaded on-demand from disk. Eventually, once a user-definable size limit is reached, old data are discarded from the cache to make room for new.
Caching Demo
With its caching functionality, MSG can in principle work with grids that are much larger on disk than available computer memory. The following Python code fragments showcase this capability, while also illustrating the basics of inspecting and controlling a grid cache.
To get started, download the sg-CAP18-coarse.h5 spectroscopic grid from the CAP18 Grids appendix, storing it your working directory. Then, initialize MSG and load the grid:
[1]:
# Import standard modules
import sys
import os
import numpy as np
import matplotlib.pyplot as plt
import random as rn
import time
# Import pymsg
import pymsg
# Set plot parameters
%matplotlib inline
plt.rcParams.update({'font.size': 16})
# Load the SpecGrid
specgrid = pymsg.SpecGrid('sg-CAP18-coarse.h5')
Next, inspect the state of the cache attached to the specgrid
object:
[11]:
# Inspect the cache
print(f'cache usage: {specgrid.cache_usage} MB')
print(f'cache limit: {specgrid.cache_limit} MB')
cache usage: 128 MB
cache limit: 128 MB
The pymsg.SpecGrid.cache_usage
property reports the current memory usage (in megabytes) of the cache. Because we’ve just loaded specgrid
, this usage is zero — no data have (yet) been read into memory. The pymsg.SpecGrid.cache_limit
property likewise specifies the maximum memory usage (again, in megabytes) of the cache. By default, this limit is set to 128 MB.
Let’s define a function that interpolates the visible-spectrum flux for a randomly-chosen set of photospheric parameters. Since we’re interested primarily in the behavior of the grid cache, the function returns its execution time rather than interpolation result.
[3]:
# Random flux routine
def random_flux():
start_time = time.process_time()
# Set up the wavelength abscissa (visible-spectrum)
lam_min = 3000.
lam_max = 7000.
lam = np.linspace(lam_min, lam_max, 1000)
# Loop until a valid set of photospheric parameters is found
while True:
# Choose random photospheric parameters
x = {}
for label in specgrid.axis_labels:
x[label] = rn.uniform(specgrid.axis_x_min[label], specgrid.axis_x_max[label])
# Interpolate the flux, allowing for the fact that the
# photospheric parameters may fall in a grid void
try:
F = specgrid.flux(x, lam)
break
except LookupError:
pass
end_time = time.process_time()
return end_time-start_time
Running this function a few times (for now ignoring the return value), we can watch the cache fill up to the 128-MB limit.
[4]:
# Run random_flux three times
for i in range(3):
random_flux()
print(f'cache usage: {specgrid.cache_usage} MB')
cache usage: 86 MB
cache usage: 128 MB
cache usage: 128 MB
If needed this limit can be increased by changing the pymsg.SpecGrid.cache_limit
property, as this example demonstrates:
[5]:
# Increase the cache limit to 256 MB
specgrid.cache_limit = 256
# Run random_flux three times
for i in range(3):
random_flux()
print(f'cache usage: {specgrid.cache_usage} MB')
cache usage: 215 MB
cache usage: 256 MB
cache usage: 256 MB
Finally, the contents of the cache can be flushed (and the memory freed up) with a call to the pymsg.SpecGrid.flush_cache()
method:
[6]:
# Flush the cache
print(f'cache usage: {specgrid.cache_usage} MB')
specgrid.flush_cache()
print(f'cache usage: {specgrid.cache_usage} MB')
cache usage: 256 MB
cache usage: 0 MB
Wavelength Subsetting
When working with spectroscopic grids, it’s often the case one is interested only in a subset of the wavelength range covered by the grid. As an example, the random_flux()
function defined above focuses on just the visible part of the electromagnetic spectrum. In such cases, you can direct MSG to cache only the data within the subset by setting the pymsg.SpecGrid.cache_lam_min
and pymsg.SpecGrid.cache_lam_min
properties. This has the benefit of slowing the rate at which the cache fills up, as the following example demonstrates:
[7]:
# Subset to the visible part of the spectrum
specgrid.cache_lam_min = 3000.
specgrid.cache_lam_max = 7000.
# Run random_flux three times
for i in range(3):
random_flux()
print(f'cache usage: {specgrid.cache_usage} MB')
cache usage: 18 MB
cache usage: 30 MB
cache usage: 55 MB
(compare this against the cache growth in the preceding examples).
Performance Impact
Let’s now explore how caching can have a significant impact on MSG’s performance. Run the following code, which expands the cache limit further to 512 MB, and then calls random_spectrum()
300 times, storing the usage and execution time after each call:
[8]:
# Flush the cache and set the limit to 512 MB
specgrid.flush_cache()
specgrid.cache_limit = 512
# Allocate usage & timing arrays
n = 300
cache_usages_512 = np.empty(n, dtype=int)
exec_timings_512 = np.empty(n)
# Call random_flux
for i in range(n):
exec_timings_512[i] = random_flux()
cache_usages_512[i] = specgrid.cache_usage
Then, plot the results:
[9]:
# Plot cache usage and execution time
fig, ax = plt.subplots(ncols=2, figsize=[12,6])
ax[0].plot(cache_usages_512)
ax[0].set_xlabel('call #')
ax[0].set_ylabel('cache usage (MB)')
ax[1].plot(exec_timings_512)
ax[1].set_xlabel('call #')
ax[1].set_ylabel('exution time (s)')
[9]:
Text(0, 0.5, 'exution time (s)')

In the left-hand panel we see the cache usage grow and eventually asymptote toward the point where the entire grid (around 400 MB) is resident in memory. In the right-hand panel the execution time rapidly drops from ~0.2s down to ~0.01s, as more and more data can be read from memory (fast) rather than disk (slow).
Let’s now repeat the exercise, but for a cache limit reduced back down to 128 MB:
[10]:
# Flush the cache and set the limit to 128 MB
specgrid.flush_cache()
specgrid.cache_limit = 128
# Allocate usage & timing arrays
n = 300
cache_usages_128 = np.empty(n, dtype=int)
exec_timings_128 = np.empty(n)
# Call random_flux
for i in range(n):
exec_timings_128[i] = random_flux()
cache_usages_128[i] = specgrid.cache_usage
# Plot cache usage and execution time
fig, ax = plt.subplots(ncols=2, figsize=[12,6])
ax[0].plot(cache_usages_512, label='512 MB limit')
ax[0].plot(cache_usages_128, label='128 MB limit')
ax[0].set_xlabel('call #')
ax[0].set_ylabel('cache usage (MB)')
ax[0].legend(loc=4)
ax[1].plot(exec_timings_512)
ax[1].plot(exec_timings_128)
ax[1].set_xlabel('call #')
ax[1].set_ylabel('exution time (s)')
[10]:
Text(0, 0.5, 'exution time (s)')

The left-hand panel reveals that, as expected, the growth of the cache usage stops once the limit is reached. Because not all the data can be held in memory, the performance gains from the cache are more modest than in the previous case, with the execution times in the right-hand panel averaging 0.15s.
Choosing Cache Settings
As the demonstration above makes clear, the choice of the pymsg.SpecGrid.cache_usage
property can have a significant impact on the performance of MSG. Ideally, this property should be set to exceed the total size of the grid 1; but if that’s not possible due to limited computer memory, try some of the following measures:
set the
pymsg.SpecGrid.cache_usage
property to at least \(4^{N}\) times the size of a single node, where \(N\) is the number of dimensions in a grid (see the Tensor Product Interpolation appendix).if undertaking a sequence of interpolations, reorder them so that ones with similar photospheric parameters are grouped together.
for spectroscopic grids, use the
pymsg.SpecGrid.cache_lam_min
and/orpymsg.SpecGrid.cache_lam_max
properties, as discussed above in the Wavelength Subsetting section.
Technical Details
The basic data unit of MSG’s caches is a grid node, representing the spectroscopic or photometric intensity for a single set of photospheric parameters. During an interpolation the nodes required for the calculation are accessed through the cache, which reads them from disk as necessary. After every node access, the current cache usage is compared against the user-specified limit. If it exceeds this limit, then one or more nodes are evicted from the cache using a least-recently used (LRU) algorithm.
This implementation — with eviction taking place after cache access — means that in principle the cache limit can be zero. For performance reasons, however, this is not recommended.
Footnotes
- 1
Note that the in-memory size of a grid is typically larger than its on-disk size, due to compression.
Exception Handling
This chapter discusses how MSG handles exceptions. These can arise in a variety of circumstances:
Attempts to load grids from files that are missing or corrupt.
Attempts to interpolate at locations outside the bounds of a grid.
Attempts to interpolate at locations where grid data are missing (so-called grid voids).
Attempts to interpolate with incomplete specification of photospheric parameters.
Attempts to interpolate for invalid wavelength and/or angle parameters.
When an exception occurs, how it’s signaled depends on the language interface being used.
Python
Using the Python interface, exceptions are signaled using the language’s built-in exception handling capabilities. The list of exceptions that can be thrown by each function is provided in the Python Interface chapter.
Fortran
Using the Fortran interface, exceptions are signaled through the
optional procedure argument stat
. If this argument is
present, then on return it is set to one of the status code values
defined in the Fortran parameters section. The
value STAT_OK
indicates that no problem was encountered;
other values signal an exception. If stat
is not present when
an exception occurs, then execution halts with an error message
printed to standard output.
C
Using the C interface, exceptions are signaled through the pointer
function argument stat
. If this pointer is non-null, then on
return the pointer target *stat
is set to one of the status code values
defined in the C enums section. The value
STAT_OK
indicates that no problem was encountered; other
values signal an exception. If stat
is null when an exception
occurs, then execution halts with an error message printed to standard
output.
Troubleshooting
Missing forum.inc
During compilation, if the error include file 'forum.inc' not
found
arises then it’s likely you have an incomplete copy of the
source code. Verify that when you checked out
the code from GitHub, you used the --recurse-submodules
option.
Other Problems
If you encounter something else that doesn’t work as it should, please open an issue on GitHub. In your issue, please specify the following:
The release of MSG you are using (e.g. 1.1).
The operating system and architecture you are using (e.g., Mac OS 13.1 on Intel).
A brief description of the problem.
A narrative of the steps to reproduce the problem.
Contributing
Code
To contribute code to the MSG project, from bug fixes through to major new features, the preferred method is the standard fork—branch—pull-request paradigm described here. If your pull request languishes for more than a couple of weeks without receiving a response, consider opening an issue to give the developers a nudge.
Documentation
To contribute documentation toward the project, follow the same
approach as above for code. (The ReStructured Text source files for
the documentation reside in the docs/source
subdirectory.)
Data
To contribute data toward the project (in the form of spectral intensity or flux grids), a number of options are available. The Importing Data section explains how to import an existing grid; but if this lies beyond your technical expertise, then the MSG team will be happy to assist. Either way, if you want to make your data publicly available and listed in the Grid Files chapter, then please open an issue.
Installation
This chapter discusses MSG installation in detail. If you just want to get up and running, have a look at the Quick Start chapter.
Pre-Requisites
To compile and use MSG, you’ll need the following software components:
A modern (2008+) Fortran compiler
The LAPACK linear algebra library
The LAPACK95 Fortran 95 interface to LAPACK
The fypp Fortran pre-processor
The HDF5 data management library
On Linux and MacOS platforms, these components are bundled together in the MESA Software Development Kit (SDK). Using this SDK is strongly recommended.
If you’re planning on using the pymsg
Python module, then
you’ll also need the following components:
Python 3.7 (or more recent)
NumPy 1.15 (or more recent)
Building MSG
Download the MSG source code, and unpack it from the command line using the tar utility:
tar xf msg-1.1.2.tar.gz
Set the MSG_DIR
environment variable with the path to the
newly created source directory; this can be achieved, e.g., using the
realpath command1:
export MSG_DIR=$(realpath msg-1.1.2)
Finally, compile MSG using the make utility:
make -j -C $MSG_DIR
(the -j flags tells make to use multiple cores, speeding up the build). If things go awry, consult the Troubleshooting chapter.
Testing MSG
To test MSG, use the command
make -C $MSG_DIR test
This runs unit tests for the various Fortran modules that together compose the MSG library. At the end of the test sequence, a summary of the number of tests passed and failed is printed. All tests should pass; if one or more fails, then please open an issue to report the problem.
Installing the pymsg
Module
To install the pymsg
Python module, use the pip tool:
pip install $MSG_DIR/python
You can alternatively add the $MSG_DIR/python/src
directory to
the PYTHONPATH
environment variable. Note that in order for
pymsg
to function correctly, the MSG_DIR
environment variable must be set at Python runtime (this variable
allows the module to find the Python extension that interfaces to the
back-end).
Custom Builds
Custom builds of MSG can be created by setting certain environment
variables, and/or variables in the file
$MSG_DIR/src/build/Makefile
, to the value yes
. The
following variables are currently supported:
- DEBUG
Enable debugging mode (default
no
)- FPE
Enable floating point exception checks (default
yes
)- OMP
Enable OpenMP parallelization (default
yes
)- PYTHON
Enable building of the Python extension (default
yes
)- TEST
Enable building of testing tools (default
yes
)- TOOLS
Enable building of development tools (default
yes
)
If a variable is not set, then its default value is assumed.
GitHub Access
Sometimes, you’ll want to try out new features in MSG that haven’t yet made it into a formal release. In such cases, you can check out MSG directly from the rhdtownsend/msg git repository on GitHub:
git clone --recurse-submodules https://github.com/rhdtownsend/msg.git
However, a word of caution: MSG is under constant development, and
features in the main
branch can change without warning.
footnote
Python Interface
The Python interface is provided through the pymsg
module.
API Specification
Classes
- class pymsg.SpecGrid(file_name)
The SpecGrid class represents a grid of spectroscopic intensity data.
This grid may be used to interpolate the intensity (or related quantities) across a wavelength abscissa and for a set of photospheric parameter values.
- __init__(file_name)
SpecGrid constructor (via loading data from a specgrid file).
- Parameters
file_name (string) – Name of the file
- Returns
Constructed object.
- Return type
- Raises
FileNotFound – If the the file cannot be found.
TypeError – If the file contains an incorrect datatype.
- property rank
Number of dimensions in grid.
- Type
int
- property axis_labels
Photospheric parameter axis labels.
- Type
list
- property axis_x_min
Photospheric parameter axis minima.
- Type
dict
- property axis_x_max
Photospheric parameter axis maxima.
- Type
dict
- property lam_min
Minimum wavelength of grid (Å).
- Type
double
- property lam_max
Maximum wavelength of grid (Å).
- Type
double
- property cache_lam_min
Minimum wavelength of grid cache (Å).
- Type
double
- property cache_lam_max
Maximum wavelength of grid cache (Å).
- Type
double
- property cache_usage
Current memory usage of grid cache (MB).
- Type
int
- property cache_limit
Maximum memory usage of grid cache (MB).
- Type
int
- flush_cache()
Flush the grid cache
- intensity(x, mu, lam, deriv=None)
Evaluate the spectroscopic intensity.
- Parameters
x (dict) – Photospheric parameters; keys must match axis_labels property, values must be double.b
mu (double) – Cosine of angle of emergence relative to surface normal.
lam (numpy.ndarray) – Wavelength abscissa (Å).
deriv (dict, optional) – Flags indicating whether to evaluate derivative with respect to each photospheric parameter; keys must match the axis_labels property, values must be boolean.
- Returns
Spectroscopic intensity (erg/cm^2/s/Å/sr) in bins delineated by lam; length len(lam)-1.
- Return type
- Raises
KeyError – If x does not define all keys appearing in the axis_labels property.
ValueError – If x, mu, or any part of the wavelength abscissa falls outside the bounds of the grid.
LookupError – If x falls in a grid void.
- E_moment(x, k, lam, deriv=None)
Evaluate the spectroscopic intensity E-moment.
- Parameters
x (dict) – Photospheric parameters; keys must match axis_labels property, values must be double.
k (int) – Degree of moment.
lam (numpy.ndarray) – Wavelength abscissa (Å).
deriv (dict, optional) – Flags indicating whether to evaluate derivative with respect to each photospheric parameter; keys must match the axis_labels property, values must be boolean.
- Returns
Spectroscopic intensity E-moment (erg/cm^2/s/Å) in bins delineated by lam; length len(lam)-1.
- Return type
- Raises
KeyError – If x does not define all keys appearing in the axis_labels property.
ValueError – If x, k, or any part of the wavelength abscissa falls outside the bounds of the grid.
LookupError – If x falls in a grid void.
- D_moment(x, l, lam, deriv=None)
Evaluate the spectroscopic intensity D-moment.
- Parameters
x (dict) – Photospheric parameters; keys must match axis_labels property, values must be double.
l (int) – Harmonic degree of moment.
lam (numpy.ndarray) – Wavelength abscissa (Å).
deriv (dict, optional) – Flags indicating whether to evaluate derivative with respect to each photospheric parameter; keys must match the axis_labels property, values must be boolean.
- Returns
Spectroscopic intensity D-moment (erg/cm^2/s/Å) in bins delineated by lam; length len(lam)-1.
- Return type
- Raises
KeyError – If x does not define all keys appearing in the axis_labels property.
ValueError – If x, l, or any part of the wavelength abscissa falls outside the bounds of the grid.
LookupError – If x falls in a grid void.
- flux(x, lam, deriv=None)
Evaluate the spectroscopic flux.
- Parameters
x (dict) – Photospheric parameters; keys must match axis_labels property, values must be double.
lam (numpy.ndarray) – Wavelength abscissa (Å)
deriv (dict, optional) – Flags indicating whether to evaluate derivative with respect to each photospheric parameter; keys must match the axis_labels property, values must be boolean.
- Returns
Spectroscopic flux (erg/cm^2/s/Å) in bins delineated by lam; length len(lam)-1.
- Return type
- Raises
KeyError – If x does not define all keys appearing in the axis_labels property.
ValueError – If x or any part of the wavelength abscissa falls outside the bounds of the grid.
LookupError – If x falls in a grid void.
- class pymsg.PhotGrid(file_name, passband_file_name=None)
The PhotGrid class represents a grid of photometric intensity data.
This grid may be used to interpolate the intensity (or related quantities) for a set of photospheric parameter values.
- __init__(file_name, passband_file_name=None)
- PhotGrid constructor (via loading data from a photgrid file, or
from a specgrid file together with a passband file).
- Parameters
file_name (string) – Name of grid file.
passband_file_name (string) – Name of passband file (if file_name corresponds to a specgrid file)
- Returns
Constructed object.
- Return type
- Raises
FileNotFound – If either file cannot be found.
TypeError – If either file contains an incorrect datatype.
- property rank
Number of dimensions in grid.
- Type
int
- property axis_labels
Photospheric parameter axis labels.
- Type
list
- property axis_x_min
Photospheric parameter axis minima.
- Type
dict
- property axis_x_max
Photospheric parameter axis maxima.
- Type
dict
- property cache_usage
Current memory usage of grid cache (MB).
- Type
int
- property cache_limit
Maximum memory usage of grid cache (MB).
- Type
int
- flush_cache()
Flush the grid cache
- intensity(x, mu, deriv=None)
Evaluate the photometric intensity, normalized to the zero- point flux.
- Parameters
x (dict) – Photospheric parameters; keys must match axis_labels property, values must be double.
mu (double) – Cosine of angle of emergence relative to surface normal.
deriv (dict, optional) – Flags indicating whether to evaluate derivative with respect to each photospheric parameter; keys must match the axis_labels property, values must be boolean.
- Returns
photometric intensity (/sr).
- Return type
double
- Raises
KeyError – If x does not define all keys appearing in the axis_labels property.
ValueError – If x or mu falls outside the bounds of the grid.
LookupError – If x falls in a grid void.
- E_moment(x, k, deriv=None)
Evaluate the photometric intensity E-moment, normalized to the zero-point flux.
- Parameters
x (dict) – Photospheric parameters; keys must match axis_labels property, values must be double.
k (int) – Degree of moment.
deriv (dict, optional) – Flags indicating whether to evaluate derivative with respect to each photospheric parameter; keys must match the axis_labels property, values must be boolean.
- Returns
photometric intensity E-moment.
- Return type
double
- Raises
KeyError – If x does not define all keys appearing in the axis_labels property.
ValueError – If x or k falls outside the bounds of the grid.
LookupError – If x falls in a grid void.
- D_moment(x, l, deriv=None)
Evaluate the photometric intensity D-moment, normalized to the zero-point flux.
- Parameters
x (dict) – Photospheric parameters; keys must match axis_labels property, values must be double.
l (int) – Harmonic degree of moment.
deriv (dict, optional) – Flags indicating whether to evaluate derivative with respect to each photospheric parameter; keys must match the axis_labels property, values must be boolean.
- Returns
photometric intensity D-moment.
- Return type
double
- Raises
KeyError – If x does not define all keys appearing in the axis_labels property.
ValueError – If x or l falls outside the bounds of the grid.
LookupError – If x falls in a grid void.
- flux(x, deriv=None)
Evaluate the photometric flux, normalized to the zero-point flux.
- Parameters
x (dict) – Photospheric parameters; keys must match axis_labels property, values must be double.
deriv (dict, optional) – Flags indicating whether to evaluate derivative with respect to each photospheric parameter; keys must match the axis_labels property, values must be boolean.
- Returns
photometric flux.
- Return type
double
- Raises
KeyError – If x does not define all keys appearing in the axis_labels property.
ValueError – If x or l falls outside the bounds of the grid.
LookupError – If x falls in a grid void.
Fortran Interface
The Fortran interface is provided through the fmsg_m
module,
which defines derived types together with supporting parameters and
procedures.
API Specification
Derived Types
- type specgrid_t
The specgrid_t type represents a grid of spectroscopic data.
This grid may be used to interpolate the intensity (or related quantities) across a wavelength abscissa and for a set of photospheric parameter values.
- subroutine get_rank(rank)
Get the rank (dimension) of the grid.
- Parameters
rank [integer ,out] :: Rank.
- subroutine get_axis(i, axis)
Get an axis of the grid.
- Parameters
i [integer ,in] :: Index of axis (between 1 and rank).
axis [axis_t ,out] :: Axis.
- subroutine get_lam_min(lam_min)
Get the minimum wavelength of the grid.
- Parameters
lam_min [real(RD),out] :: Minimum wavelength (Å).
- subroutine get_lam_max(lam_max)
Get the maximum wavelength of the grid.
- Parameters
lam_max [real(RD),out] :: Maximum wavelength (Å).
- subroutine get_cache_lam_min(cache_lam_min)
Get the minimum wavelength of the grid cache.
- Parameters
lam_min [real(RD),out] :: Minimum wavelength (Å).
- subroutine get_cache_lam_max(cache_lam_max)
Get the maximum wavelength of the grid cache.
- Parameters
cache_lam_max [real(RD),out] :: Maximum wavelength (Å).
- subroutine get_cache_limit(cache_limit)
Get the maximum memory usage of the grid cache.
- Parameters
cache_limit [integer ,out] :: Maximum memory usage (MB).
- subroutine get_cache_usage(cache_usage)
Get the current memory usage of the grid cache.
- Parameters
cache_usage [integer ,out] :: Current memory usage (MB)
- subroutine set_cache_lam_min(cache_lam_min, stat)
Set the minimum wavelength of the grid cache.
- Parameters
lam_min [real(RD),in] :: Minimum wavelength (Å).
- Options
stat [integer ,out] :: Status code.
- subroutine set_cache_lam_max(cache_lam_max, stat)
Set the maximum wavelength of the grid cache.
- Parameters
cache_lam_max [real(RD),in] :: Maximum wavelength (Å).
- Options
stat [integer ,out] :: Status code.
- subroutine set_cache_limit(cache_limit, stat)
Set the maximum memory usage of the grid cache.
- Parameters
cache_limit [integer ,in] :: Maximum memory usage (MB).
- Options
stat [integer ,out] :: Status code.
- subroutine interp_intensity(x_vec, mu, lam, I, stat, deriv_vec)
Interpolate the spectroscopic intensity.
- Parameters
x_vec (*) [real(RD),in] :: Photospheric parameter values.
mu [real(RD),in] :: Cosine of angle of emergence relative to surface normal.
lam (*) [real(RD),in] :: Wavelength abscissa (Å).
I (*) [real(RD),out] :: Spectroscopic intensity (erg/cm^2/s/Å/sr) in bins delineated by lam; length LEN(lam)-1.
- Options
stat [integer ,out] :: Status code.
deriv_vec (*) [logical ,in] :: Derivative flags.
- subroutine interp_E_moment(x_vec, k, lam, E, stat, deriv_vec)
Interpolate the spectroscopic intensity E-moment.
- Parameters
x_vec (*) [real(RD),in] :: Photospheric parameter values.
k [integer ,in] :: Degree of moment.
lam (*) [real(RD),in] :: Wavelength abscissa (Å).
E (*) [real(RD),out] :: Spectroscopic intensity E-moment (erg/cm^2/s/Å) in bins delineated by lam; length LEN(lam)-1.
- Options
stat [integer ,out] :: Status code.
deriv_vec (*) [logical ,in] :: Derivative flags.
- subroutine interp_D_moment(x_vec, l, lam, D, stat, deriv_vec)
Interpolate the spectroscopic intensity D-moment.
- Parameters
x_vec (*) [real(RD),in] :: Photospheric parameter values.
l [integer ,in] :: Harmonic degree of moment.
lam (*) [real(RD),in] :: Wavelength abscissa (Å).
D (*) [real(RD),out] :: Spectroscopic intensity D-moment (erg/cm^2/s/Å) in bins delineated by lam; length LEN(lam)-1.
- Options
stat [integer ,out] :: Status code.
deriv_vec (*) [logical ,in] :: Derivative flags.
- subroutine interp_flux(x_vec, lam, I, stat, deriv_vec)
Interpolate the spectroscopic flux.
- Parameters
x_vec (*) [real(RD),in] :: Photospheric parameter values.
lam (*) [real(RD),in] :: Wavelength abscissa (Å).
F (*) [real(RD),out] :: Spectroscopic flux (erg/cm^2/s/Å) in bins delineated by lam; length LEN(lam)-1.
- Options
stat [integer ,out] :: Status code.
deriv_vec (*) [logical ,in] :: Derivative flags.
- type photgrid_t
The photgrid_t type represents a grid of photometric data.
This grid may be used to interpolate the intensity (or related quantities) for a set of photospheric parameter values.
- subroutine get_rank(rank)
Get the rank (dimension) of the grid.
- Parameters
rank [integer ,out] :: Returned rank.
- subroutine get_axis(i, axis)
Get an axis of the grid.
- Parameters
i [integer ,in] :: Index of axis (between 1 and rank)
axis [axis_t ,out] :: Returned axis.
- subroutine get_cache_limit(cache_limit)
Get the maximum memory usage of the cache.
- Parameters
cache_limit [integer ,out] :: Maximum memory usage (MB).
- subroutine get_cache_usage(cache_usage)
Get the current memory usage of the cache.
- Parameters
cache_usage [integer ,out] :: Current memory usage (MB)
- subroutine set_cache_limit(cache_limit, stat)
Set the maximum memory usage of the cache.
- Parameters
cache_limit [integer ,in] :: Maximum memory usage (MB).
- Options
stat [integer ,out] :: Status code.
- subroutine interp_intensity(x_vec, mu, I, stat, deriv_vec)
Interpolate the photometric intensity, normalized to the zero-point flux.
- Parameters
x_vec (*) [real(RD),in] :: Photospheric parameter values.
mu [real(RD),in] :: Cosine of angle of emergence relative to surface normal.
I [real(RD),out] :: Photometric intensity (/sr).
- Options
stat [integer ,out] :: Status code.
deriv_vec (*) [logical ,in] :: Derivative flags.
- subroutine interp_E_moment(x_vec, k, E, stat, deriv_vec)
Interpolate the photometric E-moment, normalized to the zero-point flux.
- Parameters
x_vec (*) [real(RD),in] :: Photospheric parameter values.
k [integer ,in] :: Degree of of moment.
E [real(RD),out] :: Photometric E-moment.
- Options
stat [integer ,out] :: Status code.
deriv_vec (*) [logical ,in] :: Derivative flags.
- subroutine interp_D_moment(x_vec, l, D, stat, deriv_vec)
Interpolate the photometric D-moment, normalized to the zero-point flux.
- Parameters
x_vec (*) [real(RD),in] :: Photospheric parameter values.
l [integer ,in] :: Harmonic degree of moment.
D [real(RD),out] :: Photometric D-moment.
- Options
stat [integer ,out] :: Status code.
deriv_vec (*) [logical ,in] :: Derivative flags.
- subroutine interp_flux(x_vec, F, stat, deriv_vec)
Interpolate the photometric flux, normalized to the zero-point flux.
- Parameters
x_vec (*) [real(RD),in] :: Photospheric parameter values.
F [real(RD),out] :: Photometric flux.
- Options
stat [integer ,out] :: Status code.
deriv_vec (*) [logical ,in] :: Derivative flags.
- type axis_t
The axis_t type represents a grid axis.
- subroutine get_n(n)
Get the number of points making up the axis.
- Parameters
rank [integer ,out] :: Number of points.
- subroutine get_x_min(x_min)
Get the minimum value of the axis.
- Parameters
x_min [real(RD),out] :: Minimum value.
- subroutine get_x_max(x_max)
Get the maximum value of the axis.
- Parameters
x_max [real(RD),out] :: Maximum value.
- subroutine get_label(label)
Get the axis label.
- Parameters
character (*) :: Label.
- subroutine fetch(i, x, stat)
Fetch an axis value.
- Parameters
i [integer ,in] :: Index of value (from 1 to n).
x [real(RD),out] :: Axis value.
- Options
stat [integer ,out] :: Status code.
- subroutine locate(x, i)
Locate where along the axis a value falls.
- Parameters
x [real(RD),out] :: Value to locate.
i [integer ,out] :: Location index.
Parameters
- int RD
Kind type parameter for reals.
- int LABEL_LEN
Length of
axis_t
labels.
- int ~fmsg_m/STAT_OK
Status code indicating procedure call completed without error.
- int STAT_OUT_OF_BOUNDS_RANGE_LO
Status code indicating procedure call encountered an out-of-bounds reference, below the range minimum.
- int STAT_OUT_OF_BOUNDS_RANGE_HI
Status code indicating procedure call encountered an out-of-bounds reference, above the range maximum.
- int STAT_OUT_OF_BOUNDS_AXIS_LO
Status code indicating procedure call encountered an out-of-bounds reference, below the axis minimum.
- int STAT_OUT_OF_BOUNDS_AXIS_HI
Status code indicating procedure call encountered an out-of-bounds reference, above the axis maximum.
- int STAT_OUT_OF_BOUNDS_LAM_LO
Status code indicating procedure call encountered an out-of-bounds reference, below the wavelength minimum.
- int STAT_OUT_OF_BOUNDS_LAM_HI
Status code indicating procedure call encountered an out-of-bounds reference, above the wavelength maximum.
- int STAT_OUT_OF_BOUNDS_MU_LO
Status code indicating procedure call encountered an out-of-bounds reference, below the emergence cosine minimum.
- int STAT_OUT_OF_BOUNDS_MU_HI
Status code indicating procedure call encountered an out-of-bounds reference, above the emergence cosine maximum.
- int STAT_UNAVAILABLE_DATA
Status code indicating procedure call encountered unavailable data.
- int STAT_INVALID_ARGUMENT
Status code indicating procedure call encountered an invalid argument.
- int STAT_FILE_NOT_FOUND
Status code indicating procedure call encountered a file that could not be found.
- int STAT_INVALID_FILE_TYPE
Status code indicating procedure call encountered a file with an invalid type.
- int STAT_INVALID_GROUP_TYPE
Status code indicating procedure call encountered a file group with an invalid type.
- int STAT_INVALID_GROUP_REVISION
Status code indicating procedure call encountered a file group with an invalid revision number.
Procedures
- subroutine load_specgrid(specgrid_file_name, specgrid, stat)
Create a new
specgrid_t
by loading data from a specgrid file.- Parameters
specgrid_file_name [character(*),in] :: Name of the specgrid file.
specgrid [specgrid_t ,out] :: Grid object.
- Options
stat [integer ,out] :: Status code.
- subroutine load_photgrid(photgrid_file_name, photgrid, stat)
Create a new
photgrid_t
by loading data from a photgrid file.- Parameters
photgrid_file_name [character(*),in] :: Name of the photgrid file.
photgrid [specgrid_t ,out] :: Grid object.
- Options
stat [integer ,out] :: Status code.
- subroutine load_photgrid_from_specgrid(specgrid_file_name, passband_file_name, photgrid, stat)
Create a new
photgrid_t
by loading data from a specgrid file, and convolving on-the-fly with a response function loaded from a passband file.
- subroutine get_version(version)
Get the version of the MSG library
- Parameters
character (*) :: Version string.
Compiling/Linking
The module file fmsg_m
for the Fortran interface is provided
in the directory $MSG_DIR/include
, and executables should be
linked against $MSG_DIR/lib/libfmsg.so
(Linux) or
$MSG_DIR/lib/libfmsg.dylib
(MacOS). To simplify this process,
a script $MSG_DIR/scripts/fmsg_link
is provided that writes
the appropriate linker commands to standard output. This script can be
used to compile/link a program with gfortran as follows:
gfortran -o myprogram myprogram.f90 -I $MSG_DIR/include `$MSG_DIR/scripts/fmsg_link`
C Interface
The C interface is provided through the cmsg.h
header file,
which defines typedefs, enums and functions.
API Specification
Typedefs
-
type SpecGrid
A void * pointer to a Fortran
specgrid_t
spectroscopic grid object.
-
type PhotGrid
A void * pointer to a Fortran
specgrid_t
photometric grid object.
Enums
-
enum Stat
-
enumerator STAT_OK
Status code indicating function call completed without error.
-
enumerator STAT_OUT_OF_BOUNDS_RANGE_LO
Status code indicating function call encountered an out-of-bounds reference, below the range minimum.
-
enumerator STAT_OUT_OF_BOUNDS_RANGE_HI
Status code indicating function call encountered an out-of-bounds reference, above the range maximum.
-
enumerator STAT_OUT_OF_BOUNDS_AXIS_LO
Status code indicating function call encountered an out-of-bounds reference, below the axis minimum.
-
enumerator STAT_OUT_OF_BOUNDS_AXIS_HI
Status code indicating function call encountered an out-of-bounds reference, above the axis maximum.
-
enumerator STAT_OUT_OF_BOUNDS_LAM_LO
Status code indicating function call encountered an out-of-bounds reference, below the wavelength minimum.
-
enumerator STAT_OUT_OF_BOUNDS_LAM_HI
Status code indicating function call encountered an out-of-bounds reference, above the wavelength maximum.
-
enumerator STAT_OUT_OF_BOUNDS_MU_LO
Status code indicating function call encountered an out-of-bounds reference, below the emergence cosine minimum.
-
enumerator STAT_OUT_OF_BOUNDS_MU_HI
Status code indicating function call encountered an out-of-bounds reference, above the emergence cosine maximum.
-
enumerator STAT_UNAVAILABLE_DATA
Status code indicating function call encountered unavailable data.
-
enumerator STAT_INVALID_ARGUMENT
Status code indicating function call encountered an invalid argument.
-
enumerator STAT_FILE_NOT_FOUND
Status code indicating function call encountered a file that could not be found.
-
enumerator STAT_INVALID_FILE_TYPE
Status code indicating function call encountered a file with an invalid type.
-
enumerator STAT_INVALID_GROUP_TYPE
Status code indicating function call encountered a file group with an invalid type.
-
enumerator STAT_INVALID_GROUP_REVISION
Status code indicating function call encountered a file group with an invalid revision number.
-
enumerator STAT_OK
Functions
SpecGrid
-
void load_specgrid(const char *specgrid_file_name, SpecGrid *specgrid, Stat *stat)
Create a new
SpecGrid
by 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
specgrid
grid, 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_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, double x_vec[], double mu, int n, double lam[], double I[], Stat *stat, bool deriv_vec[])
Interpolate the spectroscopic intensity.
- Parameters
specgrid – Grid object.
x_vec – Photospheric parameter values.
mu – Cosine of angle of emergence relative to surface normal.
n – Number of points in wavelength abscissa.
lam[n] – Wavelength abscissa (Å).
I[n-1] – Spectroscopic intensity (erg/cm^2/s/Å/sr) in bins delineated by lam
stat – Status code (set to NULL if not required).
deriv_vec – Derivative flags (set to NULL if not required).
-
void interp_specgrid_E_moment(SpecGrid specgrid, double x_vec[], int k, int n, double lam[], double E[], Stat *stat, bool deriv_vec[])
Interpolate the spectroscopic intensity E-moment.
- Parameters
specgrid – Grid object.
x_vec – Photospheric parameter values.
k – Degree of moment.
n – Number of points in wavelength abscissa.
lam[n] – Wavelength abscissa (Å).
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 – Derivative flags (set to NULL if not required).
-
void interp_specgrid_D_moment(SpecGrid specgrid, double x_vec[], int l, int n, double lam[], double D[], Stat *stat, bool deriv_vec[])
Interpolate the spectroscopic intensity D-moment.
- Parameters
specgrid – Grid object.
x_vec – Photospheric parameter values.
l – Harmonic degree of moment.
n – Number of points in wavelength abscissa.
lam[n] – Wavelength abscissa (Å).
D[n-1] – Spectroscopic intensity D-moment (erg/cm^2/s/Å) in bins delineated by lam
stat – Status code (set to NULL if not required).
deriv_vec – Derivative flags (set to NULL if not required).
-
void interp_specgrid_flux(SpecGrid specgrid, double x_vec[], int n, double lam[], double F[], Stat *stat, bool deriv_vec[])
Interpolate the spectroscopic flux.
- Parameters
specgrid – Grid object.
x_vec – Photospheric parameter values.
n – Number of points in wavelength abscissa.
lam[n] – Wavelength abscissa (Å).
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 – Derivative flags (set to NULL if not required).
PhotGrid
-
void load_photgrid(const char *photgrid_file_name, PhotGrid *photgrid, int *stat)
Create a new
PhotGrid
by 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
PhotGrid
by 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_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, double x_vec[], double mu, double *I, int *stat, bool deriv_vec[])
Interpolate the photometric intensity, normalized to the zero-point flux.
- Parameters
photgrid – Grid object.
x_vec – Photospheric parameter values.
mu – Cosine of angle of emergence relative to surface normal.
I – Photometric intensity (/sr).
stat – Status code (set to NULL if not required).
deriv_vec – Derivative flags (set to NULL if not required).
-
void interp_photgrid_E_moment(PhotGrid photgrid, double x_vec[], int k, double *E, int *stat, bool deriv_vec[])
Interpolate the photometric intensity E-moment, normalized to the zero-point flux.
- Parameters
photgrid – Grid object.
x_vec – Photospheric parameter values.
k – Degree of moment.
D – Photometric intensity E-moment.
stat – Status code (set to NULL if not required).
deriv_vec – Derivative flags (set to NULL if not required).
-
void interp_photgrid_D_moment(PhotGrid photgrid, double x_vec[], int l, double *D, int *stat, bool deriv_vec[])
Interpolate the photometric intensity D-moment, normalized to the zero-point flux.
- Parameters
photgrid – Grid object.
x_vec – Photospheric parameter values.
l – Harmonic degree of moment.
D – Photometric intensity D-moment.
stat – Status code (set to NULL if not required).
deriv_vec – Derivative flags (set to NULL if not required).
-
void interp_photgrid_flux(PhotGrid photgrid, double x_vec[], double *F, int *stat, bool deriv_vec[])
Interpolate the photometric flux, normalized to the zero-point flux.
- Parameters
PhotGrid – Grid object.
x_vec – Photospheric parameter values.
F – Photometric flux.
stat – Status code (set to NULL if not required).
deriv_vec – Derivative flags (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).
Compiling/Linking
Headers for the C interface are provided in the header file
$MSG_DIR/include/cmsg.h
, and executables should be linked
against $MSG_DIR/lib/libcmsg.so
(Linux) or
$MSG_DIR/lib/libcmsg.dylib
(MacOS). To simplify this process,
a script $MSG_DIR/scripts/cmsg_link
is provided that writes
the appropriate linker commands to standard output. This script can be
used to compile/link a program with gcc as follows:
gcc -o myprogram myprogram.c -I $MSG_DIR/include `$MSG_DIR/scripts/cmsg_link`
Data Schema
This chapter specifies the data schema used by MSG data files. These files are stored on disk in HDF5 format, and the schema describes what data appear in each type of file.
Conventions
All MSG data files adhere to the following conventions:
HDF5 groups are used to store nested data structures that correspond reasonably closely to the derived data types (e.g.,
specgrid_t
) of the Fortran interface.Real values are written with HDF5 type H5T_IEEE_F64LE, or H5T_IEEE_F32LE when reduced precision is permitted.
Integer values are written with HDF5 type H5T_STD_I32LE.
Logical (boolean) values are written with HDF5 type H5T_STD_I32LE, with 1 corresponding to true, and 0 corresponding to false.
Character values are written with HDF5 type H5T_NATIVE_CHARACTER.
Files
- passband
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
group |
passband. |
- photgrid
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
group |
photgrid. |
- photint
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
group |
specint |
|
|
attribute |
real |
label value. |
- specgrid
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
group |
specgrid. |
- specint
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
group |
specint |
|
|
attribute |
real |
label value. |
Groups
- axis
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
character |
label. |
|
dataset |
real(:) |
abscissa values. |
- comp_range
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
dataset |
integer |
number of sub-ranges. |
|
group |
sub-ranges ( |
- cubint
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
integer |
number of points. |
|
dataset |
real( |
abscissa values. |
|
dataset |
real( |
ordinate values. |
|
dataset |
real( |
derivative values. |
- limb
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
character |
limb-darkening law (see here for a list of options). |
- limb_photint
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
logical |
precision of data (.TRUE. for 64-bit, .FALSE. for 32-bit). |
|
dataset |
real(:) |
intensity coefficients (/sr). |
|
group |
limb-darkening law. |
- limb_specint
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
real |
minimum wavelength (Å). |
|
attribute |
real |
maximum wavelength (Å). |
|
attribute |
logical |
storage precision of |
|
dataset |
real(:,:) |
intensity coefficients (erg/cm^2/s/Å/sr). |
|
group |
wavelength grid (Å). |
|
|
group |
limb-darkening law. |
- lin_range
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
real |
abscissa start value. |
|
attribute |
real |
abscissa spacing. |
|
attribute |
integer |
number of points. |
- log_range
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
real |
logarithm of abscissa start value. |
|
attribute |
real |
abscissa logarithmic spacing. |
|
attribute |
integer |
number of points. |
- passband
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
real |
normalizing flux (erg/cm^2/s). |
|
group |
filter interpolant. |
- photgrid
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
integer |
number of photints. |
|
group |
grid photints ( |
|
|
group |
virtual grid. |
- specgrid
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
integer |
number of specints. |
|
attribute |
real |
minimum wavelength (Å). |
|
attribute |
real |
maximum wavelength (Å). |
|
group |
grid specints ( |
|
|
group |
virtual grid. |
- tab_range
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
dataset |
real(:) |
abscissa values. |
- vgrid
Item |
Object Type |
Data Type |
Definition |
---|---|---|---|
|
attribute |
character |
literal |
|
attribute |
integer |
literal |
|
attribute |
integer( |
shape of grid. |
|
attribute |
integer |
rank of grid. |
|
group |
grid axes ( |
Grid Files
Due to their large size and gradually evolving content (as
improvements are made), specgrid files are
not shipped as part of the rhdtownsend/msg git repository; they
must be downloaded separately (the sole exception is the demo grid,
$MSG_DIR/data/grids/sg-demo.h5
). This chapter describes the
various grids currently available.
OSTAR2002 Grids
The OSTAR2002 grids are based on the line-blanketed O-star model atmospheres published in Lanz & Hubeny (2003). Spectra are calculated from these atmospheres using SYNSPEC, and their angle dependence parameterized with the CLARET limb-darkening law.
Specgrid File |
Resolution |
\(\lambda\) Range |
Photospheric Parameters |
---|---|---|---|
sg-OSTAR2002-low.h5 (50MB) |
\(\mathcal{R}=1000\) |
\([880\,\angstrom, 5.0\,\um]\) |
|
sg-OSTAR2002-medium.h5 (447MB) |
\(\mathcal{R}=10000\) |
\([880\,\angstrom, 5.0\,\um]\) |
|
sg-OSTAR2002-high.h5 (4.2GB) |
\(\mathcal{R}=100000\) |
\([880\,\angstrom, 5.0\,\um]\) |
|
The definitions and ranges of the photospheric parameters are as follows:
effective temperature
Teff
\(= \Teff/\kelvin \in [27500, 55000]\)surface gravity
log(g)
\(= \log_{10} (g/\cm\,\second^{-2}) \in [3.0, 4.75]\)metallicity
Z/Zo
\(= Z/Z_{\odot} \in [0.02, 2]\)
BSTAR2006 Grids
The BSTAR2006 grids are based on the line-blanketed B-star model atmospheres published in Lanz & Hubeny (2007). Spectra are calculated from these atmospheres using SYNSPEC, and their angle dependence parameterized with the CLARET limb-darkening law.
Specgrid File |
Resolution |
\(\lambda\) Range |
Photospheric Parameters |
---|---|---|---|
sg-BSTAR2006-low.h5 (77MB) |
\(\mathcal{R}=1000\) |
\([880\,\angstrom, 5\,\um]\) |
|
sg-BSTAR2006-medium.h5 (693MB) |
\(\mathcal{R}=10000\) |
\([880\,\angstrom, 5\,\um]\) |
|
sg-BSTAR2006-high.h5 (6.5GB) |
\(\mathcal{R}=100000\) |
\([880\,\angstrom, 5\,\um]\) |
|
The definitions and ranges of the photospheric parameters are as follows:
effective temperature
Teff
\(= \Teff/\kelvin \in [15000, 30000]\)surface gravity
log(g)
\(= \log_{10} (g/\cm\,\second^{-2}) \in [1.75, 4.75]\)metallicity
Z/Zo
\(= Z/Z_{\odot} \in [0, 2]\)
CAP18 Grids
The CAP18 grids are based on the data published in Allende Prieto et al. (2018) (the letters ‘CAP’ are the initials of the first author). The angle dependence of spectra is parameterized with the CONST limb-darkening law.
Specgrid File |
Resolution |
\(\lambda\) Range |
Photospheric Parameters |
---|---|---|---|
sg-CAP18-large.h5 (73GB) |
\(\mathcal{R}=10000\) |
\([1300\,\angstrom, 6.5\,\um]\) |
|
sg-CAP18-coarse.h5 (339MB) |
\(\mathcal{R}=10000\) |
\([1300\,\angstrom, 6.5\,\um]\) |
|
sg-CAP18-high.h5 (2.9GB) |
\(\mathcal{R}=100000\) |
\([1300\,\angstrom, 6.5\,\um]\) |
|
sg-CAP18-ultra.h5 (5.2GB) |
\(\mathcal{R}=300000\) |
\([1300\,\angstrom, 6.5\,\um]\) |
|
The definitions and ranges of the photospheric parameters are as follows:
effective temperature
Teff
\(= \Teff/\kelvin \in [3500, 30000]\)surface gravity
log(g)
\(= \log_{10} (g/\cm\,\second^{-2}) \in [0.0, 5.0]\)metallicity
[Fe/H]
\(= \log_{10}[ (\mathrm{Fe}/\mathrm{H}) / (\mathrm{Fe}/\mathrm{H})_{\odot} ] \in [-5.0, 0.5]\)alpha enhancement
[alpha/Fe]
\(= \log_{10}[ (\alpha/\mathrm{Fe}) / (\alpha/\mathrm{Fe})_{\odot} ] \in [-1.0, 1.0]\)microturbulent velocity
log(xi)
\(= \log_{10} (\xi/\cm\,\second^{-1}) \in [-0.301,0.903]\)
Göttingen Grids
The Göttingen grids are based on the data described in Husser et al. (2013) and available for download from phoenix.astro.physik.uni-goettingen.de. The angle dependence of spectra is parameterized with the CONST limb-darkening law.
Specgrid File |
Resolution |
\(\lambda\) Range |
Photospheric Parameters |
---|---|---|---|
sg-Goettingen-HiRes.h5 (116GB) |
variable |
\([500\,\angstrom, 5.5\um]\) |
|
sg-Goettingen-MedRes-A.h5 (6.0GB) |
\(\Delta \lambda = 1\,\angstrom\) |
\([3000\,\angstrom, 1.0\um]\) |
|
sg-Goettingen-MedRes-R.h5 (18GB) |
\(\mathcal{R}=10,000\) |
\([3000\,\angstrom, 2.5\,\um]\) |
|
The definitions and ranges of the photospheric parameters are as follows:
effective temperature
Teff
\(= \Teff/\kelvin \in [2\,300, 12\,000]\)surface gravity
log(g)
\(= \log_{10}(g/\cm\,\second^{-2}) \in [0.0, 6.0]\)metallicity
[Fe/H]
\(= \log_{10}[ (\mathrm{Fe}/\mathrm{H}) / (\mathrm{Fe}/\mathrm{H})_{\odot} ] \in [-4.0, 1.0]\)alpha enhancement
[alpha/H]
\(= \log_{10}[ (\alpha/\mathrm{Fe}) / (\alpha/\mathrm{Fe})_{\odot} ] \in [-0.2,1.2]\)
Passband Files
As with the Grid Files, the passband files used by MSG are not shipped as part of the rhdtownsend/msg git repository; they must be downloaded separately. This chapter describes the various files available for download.
SVO Filter Service Passbands
The Spanish Virtual Observatory (SVO) offers a database of filters and calibrations for a wide array of photometric systems. The table below provides links to tar archives of passband files created for selected systems in the SVO database using the make_passband tool. For each filter of each system, separate passband files for the Vega, AB and ST magnitude systems are provided.
Passband File |
Facility |
Instrument |
Filters |
---|---|---|---|
Corot |
— |
AS, PF |
|
GAIA |
Gaia2m |
Gbp_faint, Gbp_bright, G, Grp |
|
GAIA |
GAIA2r |
Gbp, G, Grp |
|
GAIA |
GAIA2 |
Gbp, G, Grp |
|
GAIA |
GAIA3 |
Gbp, G, Grp |
|
GAIA |
GAIA0 |
Gbp, G, Grp |
|
GALEX |
— |
FUV, NUV |
|
Generic |
Johnson |
U, B, V, R, I, J, M |
|
Geneva |
— |
U, B1, B, B2, V2, V, G |
|
Hipparcos |
— |
Hp, Hp_bes, Hp_MvB |
|
HST |
STIS_FUV |
F25LYA, F25LYA_G140L, F25LYA_G140M, 25MAMA_G140M, F25ND5_G140M, 25MAMA_G140L, F25ND3_G140M, F25ND5_G140L, 25MAMA, F25ND3_G140L, F25ND5, F25ND3, F25NDQ2_G140M, F25NDQ3_G140M, F25NDQ2_G140L, F25NDQ1_G140M, F25NDQ3_G140L, F25NDQ1_G140L, F25NDQ2, F25NDQ4_G140M, F25NDQ3, F25NDQ1, F25SRF2_G140M, F25NDQ4_G140L, F25SRF2_G140L, F25NDQ4, F25SRF2, F25QTZ_G140M, F25QTZ_G140L, F25QTZ |
|
HST |
ACS_SBC |
F122M, F115LP, PR130L, PR110L, F125LP, F140LP, F150LP, F165LP |
|
HST |
WFPC1-WF |
G200M2, F157W, F194W, F230W, F122M, F284W, F336W, F368M, F375N, F413M, F439W, F437N, G200, F469N, F487N, F492M, G450, F502N, F517N, F555W, F547M, F569W, F588N, F606W, F8ND, POL0, POL120, POL60, F128LP, F622W, F631N, F648M, F656N, F658N, F664N, F673N, F675W, F702W, F718M, G800, F791W, F814W, F725LP, F875M, F889N, F785LP, F850LP, F1042M, F1083N |
|
HST |
WFPC1-PC |
G200M2, F157W, F194W, F230W, F284W, F336W, F368M, F375N, F413M, F439W, F437N, F469N, F122M, F487N, F492M, F502N, G450, G200, F517N, F555W, F547M, F569W, F588N, F606W, F622W, F8ND, POL0, POL120, POL60, F128LP, F631N, F648M, F656N, F658N, F664N, F673N, F675W, F702W, F718M, G800, F791W, F814W, F725LP, F875M, F889N, F785LP, F850LP, F1042M, F1083N |
|
HST |
WFPC2-WF |
F122M, F160BW, F157W, F170W, F185W, F218W, F255W, F300W, F336W, F343N, F375N, FQUVN33, FQUVN_B, F390N, FQUVN_C, FQUVN_D, F380W, F410M, F439W, F437N, F450W, F467M, F469N, F487N, F502N, F555W, FQCH4N_D, F547M, F569W, F588N, F606W, FQCH4N33, F622W, F631N, F165LP, F130LP, F656N, F658N, F673N, F675W, F702W, POLQ, POLQ_90, POLQ_45, FQCH4N_B, F814W, F791W, F785LP, FQCH4N_C, F850LP, F953N, F1042M |
|
HST |
WFPC2-PC |
F122M, F160BW, F157W, F170W, F185W, F218W, F255W, F300W, F336W, F343N, F375N, FQUVN, F390N, F380W, F410M, F439W, F437N, F450W, F467M, F469N, F487N, F502N, F555W, FQCH4N, F547M, F569W, F588N, F606W, F622W, F631N, F165LP, F130LP, F656N, F658N, F673N, F675W, F702W, POLQ, F814W, F791W, F785LP, F850LP, F953N, F1042M |
|
HST |
HSP_UV1 |
F122M_B, F122M_A, F135W_A, F135W_B, F145M_A, F145M_B, PRISM_BLUE, F152M_A, F152M_B, F184W_A, F184W_B, F218M_A, F218M_B, F220W_A, F220W_B, F140LP_A, F140LP_B, F240W_A, F240W_B, F248M_A, F248M_B, PRISM_RED, F278N_A, F278N_B |
|
HST |
HSP_UV2 |
F122M_B, F122M_A, F145M_A, F145M_B, PRISM_BLUE, F152M_A, F152M_B, F179M_A, F179M_B, F184W_A, F184W_B, F218M_A, F218M_B, F140LP_A, F160LP_A, F140LP_B, F160LP_B, F248M_A, F248M_B, F262M_A, F262M_B, PRISM_RED, F278N_A, F278N_B, F284M_A, F284M_B |
|
HST |
FOC_F48 |
F140W, F150W, F175W, F195W, F220W, F275W, PRISM3, F130LP, F342W, PRISM1, F180LP, PRISM2, F305LP, F430W |
|
HST |
FOC_F96 |
F120M, F130M, F140M, F140W, F152M, F170M, F175W, F165W, F190M, F195W, F210M, F220W, F231M, F253M, F275W, F278M, F307M, F320W, F6ND, F2ND, F342W, F1ND, PRISM1, F346M, F130LP, POL120, POL0, F4ND, PRISM2, POL60, F8ND, F372M, F410M, F430W, F370LP, F437M, F470M, F486N, F502M, F501N, F480LP, F550M, F600M, F630M |
|
HST |
HSP_VIS |
F184W_A, F184W_B, PRISM_BLUE, F240W_A, F240W_B, F262M_A, F262M_B, F355M_A, F355M_B, F160LP_A, F160LP_B, F419N_A, F419N_B, F450W_A, F450W_B, F400LP_A, F400LP_B, F551W_B, F551W_A, PRISM_RED, F620W_B, F620W_A |
|
HST |
STIS_NUV |
F25CN182, F25CN182_PRISM, F25CIII_PRISM, F25CIII, F25CIII_G230L, F25CIII_G230M, F25CN182_G230L, F25CN182_G230M, 25MAMA, F25NDQ1, 25MAMA_PRISM, F25QTZ, F25SRF2, F25QTZ_PRISM, F25SRF2_PRISM, F25NDQ1_PRISM, 25MAMA_G230L, F25QTZ_G230L, F25SRF2_G230L, 25MAMA_G230M, F25NDQ1_G230L, F25QTZ_G230M, F25SRF2_G230M, F25NDQ2, F25NDQ1_G230M, F25NDQ2_PRISM, F25NDQ2_G230L, F25ND3, F25NDQ2_G230M, F25ND3_PRISM, F25ND3_G230L, F25ND3_G230M, F25NDQ3, F25NDQ3_G230L, F25CN270_G230L, F25CN270, F25CN270_PRISM, F25CN270_G230M, F25NDQ3_PRISM, F25NDQ3_G230M, F25NDQ4, F25NDQ4_PRISM, F25MGII, F25MGII_PRISM, F25MGII_G230L, F25MGII_G230M, F25NDQ4_G230L, F25NDQ4_G230M, F25ND5, F25ND5_PRISM, F25ND5_G230L, F25ND5_G230M |
|
HST |
WFC3_UVIS2 |
F218W, FQ232N, F225W, FQ243N, G280, F275W, F300X, F280N, F336W, F343N, F373N, FQ378N, FQ387N, F390M, F390W, F395N, F410M, FQ422M, F438W, FQ436N, FQ437N, F467M, F469N, F475W, F487N, F200LP, F475X, FQ492N, F502N, FQ508N, F555W, F547M, FQ575N, F350LP, F606W, FQ619N, F621M, F625W, F631N, FQ634N, F645N, F656N, F657N, F658N, F665N, FQ672N, FQ674N, F673N, F680N, F689M, F600LP, FQ727N, FQ750N, F763M, F775W, F814W, F845M, FQ889N, FQ906N, F850LP, FQ924N, FQ937N, F953N |
|
HST |
WFC3_UVIS1 |
F218W, FQ232N, F225W, FQ243N, F275W, G280, F300X, F280N, F336W, F343N, F373N, FQ378N, FQ387N, F390M, F390W, F395N, F410M, FQ422M, F438W, FQ436N, FQ437N, F467M, F469N, F475W, F487N, F475X, FQ492N, F502N, FQ508N, F200LP, F555W, F547M, FQ575N, F350LP, F606W, FQ619N, F621M, F625W, F631N, FQ634N, F645N, F656N, F657N, F658N, F665N, FQ672N, FQ674N, F673N, F680N, F689M, F600LP, FQ727N, FQ750N, F763M, F775W, F814W, F845M, FQ889N, FQ906N, F850LP, FQ924N, FQ937N, F953N |
|
HST |
ACS_HRC |
F220W, F250W, F330W, F344N, FR388N, F435W, FR459M, F475W, F502N, FR505N, F555W, F550M, F606W, PR200L, F625W, FR656N, F658N, F660N, POL_UV, POL_V, G800L, F775W, F814W, F892N, F850LP, FR914M |
|
HST |
HSP_POL |
F216M_0, F237M_0, F277M_0, F327M_0, F160LP_T, F160LP_A |
|
HST |
STIS_CCD |
F28X50LP_G230LB, 50CCD_G230LB, 50CORON_G230LB, F28X50LP_G230MB, 50CCD_G230MB, 50CORON_G230MB, F28X50OII_G430L, F28X50OII_G430M, F28X50OII, 50CCD_G430M, 50CORON_G430M, 50CCD_G430L, 50CORON_G430L, F28X50OIII, F28X50OIII_G430L, F28X50OIII_G430M, F28X50LP_G430M, F28X50LP_G430L, 50CCD, 50CORON, 50CCD_G750L, 50CORON_G750L, F28X50LP, 50CCD_G750M, 50CORON_G750M, F28X50LP_G750L, F28X50LP_G750M |
|
HST |
FOS_BLUE |
G130H, MIRROR, G190H, G160L, G400H, G270H, PRISM |
|
HST |
ACS_WFC |
FR388N, FR423N, F435W, FR459M, FR462N, F475W, F502N, FR505N, F555W, FR551N, F550M, F606W, FR601N, F625W, FR647M, FR656N, F658N, F660N, POL_UV, POL_V, FR716N, G800L, F775W, FR782N, F814W, FR853N, F850LP, F892N, FR914M, FR931N, FR1016N |
|
HST |
FOS_RED |
G780H, G160L, G190H, MIRROR, G400H, G270H, PRISM, G570H, G650L |
|
HST |
FGS |
ND5, PUPIL, F583W, F605W, F550W, F650W |
|
HST |
NICMOS1 |
F090M, F095N, F097N, POL0S, POL240S, POL120S, F108N, F110M, F113N, F110W, F140W, F145M, F160W, F165M, F164N, F166N, F170M, F187N, F190N |
|
HST |
WFC3_IR |
F098M, G102, F105W, F110W, F125W, F126N, F127M, F128N, F130N, F132N, F139M, G141, F140W, F153M, F160W, F164N, F167N |
|
HST |
NICMOS3 |
G096, F108N, F110W, F113N, F150W, G141, F160W, F164N, F166N, F175W, F187N, F190N, F196N, F200N, F212N, F215N, G206, F222M, F240M |
|
HST |
NICMOS2 |
F110W, F160W, F165M, F171M, F180M, F187W, F187N, F190N, POL120L, POL240L, POL0L, F204M, F205W, F207M, F212N, F215N, F216N, F222M, F237M |
|
JWST |
NIRCam |
F070W, F090W, F115W, F140M, F150W, F162M, F164N, F150W2, F182M, F187N, F200W, F210M, F212N, F250M, F277W, F300M, F323N, F322W2, F335M, F356W, F360M, F405N, F410M, F430M, F444W, F460M, F466N, F470N, F480M |
|
JWST |
NIRISS |
F090W, F115W, F140M, F150W, F158M, F200W, F277W, F356W, F380M, F430M, F444W, F480M |
|
JWST |
MIRI |
F560W, F770W, F1000W, F1065C, F1130W, F1140C, F1280W, F1500W, F1550C, F1800W, F2100W, F2300C, F2550W |
|
Kepler |
— |
K |
|
LSST |
— |
u_filter, u, g_filter, g, r, r_filter, i, i_filter, z, z_filter, y, y_filter |
|
MOST |
— |
band |
|
PAN-STARRS |
PS1 |
g, r, w, open, i, z, y |
|
Spitzer |
IRAC |
I1, I2, I3, I4 |
|
Spitzer |
IRS |
Blue, Red |
|
Spitzer |
MIPS |
24mu, 70mu, 160mu |
|
Generic |
Strömgren |
u, v, b, y |
|
TESS |
— |
Red |
|
WFIRST |
WFI |
R062, Z087, Y106, Prism, J129, W146, Grism, H158, F184 |
Grid Tools
This appendix outlines the set of tools provided with MSG to assist in
creating and managing custom grids. These tools are built by default
during compilation, and can be found in the $MSG_DIR/bin
directory.
Converting Spectra
The following tools convert existing spectra and/or spectroscopic grids into one or more specint files:
synspec_to_specint
The synspec_to_specint tool extracts an single intensity
spectrum from a fort.18
data file produced by the SYNSPEC
spectral synthesis package (Lanz & Hubeny, 2003), and writes it to a
specint file. It accepts the following
command-line arguments:
- <synspec_file_name>
Name of input SYNSPEC file.
- <n_mu>
Number of \(\mu\) values in input file (as specified by the
nmu
parameter in thefort.55
SYNSPEC control file).
- <mu_0>
Minimum \(\mu\) value in input file (as specified by the
ang0
parameter in thefort.55
SYNSPEC control file).
- <lam_min>
Minimum wavelength of output file (\(\angstrom\)).
- <lam_max>
Maximum wavelength in output file (\(\angstrom\)).
- <R>
Resolution \(\mathcal{R}=\lambda/\Delta\lambda\) in output file.
- <label> (optional)
Label of atmosphere parameter (must be accompanied by a corresponding
<value>
argument).
- <value> (optional)
Value of atmosphere parameter (must be accompanied by a corresponding
<label>
argument).
Note that <label>
and <value>
parameters must be
paired; and that there can be multiple of these pairs. For the law
selected by the <law_str>
option, the tool calculates the
limb-darkening coefficients at each wavelength via a least-squares fit
to the function
ferre_to_specint
The ferre_to_specint tool extracts a series of flux spectra from a data file in FERRE format, and writes them to specint files. It accepts the following command-line arguments:
- <ferre_file_name>
Name of input FERRE file.
- <ferre_file_type>
Type of input file. This determines the mapping between photospheric parameters given in the input file, and photospheric parameters written to the output file. Supported options are
CAP18
(for the Allende Prieto et al., 2018 grids).
goettingen_to_specint
The goettingen_to_specint tool extracts a flux spectrum from a data file in FITS format (with the schema described by Husser et al., 2013), and writes it to a specint file. This tool accepts the following command-line arguments:
- <fits_file_name>
Name of input FITS file.
- <wave_type>
Type of wavelength abscissa. This determines the number and distribution of points to assume for the input file. Supported options, corresponding to the different grids described by Husser et al. (2013), are
HiRes
(high-resolution),MedRes-A1
(medium-resolution, \(\Delta \lambda = 1\,\angstrom\)) andMedRes-R10000
(medium resolution, \(\mathcal{R}=10\,000\)).
Note
In order for goettingen_to_specint to build, you must
first uncomment/edit the line in $MSG_DIR/build/Makefile
that defines the FITS_LDFLAGS variable. This variable defines the
flags used to link against your system’s FITS library.
ascii_to_specint
The ascii_to_specint tool reads a generic flux spectrum from an ASCII text file, and writes it to a specint file. It accepts the following command-line arguments:
- <ascii_file_name>
Name of input ASCII text file (see below).
- <lam_units>
Units of wavelength data in input file. Possible choices are
'A'
(\(\angstrom\)).
- <flux_units>
Units of flux data in input file. Possible choices are
'erg/cm^2/s/A'
(\(\erg\,\cm^{-2}\,\second^{-1}\,\angstrom^{-1}\)).
- <label> (optional)
Label of atmosphere parameter (must be accompanied by a corresponding
<value>
argument).
- <value> (optional)
Value of atmosphere parameter (must be accompanied by a corresponding
<label>
argument).
Note that <label>
and <value>
parameters must be
paired; and that there can be multiple of these pairs. The input ASCII
text file should contain one or more lines, each consisting of a
wavelength value followed by a flux value.
specint_to_specint
The specint_to_specint tool subsets and/or rebins data in an existing specint file It accepts the following command-line arguments:
- lam_min=<value> (optional)
Subset to have a minimum wavelength of at least <value> (\(\angstrom\)).
- lam_max=<value> (optional)
Subset to have a maximum wavelength of at most <value> (\(\angstrom\)).
- R=<value> (optional)
Rebin to have a uniform resolution \(\mathcal{R}\) of <value>.
- dlam=<value> (optional)
Rebin to have a uniform wavelength spacing \(\Delta \lambda\) of <value> (\(\angstrom\)).
- just=<L|R> (optional)
Justify the new wavelength abscissa to the left (
L
) or right (R
).
If you would like to convert spectra from a format that isn’t covered by these tools, please open an issue describing the format and/or pointing to relevant literature.
Packaging Grids
The following tools package multiple specint files to create specgrid and photgrid files:
specint_to_specgrid
The specint_to_specgrid tool packages together multiple specint files to make a spectroscopic grid, writing it to a specgrid. It accepts the following command-line arguments:
- <manifest_file_name>
Name of input manifest file (see below).
- <allow_dupes> (optional)
Flag governing handling of duplicate grid nodes in the manifest file; set to
T
to allow duplicates.
The manifest file is a simple text file that lists all the
specint files (one per line) that should be included in
the grid. The axes and the topology of the grid are automatically
determined by the labels attached to each specint file. If
two files have identical labels, then the one appearing later in the
manifest file is used if <allow_dupes>
is T
; otherwise, an
error is thrown.
specgrid_to_photgrid
The specgrid_to_photgrid tool convolves a spectroscopic grid with a passband to make a photometric grid, writing it to a photgrid file. It accepts the following command-line arguments:
Note that it’s not always necessary to create photgrid files, as MSG can convolve with passbands on the fly (as discussed in the Photometric Colors section).
Creating Passband Files
The following tools create passband files:
make_passband
The make_passband tool reads response data, and writes them to a passband file. It accepts the following command-line arguments:
- <input_file_name>
Name of input file (see below).
- <F_0>
Normalizing flux \(F_{0}\) (\(\erg\,\cm^{-2}\,\second^{-1}\,\angstrom^{-1}\)).
The input file is a text file tabulating wavelength \(\lambda\) (\(\angstrom\)) and passband response function \(S'(\lambda)\) (see the Photometric Colors section).
Tensor Product Interpolation
Tensor product interpolation is a generalization of low-dimension (typically, 2-d or 3-d) multivariate interpolation schemes to an arbitrary number of dimensions. The literature on this topic is rather specialized (although see this Wikipedia section); therefore, this appendix provides a brief introduction. Interpolation in general is covered by many textbooks, so the focus here is reviewing how simple schemes can be generalized.
Univariate Interpolation
Suppose we are given a set of values \(f(x_{1}), f(x_{2}), \ldots, f(x_{n})\), representing some function \(f(x)\) evaluated at a set of \(n\) discrete grid points. Then, a piecewise-linear interpolation scheme approximates \(f(x)\) in each subinterval \(x_{i} \leq x \leq x_{i+1}\) via
Here \(f_{i} \equiv f(x_{i})\), while
is a normalized coordinate that expresses where in the subinterval we are; \(u=0\) corresponds to \(x=x_{i}\), and \(u=1\) to \(x=x_{i+1}\). The linear basis functions are defined by
This interpolation scheme is \(C^{0}\) continuous, and reproduces \(f(x)\) exactly at the grid points.
If we want a smoother representation of the function, we generally need more data. Suppose then that we’re also provided the derivatives \(\sderiv{f}{x}\) at the grid points. Then, a piecewise-cubic interpolation scheme is
where \(u\) has the same definition as before, \(\partial_{x} f_{i} \equiv (\sderiv{f}{x})_{x=x_{i}}\), and the cubic basis functions are
(these can be recognized as the basis functions for cubic Hermite splines). This new definition is \(C^{1}\) continuous, and reproduces \(f(x)\) and its first derivative exactly at the grid points.
Bivariate Interpolation
Bivariate interpolation allows us to approximate a function \(f(x,y)\) from its values on a rectilinear (but not necessarily equally spaced) grid described by the axis values \(x_{1},x_{2},\ldots,x_{n}\) and \(y_{1},y_{2},\ldots,y_{m}\). A piecewise-bilinear interpolating scheme approximates \(f(x,y)\) in each subinterval \(x_{i} \leq x \leq x_{i+1}\), \(y_{j} \leq y \leq y_{j+1}\) via
Here, \(u\) has the same definition as before, while
and \(f_{i,j} \equiv f(x_{i},y_{j})\). We can also write the scheme in the more-compact form
where the coefficients \(\mathcal{L}^{p,q}\) can be expressed as the matrix
A corresponding piecewise-bicubic interpolating scheme can be written in the form
where the coefficients \(\mathcal{C}^{h}\) can be expressed as the matrix
Constructing this matrix requires 16 values: the function at the four corners of the subinterval, the first derivatives with respect to \(x\) and with respect to \(y\) at the corners, and the cross derivatives \(\partial_{xy} f\) at the corners.
Multivariate Interpolation
Let’s now generalize the compact expressions above for piecewise-bilinear and bicubic interpolation, to piecewise interpolation in \(N\) dimensions. For linear interpolation, we have
Likewise, for cubic interpolation, we have
The coefficients \(\mathcal{L}^{p_{1},p_{2},\ldots,p_{N}}\) and \(\mathcal{C}^{p_{1},p_{2},\ldots,p_{N}}\) cannot easily be expressed in closed form, but they are relatively straightforward to construct algorithmically.
The summations in expressions above can be regarded as the contraction (over all indices) of a pair of rank-\(N\) tensors. In the cubic case, the components of the first tensor correspond to the coefficients \(\mathcal{C}^{p_{1},p_{2},\ldots,p_{N}}\), while the second tensor is formed by taking \(N\) outer products between the vectors
Hence, this kind of multivariate interpolation is known as tensor product interpolation.