pypeit.specutils.pypeit_loaders module
Data parsers for use with specutils.
Usage
To read a PypeIt spec1d file using specutils:
from pypeit.specutils import SpectrumList spec = SpectrumList.read(spec1d_file)
where
spec1d_file
is the relative or absolute path to a PypeIt spec1d file. You must always use aSpectrumList
to read a spec1d file, even if the file only contains one spectrum. The spec1d loader provides PypeIt-specific options that enable you to specify the type of extraction used and whether or not to use the flux-calibrated spectrum; seepypeit.specutils.pypeit_loaders.pypeit_spec1d_loader()
. By default, optimal extraction takes precedence over boxcar extraction, and flux-calibrated data take precedence over uncalibrated counts.To read a PypeIt
pypeit.onespec.OneSpec
file:from pypeit.specutils import Spectrum1D spec = Spectrum1D.read(onespec_file)
where
onespec_file
is the relative or absolute path to a PypeItpypeit.onespec.OneSpec
file. For these files, you can use eitherSpectrum1D
orSpectrumList
to read the file, but (by definition) the result of usingSpectrumList
will just be a list of oneSpectrum1D
object. Thepypeit.onespec.OneSpec
loader provides a PypeIt-specific option that enables you to select the uniform grid wavelength vector, instead of the contribution-weighted wavelengths; seepypeit.specutils.pypeit_loaders.pypeit_onespec_loader()
.
Note
Importing Spectrum1D
and SpectrumList
are shown as coming from the
pypeit.specutils
module, but the objects themselves are identical to the
specutils objects. The reason they are imported from within PypeIt is
that, under the hood, the import also “registers” the PypeIt-specific
loaders with the relevant specutils module. This circumvents the need to
place any pypeit specific code in a ~/.specutils
directory (as
recommended here) and
keeps the import statement to one line. That is,
from pypeit.specutils import Spectrum1D
is really just shorthand for
from specutils import Spectrum1D
from pypeit.specutils import pypeit_loaders
Examples
In addition to the pypeit_show_1dspec GUI, these specutils loaders allow you to interact with your spectra using jdaviz. To do so, use the following lines in a jupyter notebook (you currently must do this from within a notebook):
from pypeit.specutils import SpectrumList
from jdaviz import Specviz
file = 'my_spec1d_file.fits'
spec = SpectrumList.read(file)
specviz = Specviz()
specviz.load_spectrum(spec)
specviz.show()
Version History
2022-02-04: Initial Version (tbowers)
2022-09-16: Correct an import error and add module docstring (tbowers)
2023-03-09: Moved into the main pypeit repo and refactored (KBW)
2023-06-23: Added sensible error message for incorrect spec1d loading (tbowers)
- pypeit.specutils.pypeit_loaders._enforce_monotonic_wavelengths(wave, flux, ivar, strict=True)[source]
Force the spectrum to have a monotonically increasing wavelength vector.
- Parameters:
wave (numpy.ndarray) – Spectrum wavelengths
flux (numpy.ndarray) – Spectrum flux
ivar (numpy.ndarray) – Spectrum inverse variance. Can be None or the standard deviation. The only operation on this and the
flux
vector is to downselect the monotonically increasing values.strict (bool, optional) – Check that the wavelength vector is monotonically increasing. If not, raise an error (as would be done by the specutils.SpectrumList class). If False, wavelengths that are not monotonically increasing are masked in the construction of the returned specutils.SpectrumList object.
- Returns:
wave (numpy.ndarray) – Edited wavelength vector. This may be an unchanged reference to the original vector.
flux (numpy.ndarray) – Edited flux vector. This may be an unchanged reference to the original vector.
ivar (numpy.ndarray) – Edited inverse variance vector. This may be an unchanged reference to the original vector.
- pypeit.specutils.pypeit_loaders._identify_pypeit(*args, **kwargs)[source]
Check if a file is a PypeIt output file, in the most general sense.
This currently only checks if
VERSPYP
is in the primary header.
- pypeit.specutils.pypeit_loaders.identify_pypeit_onespec(origin, *args, **kwargs)[source]
Check if a file is a PypeIt file that follows the
pypeit.onespec.OneSpec
datamodel.In addition to checking if it is a PypeIt file (see
_identify_pypeit()
), this also checks that the datamodel classes are correct.
- pypeit.specutils.pypeit_loaders.identify_pypeit_spec1d(origin, *args, **kwargs)[source]
Check if a file is a PypeIt spec1d file.
In addition to checking if it is a PypeIt file (see
_identify_pypeit()
), this also checks that the datamodel classes are correct.
- pypeit.specutils.pypeit_loaders.pypeit_onespec_loader(filename, grid=False, strict=True, chk_version=True, **kwargs)[source]
Load a spectrum from a PypeIt OneSpec file into a Spectrum1D object.
- Parameters:
filename (str) – The path to the FITS file
grid (bool, optional) – Use the uniform grid wavelengths, instead of the contribution-weighted center.
strict (bool, optional) – Check that the wavelength vector is monotonically increasing. If not, raise an error (as would be done by the specutils.Spectrum1D class). If False, wavelengths that are not monotonically increasing are masked in the construction of the returned specutils.Spectrum1D object.
chk_version (
bool
, optional) – When reading in existing files written by PypeIt, perform strict version checking to ensure a valid file. If False, the code will try to keep going, but this may lead to faults and quiet failures. User beware!kwargs (dict, optional) – Ignored! Used to catch spurious arguments passed to the base class that are ignored by this function.
- Returns:
spec – Spectrum in the PypeIt OneSpec file
- Return type:
- pypeit.specutils.pypeit_loaders.pypeit_spec1d_loader(filename, extract=None, fluxed=True, strict=True, chk_version=True, **kwargs)[source]
Load spectra from a PypeIt spec1d file into a SpectrumList.
- Parameters:
filename (str) – The path to the FITS file
extract (str, optional) – The extraction used to produce the spectrum. Must be either None,
'BOX'
(for a boxcar extraction), or'OPT'
for optimal extraction. If None, the optimal extraction will be returned, if it exists, otherwise the boxcar extraction will be returned.fluxed (bool, optional) – If True, return the flux-calibrated spectrum, if it exists. If the flux calibration hasn’t been performed or
fluxed=False
, the spectrum is returned in counts.strict (bool, optional) – Check that the wavelength vector is monotonically increasing. If not, raise an error (as would be done by the specutils.SpectrumList class). If False, wavelengths that are not monotonically increasing are masked in the construction of the returned specutils.SpectrumList object.
chk_version (
bool
, optional) – When reading in existing files written by PypeIt, perform strict version checking to ensure a valid file. If False, the code will try to keep going, but this may lead to faults and quiet failures. User beware!kwargs (dict, optional) – Ignored! Used to catch spurious arguments passed to the base class that are ignored by this function.
- Returns:
spec – Contains all spectra in the PypeIt spec1d file
- Return type:
- pypeit.specutils.pypeit_loaders.pypeit_spec1d_loader_nolist(filename, extract=None, fluxed=True, **kwargs)[source]
Sensible error message if a user tries to load spectra from a PypeIt spec1d file into a Spectrum1D.
This is not allowed because spec1d files may contain mutliple spectra. This function accepts all arguments as the SpectrumList version, but only outputs a PypeIt Error with a sensible message.
This avoids receiving unhelpful error messages such as:
OSError: Could not identify column containing the wavelength, frequency or energy
- Parameters:
filename (str) – The path to the FITS file
extract (str, optional) – The extraction used to produce the spectrum. Must be either None,
'BOX'
(for a boxcar extraction), or'OPT'
for optimal extraction. If None, the optimal extraction will be returned, if it exists, otherwise the boxcar extraction will be returned.fluxed (bool, optional) – If True, return the flux-calibrated spectrum, if it exists. If the flux calibration hasn’t been performed or
fluxed=False
, the spectrum is returned in counts.