File: spectrum_collection.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 (111 lines) | stat: -rw-r--r-- 3,815 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
================================
Working With SpectrumCollections
================================

A spectrum collection is a way to keep a set of spectra data together and have
the collection behave as if it were a single spectrum object. This means that
it can be used in regular analysis functions to perform operations over entire
sets of data.

Currently, all :class:`~specutils.SpectrumCollection` items must be the same
shape. No assumptions are made about the dispersion solutions, and users are
encouraged to ensure their spectrum collections make sense either by resampling
them beforehand, or being aware that they do not share the same dispersion
solution.

.. code:: python

    >>> import numpy as np
    >>> import astropy.units as u
    >>> from astropy.nddata import StdDevUncertainty
    >>> from specutils import SpectrumCollection
    >>> from specutils.utils.wcs_utils import gwcs_from_array

    >>> flux = u.Quantity(np.random.sample((5, 10)), unit='Jy')
    >>> spectral_axis = u.Quantity(np.arange(50).reshape((5, 10)), unit='AA')
    >>> wcs = np.array([gwcs_from_array(x, x.shape) for x in spectral_axis])
    >>> uncertainty = StdDevUncertainty(np.random.sample((5, 10)), unit='Jy')
    >>> mask = np.ones((5, 10)).astype(bool)
    >>> meta = [{'test': 5, 'info': [1, 2, 3]} for i in range(5)]

    >>> spec_coll = SpectrumCollection(
    ... flux=flux, spectral_axis=spectral_axis, wcs=wcs,
    ... uncertainty=uncertainty, mask=mask, meta=meta)

    >>> spec_coll.shape
    (5,)
    >>> spec_coll.flux.unit
    Unit("Jy")
    >>> spec_coll.spectral_axis.shape
    (5, 10)
    >>> spec_coll.spectral_axis.unit
    Unit("Angstrom")

Collections from 1D spectra
---------------------------

It is also possible to create a :class:`~specutils.SpectrumCollection` from
a list of :class:`~specutils.Spectrum`:

.. code:: python

    >>> import warnings
    >>> import numpy as np
    >>> from astropy import units as u
    >>> from specutils import Spectrum, SpectrumCollection
    >>> spec = Spectrum(spectral_axis=np.linspace(0, 50, 50) * u.AA,
    ...                   flux=np.random.randn(50) * u.Jy,
    ...                   uncertainty=StdDevUncertainty(np.random.sample(50), unit='Jy'))
    >>> with warnings.catch_warnings():  # Ignore warnings
    ...     warnings.simplefilter('ignore')
    ...     spec1 = Spectrum(spectral_axis=np.linspace(20, 60, 50) * u.AA,
    ...                        flux=np.random.randn(50) * u.Jy,
    ...                        uncertainty=StdDevUncertainty(np.random.sample(50), unit='Jy'))
    ...     spec_coll = SpectrumCollection.from_spectra([spec, spec1])

    >>> spec_coll.shape
    (2,)
    >>> spec_coll.flux.unit
    Unit("Jy")
    >>> spec_coll.spectral_axis.shape
    (2, 50)
    >>> spec_coll.spectral_axis.unit
    Unit("Angstrom")

:class:`~specutils.SpectrumCollection` objects can be treated just like
:class:`~specutils.Spectrum` objects; calling a particular attribute on the
object will return an array whose type depends on the type of the attribute in
the :class:`~specutils.Spectrum` object.

.. code:: python

    >>> print(type(spec1.flux))
    <class 'astropy.units.quantity.Quantity'>
    >>> print(type(spec_coll.flux))
    <class 'astropy.units.quantity.Quantity'>

The difference is their shape. The returned array from the
:class:`~specutils.SpectrumCollection`
object will have shape ``(N, M)`` where ``N`` is the number of input spectra
and ``M`` is the length of the output dispersion grid.

.. code:: python

    >>> print(spec1.flux.shape)
    (50,)
    >>> print(spec_coll.flux.shape)
    (2, 50)


Reference/API
-------------

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

    :skip: QTable
    :skip: Spectrum
    :skip: SpectralRegion
    :skip: SpectralAxis