Working with Spectrum1Ds¶
As described in more detail in Overview of How Specutils Represents Spectra, the core data class in
specutils for a single spectrum is Spectrum1D. This object
can represent either one or many spectra, all with the same spectral_axis.
This section describes some of the basic features of this class.
Basic Spectrum Creation¶
The simplest (and most powerful) way to create a Spectrum1D is to
create it explicitly from arrays or Quantity objects:
>>> import numpy as np
>>> import astropy.units as u
>>> import matplotlib.pyplot as plt
>>> from specutils import Spectrum1D
>>> flux = np.random.randn(200)*u.Jy
>>> wavelength = np.arange(5100, 5300)*u.AA
>>> spec1d = Spectrum1D(spectral_axis=wavelength, flux=flux)
>>> ax = plt.subplots()[1]
>>> ax.plot(spec1d.spectral_axis, spec1d.flux)
>>> ax.set_xlabel("Dispersion")
>>> ax.set_ylabel("Flux")
(Source code, png, hires.png, pdf)
Note
Note that the spectral_axis can also be provided as a SpectralAxis object,
and in fact will internally convert the spectral_axis to SpectralAxis if it
is provided as an array or Quantity.
Reading from a File¶
specutils takes advantage of the Astropy IO machinery and allows loading and
writing to files. The example below shows loading a FITS file. While specutils
has some basic data loaders, for more complicated or custom files, users are
encouraged to create their own loader.
>>> from specutils import Spectrum1D
>>> spec1d = Spectrum1D.read("/path/to/file.fits")
Most of the built-in specutils default loaders can also read an existing
astropy.io.fits.HDUList object or an open file object (as resulting
from e.g. streaming a file from the internet). Note that in these cases, a
format string corresponding to an existing loader must be supplied because
these objects lack enough contextual information to automatically identify
a loader.
>>> from specutils import Spectrum1D
>>> import urllib
>>> specs = urllib.request.urlopen('https://data.sdss.org/sas/dr14/sdss/spectro/redux/26/spectra/0751/spec-0751-52251-0160.fits')
>>> Spectrum1D.read(specs, format="SDSS-III/IV spec")
<Spectrum1D(flux=<Quantity [30.596626,...]...>
Note that the same spectrum could be more conveniently downloaded via astroquery, if the user has that package installed:
>>> from astroquery.sdss import SDSS
>>> specs = SDSS.get_spectra(plate=751, mjd=52251, fiberID=160)
>>> Spectrum1D.read(specs[0], format="SDSS-III/IV spec")
<Spectrum1D(flux=<Quantity [30.596626,...]...>
List of Loaders¶
The Spectrum1D class has built-in support for various input and output formats.
A full list of the supported formats is shown in the table below.
Format |
Read |
Write |
Auto-identify |
|---|---|---|---|
6dFGS-split |
Yes |
No |
Yes |
6dFGS-tabular |
Yes |
No |
Yes |
APOGEE apStar |
Yes |
No |
Yes |
APOGEE apVisit |
Yes |
No |
Yes |
APOGEE aspcapStar |
Yes |
No |
Yes |
ASCII |
Yes |
No |
Yes |
ECSV |
Yes |
No |
Yes |
HST/COS |
Yes |
No |
Yes |
HST/STIS |
Yes |
No |
Yes |
IPAC |
Yes |
No |
Yes |
JWST s2d |
Yes |
No |
Yes |
JWST s3d |
Yes |
No |
Yes |
JWST x1d |
Yes |
No |
Yes |
MUSCLES SED |
Yes |
No |
Yes |
SDSS-I/II spSpec |
Yes |
No |
Yes |
SDSS-III/IV spec |
Yes |
No |
Yes |
Subaru-pfsObject |
Yes |
No |
Yes |
iraf |
Yes |
No |
Yes |
tabular-fits |
Yes |
Yes |
Yes |
wcs1d-fits |
Yes |
No |
Yes |
Including Uncertainties¶
The Spectrum1D class supports uncertainties, and
arithmetic operations performed with Spectrum1D
objects will propagate uncertainties.
Uncertainties are a special subclass of NDData, and their
propagation rules are implemented at the class level. Therefore, users must
specify the uncertainty type at creation time
>>> from specutils import Spectrum1D
>>> from astropy.nddata import StdDevUncertainty
>>> spec = Spectrum1D(spectral_axis=np.arange(5000, 5010)*u.AA, flux=np.random.sample(10)*u.Jy, uncertainty=StdDevUncertainty(np.random.sample(10) * 0.1))
Warning
Not defining an uncertainty class will result in an
UnknownUncertainty object which will not
propagate uncertainties in arithmetic operations.
Defining WCS¶
Specutils always maintains a WCS object whether it is passed explicitly by the
user, or is created dynamically by specutils itself. In the latter case, the
user need not be aware that the WCS object is being used, and can interact
with the Spectrum1D object as if it were only a simple
data container.
Currently, specutils understands two WCS formats: FITS WCS and GWCS. When a user does not explicitly supply a WCS object, specutils will fallback on an internal GWCS object it will create.
Note
To create a custom adapter for a different WCS class (i.e. aside from FITSWCS or GWCS), please see the documentation on WCS Adapter classes.
Providing a FITS-style WCS¶
>>> from specutils.spectra import Spectrum1D
>>> import astropy.wcs as fitswcs
>>> import astropy.units as u
>>> import numpy as np
>>> my_wcs = fitswcs.WCS(header={'CDELT1': 1, 'CRVAL1': 6562.8, 'CUNIT1': 'Angstrom', 'CTYPE1': 'WAVE', 'RESTFRQ': 1400000000, 'CRPIX1': 25})
>>> spec = Spectrum1D(flux=[5,6,7] * u.Jy, wcs=my_wcs)
>>> spec.spectral_axis
<Quantity [ 6538.8, 6539.8, 6540.8] Angstrom>
>>> spec.wcs.pixel_to_world(np.arange(3))
array([6.5388e-07, 6.5398e-07, 6.5408e-07])
Multi-dimensional Data Sets¶
Spectrum1D also supports the multidimensional case where you
have, say, an (n_spectra, n_pix)
shaped data set where each n_spectra element provides a different flux
data array and so flux and uncertainty may be multidimensional as
long as the last dimension matches the shape of spectral_axis This is meant
to allow fast operations on collections of spectra that share the same
spectral_axis. While it may seem to conflict with the “1D” in the class
name, this name scheme is meant to communicate the presence of a single
common spectral axis.
Note
The case where each flux data array is related to a different spectral
axis is encapsulated in the SpectrumCollection
object described in the related docs.
>>> from specutils import Spectrum1D
>>> spec = Spectrum1D(spectral_axis=np.arange(5000, 5010)*u.AA, flux=np.random.sample((5, 10))*u.Jy)
>>> spec_slice = spec[0]
>>> spec_slice.spectral_axis
<Quantity [5000., 5001., 5002., 5003., 5004., 5005., 5006., 5007., 5008., 5009.] Angstrom>
>>> spec_slice.flux
<Quantity [0.72722821, 0.32147784, 0.70256482, 0.04445197, 0.03390352,
0.50835299, 0.87581725, 0.50270413, 0.08556376, 0.53713355] Jy>
While the above example only shows two dimensions, this concept generalizes to
any number of dimensions for Spectrum1D, as long as the spectral
axis is always the last.
Reference/API¶
Classes¶
Configuration parameters for specutils. |
|
Coordinate object representing spectral values corresponding to a specific spectrum. |
|
|
Spectrum container for 1D spectral data. |
A list that is used to hold a list of Spectrum1D objects |