File: spectral_regions.rst

package info (click to toggle)
specutils 2.2.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 12,460 kB
  • sloc: python: 13,497; makefile: 111
file content (355 lines) | stat: -rw-r--r-- 13,396 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
================
Spectral Regions
================

A spectral region may be defined and may encompass one, or more,
sub-regions. They are defined independently of a `~specutils.Spectrum`
object in the sense that spectral regions like "near the Halpha line rest
wavelength" have meaning independent of the details of a particular spectrum.

Spectral regions can be defined either as a single region by passing
two `~astropy.units.Quantity`'s or by passing a list of 2-tuples. Note that
the units of these quantites can be any valid spectral unit *or* ``u.pixel``
(which indicates to use indexing directly).

.. code-block:: python

    >>> from astropy import units as u
    >>> from specutils.spectra import SpectralRegion

    >>> sr = SpectralRegion(0.45*u.um, 0.6*u.um)
    >>> sr_two = SpectralRegion([(0.45*u.um, 0.6*u.um), (0.8*u.um, 0.9*u.um)])

`~specutils.SpectralRegion` can be combined by using the '+' operator:

.. code-block:: python

    >>> from astropy import units as u
    >>> from specutils.spectra import SpectralRegion
    >>> sr = SpectralRegion(0.45*u.um, 0.6*u.um) + SpectralRegion(0.8*u.um, 0.9*u.um)

Regions can also be added in place:

.. code-block:: python

    >>> from astropy import units as u
    >>> from specutils.spectra import SpectralRegion

    >>> sr1 = SpectralRegion(0.45*u.um, 0.6*u.um)
    >>> sr2 = SpectralRegion(0.8*u.um, 0.9*u.um)
    >>> sr1 += sr2

Regions can be sliced by indexing by an integer or by a range:

.. code-block:: python

    >>> from astropy import units as u
    >>> from specutils.spectra import SpectralRegion

    >>> sr = SpectralRegion(0.15*u.um, 0.2*u.um) + SpectralRegion(0.3*u.um, 0.4*u.um) +\
    ...      SpectralRegion(0.45*u.um, 0.6*u.um) + SpectralRegion(0.8*u.um, 0.9*u.um) +\
    ...      SpectralRegion(1.0*u.um, 1.2*u.um) + SpectralRegion(1.3*u.um, 1.5*u.um)

    >>> # Get one spectral region (returns a SpectralRegion instance)
    >>> sone = sr1[0]

    >>> # Slice spectral region.
    >>> subsr = sr[3:5]
    >>> # SpectralRegion: 0.8 um - 0.9 um, 1.0 um - 1.2 um

The lower and upper bounds on a region are accessed by calling lower
or upper. The lower bound of a `~specutils.SpectralRegion` is the
minimum of the lower bounds of each sub-region and the upper bound is the
maximum of the upper bounds:

.. code-block:: python

    >>> from astropy import units as u
    >>> from specutils.spectra import SpectralRegion

    >>> sr = (SpectralRegion(0.15*u.um, 0.2*u.um) + SpectralRegion(0.3*u.um, 0.4*u.um) +
    ...       SpectralRegion(0.45*u.um, 0.6*u.um) + SpectralRegion(0.8*u.um, 0.9*u.um) +
    ...       SpectralRegion(1.0*u.um, 1.2*u.um) + SpectralRegion(1.3*u.um, 1.5*u.um))

    >>> # Bounds on the spectral region (most minimum and maximum bound)
    >>> sr.bounds
    (<Quantity 0.15 um>, <Quantity 1.5 um>)

    >>> # Lower bound on the spectral region (most minimum)
    >>> sr.lower
    <Quantity 0.15 um>

    >>> sr.upper
    <Quantity 1.5 um>

    >>> # Lower bound on one element of the spectral region.
    >>> sr[3].lower
    <Quantity 0.8 um>

One can also delete a sub-region:

.. code-block:: python

    >>> from astropy import units as u
    >>> from specutils.spectra import SpectralRegion

    >>> sr = (SpectralRegion(0.15*u.um, 0.2*u.um) + SpectralRegion(0.3*u.um, 0.4*u.um) +
    ...       SpectralRegion(0.45*u.um, 0.6*u.um) + SpectralRegion(0.8*u.um, 0.9*u.um) +
    ...       SpectralRegion(1.0*u.um, 1.2*u.um) + SpectralRegion(1.3*u.um, 1.5*u.um))

    >>> del sr[1]
    >>> sr
    Spectral Region, 5 sub-regions:
    (0.15 um, 0.2 um)   (0.45 um, 0.6 um)   (0.8 um, 0.9 um)
    (1.0 um, 1.2 um)    (1.3 um, 1.5 um)

There is also the ability to iterate:

.. code-block:: python

    >>> from astropy import units as u
    >>> from specutils.spectra import SpectralRegion

    >>> sr = (SpectralRegion(0.15*u.um, 0.2*u.um) + SpectralRegion(0.3*u.um, 0.4*u.um) +
    ...       SpectralRegion(0.45*u.um, 0.6*u.um) + SpectralRegion(0.8*u.um, 0.9*u.um) +
    ...       SpectralRegion(1.0*u.um, 1.2*u.um) + SpectralRegion(1.3*u.um, 1.5*u.um))

    >>> for s in sr:
    ...     print(s.lower)
    0.15 um
    0.3 um
    0.45 um
    0.8 um
    1.0 um
    1.3 um


And, lastly, there is the ability to invert a `~specutils.SpectralRegion` given a
lower and upper bound. For example, if a set of ranges are defined each defining a range
around lines, then calling invert will return a `~specutils.SpectralRegion` that
defines the baseline/noise regions:

.. code-block:: python

    >>> from astropy import units as u
    >>> from specutils.spectra import SpectralRegion

    >>> sr = (SpectralRegion(0.15*u.um, 0.2*u.um) + SpectralRegion(0.3*u.um, 0.4*u.um) +
    ...       SpectralRegion(0.45*u.um, 0.6*u.um) + SpectralRegion(0.8*u.um, 0.9*u.um) +
    ...       SpectralRegion(1.0*u.um, 1.2*u.um) + SpectralRegion(1.3*u.um, 1.5*u.um))

    >>> sr_inverted = sr.invert(0.05*u.um, 3*u.um)
    >>> sr_inverted
    Spectral Region, 7 sub-regions:
    (0.05 um, 0.15 um)   (0.2 um, 0.3 um)     (0.4 um, 0.45 um)
    (0.6 um, 0.8 um)     (0.9 um, 1.0 um)     (1.2 um, 1.3 um)
    (1.5 um, 3.0 um)

Region Extraction
-----------------

Given a `~specutils.SpectralRegion`, one can extract a sub-spectrum
from a `~specutils.Spectrum` object. If the `~specutils.SpectralRegion`
has multiple sub-regions then by default a list of `~specutils.Spectrum` 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.Spectrum` object instead.

An example of a single sub-region `~specutils.SpectralRegion`:

.. code-block:: python

    >>> from astropy import units as u
    >>> import numpy as np
    >>> from specutils import Spectrum, SpectralRegion
    >>> from specutils.manipulation import extract_region

    >>> region = SpectralRegion(8*u.nm, 22*u.nm)
    >>> spectrum = Spectrum(spectral_axis=np.arange(1, 50) * u.nm, flux=np.random.sample(49)*u.Jy)
    >>> sub_spectrum = extract_region(spectrum, region)
    >>> sub_spectrum.spectral_axis
    <SpectralAxis [ 8.,  9., 10., 11., 12., 13., 14., 15., 16., 17., 18., 19., 20., 21.,
       22.] nm>

Extraction also correctly interprets different kinds of spectral region units
as would be expected:

.. code-block:: python

    >>> from astropy import units as u
    >>> import numpy as np
    >>> from specutils import Spectrum, SpectralRegion
    >>> from specutils.manipulation import extract_region

    >>> spectrum = Spectrum(spectral_axis=np.arange(1, 50) * u.nm, flux=np.random.sample(49)*u.Jy)
    >>> region_angstroms = SpectralRegion(80*u.AA, 220*u.AA)
    >>> sub_spectrum = extract_region(spectrum, region_angstroms)
    >>> sub_spectrum.spectral_axis
    <SpectralAxis [ 8.,  9., 10., 11., 12., 13., 14., 15., 16., 17., 18., 19., 20., 21.,
       22.] nm>
    >>> region_pixels = SpectralRegion(7.5*u.pixel, 21.5*u.pixel)
    >>> sub_spectrum = extract_region(spectrum, region_pixels)
    >>> sub_spectrum.spectral_axis
    <SpectralAxis [ 8.,  9., 10., 11., 12., 13., 14., 15., 16., 17., 18., 19., 20., 21.,
       22.] nm>

An example of a multiple sub-region `~specutils.SpectralRegion`:

.. code-block:: python

    >>> from astropy import units as u
    >>> import numpy as np
    >>> from specutils import Spectrum, SpectralRegion
    >>> from specutils.manipulation import extract_region

    >>> region = SpectralRegion([(8*u.nm, 22*u.nm), (34*u.nm, 40*u.nm)])
    >>> spectrum = Spectrum(spectral_axis=np.arange(1, 50) * u.nm, flux=np.random.sample(49)*u.Jy)
    >>> sub_spectra = extract_region(spectrum, region)
    >>> sub_spectra[0].spectral_axis
    <SpectralAxis [ 8.,  9., 10., 11., 12., 13., 14., 15., 16., 17., 18., 19., 20., 21.,
       22.] nm>
    >>> sub_spectra[1].spectral_axis
    <SpectralAxis [34., 35., 36., 37., 38., 39., 40.] nm>

Multiple sub-regions can also be returned as a single concatenated spectrum:

.. code-block:: python

    >>> sub_spectrum = extract_region(spectrum, region, return_single_spectrum=True)
    >>> sub_spectrum.spectral_axis
    <SpectralAxis [ 8.,  9., 10., 11., 12., 13., 14., 15., 16., 17., 18., 19., 20., 21.,
        22., 34., 35., 36., 37., 38., 39., 40.] nm>

The bounding region that includes all data, including the ones that lie
in between disjointed spectral regions, can be extracted with
`specutils.manipulation.extract_bounding_spectral_region`:

.. code-block:: python

    >>> from astropy import units as u
    >>> import numpy as np
    >>> from specutils import Spectrum, SpectralRegion
    >>> from specutils.manipulation import extract_bounding_spectral_region

    >>> spectrum = Spectrum(spectral_axis=np.arange(1, 50) * u.nm,
    ...                       flux=np.random.default_rng(12345).random(49)*u.Jy)
    >>> region = SpectralRegion([(8*u.nm, 12*u.nm), (24*u.nm, 30*u.nm)])
    >>> sub_spectrum = extract_bounding_spectral_region(spectrum, region)
    >>> sub_spectrum.spectral_axis
    <SpectralAxis [ 8.,  9., 10., 11., 12., 13., 14., 15., 16., 17., 18., 19., 20., 21.,
        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
axis range to be extracted.
This function's purpose is to facilitate migration of ``spectral_cube`` functionality
into ``specutils``:

.. code-block:: python

    >>> from astropy import units as u
    >>> import numpy as np
    >>> from specutils import Spectrum, SpectralRegion
    >>> from specutils.manipulation import spectral_slab

    >>> spectrum = Spectrum(spectral_axis=np.arange(1, 50) * u.nm,
    ...                       flux=np.random.default_rng(12345).random(49)*u.Jy)
    >>> sub_spectrum = spectral_slab(spectrum, 8*u.nm, 20*u.nm)
    >>> sub_spectrum.spectral_axis
    <SpectralAxis [ 8.,  9., 10., 11., 12., 13., 14., 15., 16., 17., 18., 19., 20.] nm>


Line List Compatibility
-----------------------

`~specutils.SpectralRegion` objects can also be created from the `~astropy.table.QTable` object returned from the line
finding functions:

.. code-block:: python

    >>> from astropy import units as u
    >>> import numpy as np
    >>> from specutils import Spectrum, SpectralRegion
    >>> from astropy.modeling.models import Gaussian1D
    >>> from specutils.fitting import find_lines_derivative

    >>> g1 = Gaussian1D(1, 4.6, 0.2)
    >>> g2 = Gaussian1D(2.5, 5.5, 0.1)
    >>> g3 = Gaussian1D(-1.7, 8.2, 0.1)

    >>> x = np.linspace(0, 10, 200)
    >>> y = g1(x) + g2(x) + g3(x)
    >>> spectrum = Spectrum(flux=y * u.Jy, spectral_axis=x * u.um)

    >>> lines = find_lines_derivative(spectrum, flux_threshold=0.01)
    >>> spec_reg = SpectralRegion.from_line_list(lines)
    >>> spec_reg
    Spectral Region, 3 sub-regions:
      (4.072864321608041 um, 5.072864321608041 um)
      (4.977386934673367 um, 5.977386934673367 um)
      (7.690954773869347 um, 8.690954773869347 um)

This can be fed into the ``exclude_regions`` argument of the `~specutils.fitting.fit_generic_continuum` or
`~specutils.fitting.fit_continuum` functions to avoid fitting regions that contain line features. Note that,
by default, this uses pythonic slicing, i.e., spectral values greater than or equal to the lower bound and
less than the upper bound of the region will be excluded from the fit. For convenience in some cases, the
``exclude_region_upper_bounds`` keyword can be set to ``True`` to exlude spectral values less than *or equal*
to the upper bound instead.

Conversely, users can also invert the spectral region

.. code-block:: python

    >>> inv_spec_reg = spec_reg.invert(spectrum.spectral_axis[0], spectrum.spectral_axis[-1])
    >>> inv_spec_reg
    Spectral Region, 3 sub-regions:
      (0.0 um, 4.072864321608041 um)
      (5.977386934673367 um, 7.690954773869347 um)
      (8.690954773869347 um, 10.0 um)

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")

.. testcleanup::

    >>> import os
    >>> os.remove("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
-------------

.. automodapi:: specutils
    :no-main-docstr:
    :no-heading:
    :no-inheritance-diagram:

    :skip: QTable
    :skip: Spectrum
    :skip: SpectrumCollection
    :skip: SpectralAxis