Skip to content

Commit

Permalink
Add read, write, as_table, and from_qtable to SpectralRegion (#1133)
Browse files Browse the repository at this point in the history
* Add read, write, as_table, and from_qtable to SpectralRegion

* Clarify docstrings

* Add skips to automodapi

* Add section in docs for new functionality

* Fix doc build

* Add test

* Codestyle

* Use plain old property instead of cached_property
  • Loading branch information
rosteen committed Apr 30, 2024
1 parent bf2355b commit 4e8412a
Show file tree
Hide file tree
Showing 5 changed files with 107 additions and 5 deletions.
36 changes: 31 additions & 5 deletions docs/spectral_regions.rst
Original file line number Diff line number Diff line change
Expand Up @@ -150,9 +150,9 @@ Region Extraction

Given a `~specutils.SpectralRegion`, one can extract a sub-spectrum
from a `~specutils.Spectrum1D` object. If the `~specutils.SpectralRegion`
has multiple sub-regions then by default a list of `~specutils.Spectrum1D` objects
will be returned. If the ``return_single_spectrum`` argument is set to ``True``,
the resulting spectra will be concatenated together into a single
has multiple sub-regions then by default a list of `~specutils.Spectrum1D` objects
will be returned. If the ``return_single_spectrum`` argument is set to ``True``,
the resulting spectra will be concatenated together into a single
`~specutils.Spectrum1D` object instead.

An example of a single sub-region `~specutils.SpectralRegion`:
Expand Down Expand Up @@ -240,8 +240,8 @@ in between disjointed spectral regions, can be extracted with
22., 23., 24., 25., 26., 27., 28., 29., 30.] nm>
`~specutils.manipulation.spectral_slab` is basically an alternate entry point for
`~specutils.manipulation.extract_region`. Notice the slightly different way to input the spectral
`~specutils.manipulation.spectral_slab` is basically an alternate entry point for
`~specutils.manipulation.extract_region`. Notice the slightly different way to input the spectral
axis range to be extracted.
This function's purpose is to facilitate migration of ``spectral_cube`` functionality
into ``specutils``:
Expand Down Expand Up @@ -306,6 +306,31 @@ also invert the spectral region
and use that result as the ``exclude_regions`` argument in the `~specutils.fitting.fit_lines` function in order to avoid
attempting to fit any of the continuum region.

Reading and Writing
-------------------

`~specutils.SpectralRegion` objects can be written to ``ecsv`` files, which uses the `~astropy.table.QTable` write machinery:

.. code-block:: python
>>> spec_reg.write("spectral_region.ecsv")
This overwrites existing files by default if a duplicate filename is input. The resulting files can be read back in
to create a new `~specutils.SpectralRegion` using the ``read`` class method:

.. code-block:: python
>>> spec_reg = SpectralRegion.read("spectral_region.ecsv")
The `~astropy.table.QTable` created to write out the `~specutils.SpectralRegion` to file can also be accessed
directly with the ``as_table`` method, and a `~specutils.SpectralRegion` can be created directly from a `~astropy.table.QTable`
with the appropriate columns (minimally ``lower_bound`` and ``upper_bound``) using the ``from_qtable`` class method.

.. code-block:: python
>>> spec_reg_table = spec_reg.as_table()
>>> spec_reg_2 = SpectralRegion.from_qtable(spec_reg_table)
Reference/API
-------------

Expand All @@ -315,6 +340,7 @@ Reference/API
:no-inheritance-diagram:

:skip: test
:skip: QTable
:skip: Spectrum1D
:skip: SpectrumCollection
:skip: SpectralAxis
1 change: 1 addition & 0 deletions docs/spectrum1d.rst
Original file line number Diff line number Diff line change
Expand Up @@ -386,5 +386,6 @@ Reference/API
:headings: -~

:skip: test
:skip: QTable
:skip: SpectrumCollection
:skip: SpectralRegion
1 change: 1 addition & 0 deletions docs/spectrum_collection.rst
Original file line number Diff line number Diff line change
Expand Up @@ -106,6 +106,7 @@ Reference/API
:no-inheritance-diagram:

:skip: test
:skip: QTable
:skip: Spectrum1D
:skip: SpectralRegion
:skip: SpectralAxis
63 changes: 63 additions & 0 deletions specutils/spectra/spectral_region.py
Original file line number Diff line number Diff line change
@@ -1,6 +1,7 @@
import itertools
import sys

from astropy.table import QTable
import astropy.units as u


Expand Down Expand Up @@ -106,6 +107,42 @@ def from_line_list(cls, table, width=1):
return cls([(x - width * 0.5, x + width * 0.5)
for x in table['line_center']])

@classmethod
def from_qtable(cls, table):
"""
Generate a ``SpectralRegion`` instance from an `~astropy.table.QTable`
object has has ``lower_bound`` and ``upper_bound`` columns
Parameters
----------
table : `~astropy.table.QTable`
An `~astropy.table.QTable` object with ``lower_bound`` and ``upper_bound`` columns.
Returns
-------
`~specutils.SpectralRegion`
The spectral region based on the table of bounds.
"""
subregions = []
for row in table:
subregions.append([row["lower_bound"], row["upper_bound"]])

return cls(subregions)

@classmethod
def read(cls, filename):
"""
Create a ``SpectralRegion`` from an ecsv file output by the
``SpectralRegion.write`` method.
Parameters
----------
filename : str
The name of the ecsv file on disk to be read in as a ``SpectralRegion``.
"""
table = QTable.read(filename)
return cls.from_qtable(table)

def _info(self):
"""
Pretty print the sub-regions.
Expand Down Expand Up @@ -337,3 +374,29 @@ def invert(self, lower_bound, upper_bound):
#

return SpectralRegion([(x, y) for x, y in zip(newlist[0::2], newlist[1::2])])

@property
def _table(self):
# Return the information defining the spectral region as an astropy QTable
lower_bounds = []
upper_bounds = []
for subregion in self.subregions:
lower_bounds.append(subregion[0])
upper_bounds.append(subregion[1])

return QTable([lower_bounds, upper_bounds], names=("lower_bound", "upper_bound"))

def as_table(self):
"""Returns an `~astropy.table.QTable` with the upper and lower bound
of each subregion in the ``SpectralRegion``."""
return self._table

def write(self, filename="spectral_region.ecsv", overwrite=True):
"""
Write the SpectralRegion to an ecsv file using `~astropy.table.QTable`.
Overwrites by default.
"""
if not filename.endswith("ecsv"):
raise ValueError("SpectralRegions can only be written out to ecsv files.")

self._table.write(filename, overwrite=overwrite)
11 changes: 11 additions & 0 deletions specutils/tests/test_regions.py
Original file line number Diff line number Diff line change
Expand Up @@ -189,3 +189,14 @@ def test_from_list_list():

for i, reg in enumerate(expected):
assert_quantity_allclose(reg, (spec_reg[i].lower, spec_reg[i].upper))


def test_read_write():
sr = SpectralRegion([(0.45*u.um, 0.6*u.um), (0.8*u.um, 0.9*u.um)])
sr.write("test_sr.ecsv")
sr2 = SpectralRegion.read("test_sr.ecsv")
assert list(sr2.as_table().columns) == ["lower_bound", "upper_bound"]

sr3 = SpectralRegion.from_qtable(sr2.as_table())
assert sr3.subregions[0] == sr.subregions[0]
assert sr3.subregions[1] == sr.subregions[1]

0 comments on commit 4e8412a

Please sign in to comment.