WDmodel.WDmodel module

DA White Dwarf Atmosphere Models and SED generator.

Model grid originally from J. Holberg using I. Hubeny’s Tlusty code (v200) and custom Synspec routines, repackaged into HDF5 by G. Narayan.

class WDmodel.WDmodel.WDmodel(grid_file=None, grid_name=None, rvmodel='f99')[source]

Bases: object

DA White Dwarf Atmosphere Model and SED generator

Base class defines the routines to generate and work with DA White Dwarf model spectra. Requires a grid file of DA White Dwarf atmospheres. This grid file is included along with the package - TlustyGrids.hdf5 - and is the default.

Parameters
  • grid_file (str, optional) – Filename of the HDF5 grid file to read. See WDmodel.io.read_model_grid() for format of the grid file. Default is TlustyGrids.hdf5, included with the WDmodel package.

  • grid_name (str, optional) – Name of the HDF5 group containing the white dwarf model atmosphere grids in grid_file. Default is default

  • rvmodel ({'ccm89','od94','f99','custom'}, optional) –

    Specify parametrization of the reddening law. Default is 'f99'.

    rvmodel

    parametrization

    'ccm89'

    Cardelli, Clayton and Mathis (1989, ApJ, 345, 245)

    'od94'

    O’Donnell (1994, ApJ, 422, 158)

    'f99'

    Fitzpatrick (1999, PASP, 111, 63)

    'custom'

    Custom law from Jay Holberg (email, 20180424)

_lines

dictionary mapping Hydrogen Balmer series line names to line number, central wavelength in Angstrom, approximate line width and continuum region width around line. Used to extract Balmer lines from spectra for visualization.

Type

dict

_grid_file

Filename of the HDF5 grid file that was read.

Type

str

_grid_name

Name of the HDF5 group containing the white dwarf model atmosphere

Type

str

_wave

Array of model grid wavelengths in Angstroms, sorted in ascending order

Type

array-like

_ggrid

Array of model grid surface gravity values in dex, sorted in ascending order

Type

array-like

_tgrid

Array of model grid temperature values in Kelvin, sorted in ascending order

Type

array-like

_nwave

Size of the model grid wavelength array, _wave

Type

int

_ngrav

Size of the model grid surface gravity array, _ggrid

Type

int

_ntemp

Size of the model grid temperature array, _tgrid

Type

int

_flux

Array of model grid fluxes, shape (_nwave, _ntemp, _ngrav)

Type

array-like

_lwave

Array of model grid log10 wavelengths for interpolation

Type

array-like

_lflux

Array of model grid log10 fluxes for interpolation, shape (_ntemp, _ngrav, _nwave)

Type

array-like

_law
Type

extinction function corresponding to rvmodel

Returns

out

Return type

WDmodel.WDmodel.WDmodel instance

Raises

ValueError – If the supplied rvmodel is unknown

Notes

Virtually none of the attributes should be used directly since it is trivially possible to break the model by redefining them. Access to them is best through the functions connected to the models.

A custom user-specified grid file can be specified. See WDmodel.io.read_model_grid() for the format of the grid file.

Uses scipy.interpolate.RegularGridInterpolator to interpolate the models.

The class contains various convenience methods that begin with an underscore (_) that will not be imported by default. These are intended for internal use, and do not have the sanity checking and associated overhead of the public methods.

_WDmodel__init__rvmodel(rvmodel='f99')
_WDmodel__init__tlusty(grid_file=None, grid_name=None)
__call__(teff, logg, wave=None, log=False, strict=True)

Returns the model flux given teff and logg at wavelengths wave

Wraps WDmodel.WDmodel.WDmodel._get_model() adding checking of inputs.

Parameters
  • teff (float) – Desired model white dwarf atmosphere temperature (in Kelvin)

  • logg (float) – Desired model white dwarf atmosphere surface gravity (in dex)

  • wave (array-like, optional) – Desired wavelengths at which to compute the model atmosphere flux. If not supplied, the full model wavelength grid is returned.

  • log (bool, optional) – Return the log10 flux, rather than the flux

  • strict (bool, optional) – If strict, teff and logg out of model grid range raise a ValueError, otherwise raise a RuntimeWarning and set teff, logg to the nearest grid value.

Returns

  • wave (array-like) – Valid output wavelengths

  • flux (array-like) – Interpolated model flux at teff, logg and wavelengths wave

Raises

ValueError – If teff or logg are out of range of model grid and strict is True or if there are any invalid wavelengths, or the requested wavelengths to do not overlap with the model grid

Notes

Unlike the corresponding private methods, the public methods implement checking of the inputs and returns the wavelengths in addition to the flux. Internally, we only use the private methods as the inputs only need to be checked once, and their state is not altered anywhere after.

__init__(grid_file=None, grid_name=None, rvmodel='f99')[source]

Initialize self. See help(type(self)) for accurate signature.

_custom_extinction(wave, av, rv=3.1, unit='aa')[source]

Return the extinction for av, rv at wavelengths wave for the custom reddening law defined by J. Holberg

Mimics the interface provided by WDmodel.WDmodel.WDmodel._law to calculate the extinction as a function of wavelength (in Angstroms), \(A_{\lambda}\).

Parameters
  • wave (array-like) – Array of wavelengths in Angstrom at which to compute extinction, sorted in ascending order

  • av (float) – Extinction in the V band, \(A_V\)

  • rv (float, optional) – Fixed to 3.1 for J. Holberg’s custom reddening law

Returns

out – Extinction at wavelengths wave for av and rv

Return type

array-like

Notes

av should be >= 0.

_extract_from_indices(w, f, ZE, df=None)[source]

Extracts slices of multiple arrays for the same set of indices.

Convenience function to extract elements of wavelength w, flux f and optionally flux uncertainty df using indices ZE

Parameters
  • w (array-like) – Wavelength array from which to extract indices ZE

  • f (array-like) – Flux array from which to extract indices ZE

  • ZE (array-like) – indices to extract

  • df (None or array-like, optional) – If array-like, extracted elements of this array are also returned

Returns

  • w (array-like) – elements of input wavelength array at indices ZE

  • f (array-like) – elements of input flux array at indices ZE

  • [df] (array-like) – elements of input flux uncertainty array at indices ZE if optional input df is supplied

_extract_spectral_line(w, f, line, df=None)[source]

Extracts slices of multiple arrays corresponding to a hydrogen Balmer line

Convenience function to extract elements of wavelength w, flux f and optionally flux uncertainty df for a hydrogen Balmer line. Wraps WDmodel.WDmodel.WDmodel._get_line_indices() and WDmodel.WDmodel.WDmodel._extract_from_indices(), both of which have their own reasons for existence as well.

Parameters
  • w (array-like) – Wavelength array from which to extract elements corresponding to hydrogen Balmer line

  • f (array-like) – Flux array from which to extract elements corresponding to hydrogen Balmer line

  • line ({'alpha', 'beta', 'gamma', 'delta', 'zeta', 'eta'}) – Name of hydrogen Balmer line to extract. Properties are pre-defined in WDmodel.WDmodel.WDmodel._lines

  • df (None or array-like, optional) – If array-like, extracted elements of this array are also returned

Returns

  • w (array-like) – elements of input wavelength array for hydrogen Balmer feature line

  • f (array-like) – elements of input flux array for hydrogen Balmer feature line

  • [df] (array-like) – elements of input flux uncertainty array for hydrogen Balmer feature line if optional input df is supplied

Notes

Same as WDmodel.WDmodel.WDmodel.extract_spectral_line() without checking of inputs and therefore corresponding overhead. Used internally.

_get_full_obs_model(teff, logg, av, fwhm, wave, rv=3.1, log=False, pixel_scale=1.0)[source]

Returns the observed model flux given teff, logg, av, rv, fwhm (for Gaussian instrumental broadening) at wavelengths, wave as well as the full SED.

Convenience function that does the same thing as WDmodel.WDmodel.WDmodel._get_obs_model(), but also returns the full SED without any instrumental broadening applied, appropriate for synthetic photometry.

Uses WDmodel.WDmodel.WDmodel._get_model() to get the unreddened model, and reddens it with WDmodel.WDmodel.WDmodel.reddening() and convolves it with a Gaussian kernel using scipy.ndimage.filters.gaussian_filter1d()

Parameters
  • teff (float) – Desired model white dwarf atmosphere temperature (in Kelvin)

  • logg (float) – Desired model white dwarf atmosphere surface gravity (in dex)

  • av (float) – Extinction in the V band, \(A_V\)

  • fwhm (float) – Instrumental FWHM in Angstrom

  • wave (array-like) – Desired wavelengths at which to compute the model atmosphere flux.

  • rv (float, optional) – The reddening law parameter, \(R_V\), the ration of the V band extinction \(A_V\) to the reddening between the B and V bands, \(E(B-V)\). Default is 3.1, appropriate for stellar SEDs in the Milky Way.

  • log (bool, optional) – Return the log10 flux, rather than the flux (what’s actually interpolated)

  • pixel_scale (float, optional) – Jacobian of the transformation between wavelength in Angstrom and pixels. In principle, this should be a vector, but virtually all spectral reduction packages resample the spectrum onto a uniform wavelength scale that is close to the native pixel scale of the spectrograph. Default is 1.

Returns

  • flux (array-like) – Interpolated model flux at teff, logg with reddening parametrized by av, rv and broadened by a Gaussian kernel defined by fwhm at wavelengths wave

  • mod (numpy.recarray with dtype=[('wave', '<f8'), ('flux', '<f8')]) – Full model SED at teff, logg with reddening parametrized by av, rv

Notes

fwhm and pixel_scale must be > 0

classmethod _get_indices_in_range(wave, WA, WB, W0=None)[source]

Returns indices of wavelength between blue and red wavelength limits and the central wavelength

Parameters
  • wave (array-like) – Wavelengths array from which to extract indices

  • WA (float) – blue limit of wavelengths to extract

  • WB (float) – red limit of wavelenghts to extract

  • W0 (float or None, optional) – None or a central wavelength of range [WA, WB] to return. If None, the central wavelength is computed, else the input is simply returned.

Returns

  • W0 (float) – central wavelength of range [WA, WB]

  • ZE (array-like) – indices of wave in range [WA, WB]

_get_line_indices(wave, line)[source]

Returns the central wavelength and indices of wavelength corresponding to a hydrogen Balmer line

Parameters
  • wave (array-like) – Wavelengths array from which to extract indices

  • line ({'alpha', 'beta', 'gamma', 'delta', 'zeta', 'eta'}) – Name of hydrogen Balmer line to extract. Properties are pre-defined in WDmodel.WDmodel.WDmodel._lines

Returns

  • W0 (float) – central wavelength of line

  • ZE (array-like) – indices of wave of line

Notes

No checking of input - will throw KeyError if line is not accepted value

_get_model(teff, logg, wave=None, log=False)[source]

Returns the model flux given teff and logg at wavelengths wave

Simple 3-D interpolation of model grid. Computes unreddened, unnormalized, unconvolved, interpolated model flux. Uses scipy.interpolate.RegularGridInterpolator to generate the interpolated model. This output has been tested against WDmodel.WDmodel.WDmodel._get_model_nosp().

Parameters
  • teff (float) – Desired model white dwarf atmosphere temperature (in Kelvin)

  • logg (float) – Desired model white dwarf atmosphere surface gravity (in dex)

  • wave (array-like, optional) – Desired wavelengths at which to compute the model atmosphere flux. If not supplied, the full model wavelength grid is returned.

  • log (bool, optional) – Return the log10 flux rather than the flux.

Returns

flux – Interpolated model flux at teff, logg and wavelengths wave.

Return type

array-like

Notes

Inputs teff, logg and wave must be within the bounds of the grid. See WDmodel.WDmodel.WDmodel._wave, WDmodel.WDmodel.WDmodel._ggrid, WDmodel.WDmodel.WDmodel._tgrid, for grid locations and limits.

_get_model_nosp(teff, logg, wave=None, log=False)[source]

Returns the model flux given teff and logg at wavelengths wave

Simple 3-D interpolation of model grid. Computes unreddened, unnormalized, unconvolved, interpolated model flux. Not used, but serves as check of output of interpolation of scipy.interpolate.RegularGridInterpolator output.

Parameters
  • teff (float) – Desired model white dwarf atmosphere temperature (in Kelvin)

  • logg (float) – Desired model white dwarf atmosphere surface gravity (in dex)

  • wave (array-like, optional) – Desired wavelengths at which to compute the model atmosphere flux. If not supplied, the full model wavelength grid is returned.

  • log (bool, optional) – Return the log10 flux, rather than the flux.

Returns

flux – Interpolated model flux at teff, logg and wavelengths wave

Return type

array-like

Notes

Inputs teff, logg and wave must be within the bounds of the grid. See WDmodel.WDmodel.WDmodel._wave, WDmodel.WDmodel.WDmodel._ggrid, WDmodel.WDmodel.WDmodel._tgrid, for grid locations and limits.

This restriction is not imposed here for performance reasons, but is implicitly set by routines that call this method. The user is expected to verify this condition if this method is used outside the context of the WDmodel package. Caveat emptor.

_get_obs_model(teff, logg, av, fwhm, wave, rv=3.1, log=False, pixel_scale=1.0)[source]

Returns the observed model flux given teff, logg, av, rv, fwhm (for Gaussian instrumental broadening) and wavelengths wave

Uses WDmodel.WDmodel.WDmodel._get_model() to get the unreddened model, and reddens it with WDmodel.WDmodel.WDmodel.reddening() and convolves it with a Gaussian kernel using scipy.ndimage.filters.gaussian_filter1d()

Parameters
  • teff (float) – Desired model white dwarf atmosphere temperature (in Kelvin)

  • logg (float) – Desired model white dwarf atmosphere surface gravity (in dex)

  • av (float) – Extinction in the V band, \(A_V\)

  • fwhm (float) – Instrumental FWHM in Angstrom

  • wave (array-like) – Desired wavelengths at which to compute the model atmosphere flux.

  • rv (float, optional) – The reddening law parameter, \(R_V\), the ration of the V band extinction \(A_V\) to the reddening between the B and V bands, \(E(B-V)\). Default is 3.1, appropriate for stellar SEDs in the Milky Way.

  • log (bool, optional) – Return the log10 flux, rather than the flux (what’s actually interpolated)

  • pixel_scale (float, optional) – Jacobian of the transformation between wavelength in Angstrom and pixels. In principle, this should be a vector, but virtually all spectral reduction packages resample the spectrum onto a uniform wavelength scale that is close to the native pixel scale of the spectrograph. Default is 1.

Returns

flux – Interpolated model flux at teff, logg with reddening parametrized by av, rv and broadened by a Gaussian kernel defined by fwhm at wavelengths wave

Return type

array-like

Notes

fwhm and pixel_scale must be > 0

_get_red_model(teff, logg, av, wave, rv=3.1, log=False)[source]

Returns the reddened model flux given teff, logg, av, rv at wavelengths wave

Uses WDmodel.WDmodel.WDmodel._get_model() to get the unreddened model, and reddens it with WDmodel.WDmodel.WDmodel.reddening()

Parameters
  • teff (float) – Desired model white dwarf atmosphere temperature (in Kelvin)

  • logg (float) – Desired model white dwarf atmosphere surface gravity (in dex)

  • av (float) – Extinction in the V band, \(A_V\)

  • wave (array-like) – Desired wavelengths at which to compute the model atmosphere flux.

  • rv (float, optional) – The reddening law parameter, \(R_V\), the ration of the V band extinction \(A_V\) to the reddening between the B and V bands, \(E(B-V)\). Default is 3.1, appropriate for stellar SEDs in the Milky Way.

  • log (bool, optional) – Return the log10 flux, rather than the flux (what’s actually interpolated)

Returns

flux – Interpolated model flux at teff, logg with reddening parametrized by av, rv at wavelengths wave

Return type

array-like

classmethod _wave_test(wave)[source]

Raises an error if wavelengths are not valid

Parameters

wave (array-like) – Array of wavelengths to test for validity

Raises

ValueError – If wavelength array is empty, has negative values, or is not monotonic

extinction(wave, av, rv=3.1)[source]

Return the extinction for av, rv at wavelengths wave

Uses the extinction function corresponding to the rvmodel parametrization set as WDmodel.WDmodel.WDmodel._law to calculate the extinction as a function of wavelength (in Angstroms), \(A_{\lambda}\).

Parameters
  • wave (array-like) – Array of wavelengths in Angstrom at which to compute extinction, sorted in ascending order

  • av (float) – Extinction in the V band, \(A_V\)

  • rv (float, optional) – The reddening law parameter, \(R_V\), the ration of the V band extinction \(A_V\) to the reddening between the B and V bands, \(E(B-V)\). Default is 3.1, appropriate for stellar SEDs in the Milky Way.

Returns

out – Extinction at wavelengths wave for av and rv

Return type

array-like

Notes

av should be >= 0.

extract_spectral_line(w, f, line, df=None)[source]

Extracts slices of multiple arrays corresponding to a hydrogen Balmer line

Convenience function to extract elements of wavelength w, flux f and optionally flux uncertainty df for a hydrogen Balmer line. Wraps WDmodel.WDmodel.WDmodel._extract_spectral_line() adding checking of inputs.

Parameters
  • w (array-like) – Wavelength array from which to extract elements corresponding to hydrogen Balmer line

  • f (array-like) – Flux array from which to extract elements corresponding to hydrogen Balmer line

  • line ({'alpha', 'beta', 'gamma', 'delta', 'zeta', 'eta'}) – Name of hydrogen Balmer line to extract. Properties are pre-defined in WDmodel.WDmodel.WDmodel._lines

  • df (None or array-like, optional) – If array-like, extracted elements of this array are also returned

Returns

  • w (array-like) – elements of input wavelength array for hydrogen Balmer feature line

  • f (array-like) – elements of input flux array for hydrogen Balmer feature line

  • [df] (array-like) – elements of input flux uncertainty array for hydrogen Balmer feature line if optional input df is supplied

Raises

ValueError – If line is not one of the first six of the Balmer series or If wavelengths are invalid of If there’s a difference in shape of any of the arrays

get_model(teff, logg, wave=None, log=False, strict=True)[source]

Returns the model flux given teff and logg at wavelengths wave

Wraps WDmodel.WDmodel.WDmodel._get_model() adding checking of inputs.

Parameters
  • teff (float) – Desired model white dwarf atmosphere temperature (in Kelvin)

  • logg (float) – Desired model white dwarf atmosphere surface gravity (in dex)

  • wave (array-like, optional) – Desired wavelengths at which to compute the model atmosphere flux. If not supplied, the full model wavelength grid is returned.

  • log (bool, optional) – Return the log10 flux, rather than the flux

  • strict (bool, optional) – If strict, teff and logg out of model grid range raise a ValueError, otherwise raise a RuntimeWarning and set teff, logg to the nearest grid value.

Returns

  • wave (array-like) – Valid output wavelengths

  • flux (array-like) – Interpolated model flux at teff, logg and wavelengths wave

Raises

ValueError – If teff or logg are out of range of model grid and strict is True or if there are any invalid wavelengths, or the requested wavelengths to do not overlap with the model grid

Notes

Unlike the corresponding private methods, the public methods implement checking of the inputs and returns the wavelengths in addition to the flux. Internally, we only use the private methods as the inputs only need to be checked once, and their state is not altered anywhere after.

get_obs_model(teff, logg, av, fwhm, rv=3.1, wave=None, log=False, strict=True, pixel_scale=1.0)[source]

Returns the observed model flux given teff, logg, av, rv, fwhm (for Gaussian instrumental broadening) and wavelengths wave

Uses WDmodel.WDmodel.WDmodel.get_red_model() to get the reddened model and convolves it with a Gaussian kernel using scipy.ndimage.filters.gaussian_filter1d()

Parameters
  • teff (float) – Desired model white dwarf atmosphere temperature (in Kelvin)

  • logg (float) – Desired model white dwarf atmosphere surface gravity (in dex)

  • av (float) – Extinction in the V band, \(A_V\)

  • fwhm (float) – Instrumental FWHM in Angstrom

  • rv (float, optional) – The reddening law parameter, \(R_V\), the ration of the V band extinction \(A_V\) to the reddening between the B and V bands, \(E(B-V)\). Default is 3.1, appropriate for stellar SEDs in the Milky Way.

  • wave (array-like, optional) – Desired wavelengths at which to compute the model atmosphere flux. If not supplied, the full model wavelength grid is returned.

  • log (bool, optional) – Return the log10 flux, rather than the flux (what’s actually interpolated)

  • strict (bool, optional) – If strict, teff and logg out of model grid range raise a ValueError, otherwise raise a RuntimeWarning and set teff, logg to the nearest grid value.

  • pixel_scale (float, optional) – Jacobian of the transformation between wavelength in Angstrom and pixels. In principle, this should be a vector, but virtually all spectral reduction packages resample the spectrum onto a uniform wavelength scale that is close to the native pixel scale of the spectrograph. Default is 1.

Returns

  • wave (array-like) – Valid output wavelengths

  • flux (array-like) – Interpolated model flux at teff, logg with reddening parametrized by av, rv broadened by a Gaussian kernel defined by fwhm at wavelengths wave

Notes

fwhm and pixel_scale must be > 0

get_red_model(teff, logg, av, rv=3.1, wave=None, log=False, strict=True)[source]

Returns the reddened model flux given teff, logg, av, rv at wavelengths wave

Uses WDmodel.WDmodel.WDmodel.get_model() to get the unreddened model, and reddens it with WDmodel.WDmodel.WDmodel.reddening()

Parameters
  • teff (float) – Desired model white dwarf atmosphere temperature (in Kelvin)

  • logg (float) – Desired model white dwarf atmosphere surface gravity (in dex)

  • av (float) – Extinction in the V band, \(A_V\)

  • rv (float, optional) – The reddening law parameter, \(R_V\), the ration of the V band extinction \(A_V\) to the reddening between the B and V bands, \(E(B-V)\). Default is 3.1, appropriate for stellar SEDs in the Milky Way.

  • wave (array-like, optional) – Desired wavelengths at which to compute the model atmosphere flux. If not supplied, the full model wavelength grid is returned.

  • log (bool, optional) – Return the log10 flux, rather than the flux (what’s actually interpolated)

  • strict (bool, optional) – If strict, teff and logg out of model grid range raise a ValueError, otherwise raise a RuntimeWarning and set teff, logg to the nearest grid value.

Returns

  • wave (array-like) – Valid output wavelengths

  • flux (array-like) – Interpolated model flux at teff, logg with reddening parametrized by av, rv at wavelengths wave

Raises

ValueError – If av  < 0 or rv not in [1.7, 5.1]

reddening(wave, flux, av, rv=3.1)[source]

Redden a 1-D spectrum with extinction

Uses the extinction function corresponding to the rvmodel parametrization set in WDmodel.WDmodel.WDmodel._WDmodel__init__rvmodel() to calculate the extinction as a function of wavelength (in Angstroms), \(A_{\lambda}\).

Parameters
  • wave (array-like) – Array of wavelengths in Angstrom at which to compute extinction, sorted in ascending order

  • flux (array-like) – Array of fluxes at wave at which to apply extinction

  • av (float) – Extinction in the V band, \(A_V\)

  • rv (float, optional) – The reddening law parameter, \(R_V\), the ration of the V band extinction \(A_V\) to the reddening between the B and V bands, \(E(B-V)\). Default is 3.1, appropriate for stellar SEDs in the Milky Way.

Returns

out – The reddened spectrum

Return type

array-like

Notes

av and flux should be >= 0.