File: ccddata.rst

package info (click to toggle)
astropy 5.2.1-2
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 41,972 kB
  • sloc: python: 219,331; ansic: 147,297; javascript: 13,556; lex: 8,496; sh: 3,319; xml: 1,622; makefile: 185
file content (219 lines) | stat: -rw-r--r-- 8,710 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
.. _ccddata:


CCDData Class
=============

Getting Started
---------------

Getting Data In
+++++++++++++++

Creating a `~astropy.nddata.CCDData` object from any array-like data using
`astropy.nddata` is convenient:

    >>> import numpy as np
    >>> from astropy.nddata import CCDData
    >>> from astropy.utils.data import get_pkg_data_filename
    >>> ccd = CCDData(np.arange(10), unit="adu")

Note that behind the scenes, this creates references to (not copies of) your
data when possible, so modifying the data in ``ccd`` will modify the
underlying data.

You are **required** to provide a unit for your data. The most frequently used
units for these objects are likely to be ``adu``, ``photon``, and ``electron``,
which can be set either by providing the string name of the unit (as in the
example above) or from unit objects:

    >>> from astropy import units as u
    >>> ccd_photon = CCDData([1, 2, 3], unit=u.photon)
    >>> ccd_electron = CCDData([1, 2, 3], unit="electron")

If you prefer *not* to use the unit functionality, then use the special unit
``u.dimensionless_unscaled`` when you create your `~astropy.nddata.CCDData`
images:

    >>> ccd_unitless = CCDData(np.zeros((10, 10)),
    ...                        unit=u.dimensionless_unscaled)

A `~astropy.nddata.CCDData` object can also be initialized from a FITS filename
or URL:

    >>> ccd = CCDData.read('my_file.fits', unit="adu")  # doctest: +SKIP
    >>> ccd = CCDData.read(get_pkg_data_filename('tutorials/FITS-images/HorseHead.fits'), unit="adu", cache=True)  # doctest: +REMOTE_DATA +IGNORE_WARNINGS

If there is a unit in the FITS file (in the ``BUNIT`` keyword), that will be
used, but explicitly providing a unit in ``read`` will override any unit in the
FITS file.

There is no restriction at all on what the unit can be — any unit in
`astropy.units` or another that you create yourself will work.

In addition, the user can specify the extension in a FITS file to use:

    >>> ccd = CCDData.read('my_file.fits', hdu=1, unit="adu")  # doctest: +SKIP

If ``hdu`` is not specified, it will assume the data is in the primary
extension. If there is no data in the primary extension, the first extension
with image data will be used.

Metadata
++++++++

When initializing from a FITS file, the ``header`` property is initialized using
the header of the FITS file. Metadata is optional, and can be provided by any
dictionary or dict-like object:

    >>> ccd_simple = CCDData(np.arange(10), unit="adu")
    >>> my_meta = {'observer': 'Edwin Hubble', 'exposure': 30.0}
    >>> ccd_simple.header = my_meta  # or use ccd_simple.meta = my_meta

Whether the metadata is case-sensitive or not depends on how it is
initialized. A FITS header, for example, is not case-sensitive, but a Python
dictionary is.

Getting Data Out
++++++++++++++++

A `~astropy.nddata.CCDData` object behaves like a ``numpy`` array (masked if the
`~astropy.nddata.CCDData` mask is set) in expressions, and the underlying
data (ignoring any mask) is accessed through the ``data`` attribute:

    >>> ccd_masked = CCDData([1, 2, 3], unit="adu", mask=[0, 0, 1])
    >>> 2 * np.ones(3) * ccd_masked   # one return value will be masked
    masked_array(data=[2.0, 4.0, --],
                 mask=[False, False,  True],
           fill_value=1e+20)
    >>> 2 * np.ones(3) * ccd_masked.data   # ignores the mask  # doctest: +FLOAT_CMP
    array([2., 4., 6.])

You can force conversion to a ``numpy`` array with:

    >>> np.asarray(ccd_masked)
    array([1, 2, 3])
    >>> np.ma.array(ccd_masked.data, mask=ccd_masked.mask)
    masked_array(data=[1, 2, --],
                 mask=[False, False,  True],
           fill_value=999999)

A method for converting a `~astropy.nddata.CCDData` object to a FITS HDU list
is also available. It converts the metadata to a FITS header:

    >>> hdulist = ccd_masked.to_hdu()

You can also write directly to a FITS file:

    >>> ccd_masked.write('my_image.fits')

Masks and Flags
+++++++++++++++

Although it is not required when a `~astropy.nddata.CCDData` image is created,
you can also specify a mask and/or flags.

A mask is a boolean array the same size as the data in which a value of
``True`` indicates that a particular pixel should be masked (*i.e.*, not be
included in arithmetic operations or aggregation).

Flags are one or more additional arrays (of any type) whose shape matches the
shape of the data. One particularly useful type of flag is a bit planes; for
more details about bit planes and the functions ``astropy`` provides for
converting them to binary masks, see :ref:`bitmask_details`. For more details
on setting flags, see `~astropy.nddata.NDData`.

WCS
+++

The ``wcs`` attribute of a `~astropy.nddata.CCDData` object can be set two ways.

+ If the `~astropy.nddata.CCDData` object is created from a FITS file that has
  WCS keywords in the header, the ``wcs`` attribute is set to a
  `~astropy.wcs.WCS` object using the information in the FITS header.

+ The WCS can also be provided when the `~astropy.nddata.CCDData` object is
  constructed with the ``wcs`` argument.

Either way, the ``wcs`` attribute is kept up to date if the
`~astropy.nddata.CCDData` image is trimmed.

PSF
+++

The ``psf`` attributes of a `~astropy.nddata.CCDData` object can be set two ways.

+ If the FITS file has an image HDU extension matching the appropriate name (defaulted to ``"PSFIMAGE"``), the ``psf`` attribute is loaded from that image HDU.

+ The PSF can also be provided when the `~astropy.nddata.CCDData` object is
  constructed with the ``psf`` argument.

The ``psf`` attribute should be a normalized image representing the PSF at the center of the `~astropy.nddata.CCDData`, sized appropriately for the data; users are responsible for managing and interpreting it in context.
For more on normalizing a PSF image, see :ref:`astropy:kernel_normalization`.

The ``psf`` attribute is set to `None` in the output of an arithmetic operation, no matter the inputs. A warning message is emitted if either of the input images contain a non-`None` PSF; users are responsible for determining the appropriate thing to do in that context.

Uncertainty
-----------

You can set the uncertainty directly, either by creating a
`~astropy.nddata.StdDevUncertainty` object first:

    >>> rng = np.random.default_rng()
    >>> data = rng.normal(size=(10, 10), loc=1.0, scale=0.1)
    >>> ccd = CCDData(data, unit="electron")
    >>> from astropy.nddata.nduncertainty import StdDevUncertainty
    >>> uncertainty = 0.1 * ccd.data  # can be any array whose shape matches the data
    >>> my_uncertainty = StdDevUncertainty(uncertainty)
    >>> ccd.uncertainty = my_uncertainty

Or by providing a `~numpy.ndarray` with the same shape as the data:

    >>> ccd.uncertainty = 0.1 * ccd.data  # doctest: +ELLIPSIS
    INFO: array provided for uncertainty; assuming it is a StdDevUncertainty. [...]

In this case, the uncertainty is assumed to be
`~astropy.nddata.StdDevUncertainty`.

Two other uncertainty classes are available for which error propagation is
also supported: `~astropy.nddata.VarianceUncertainty` and
`~astropy.nddata.InverseVariance`. Using one of these three uncertainties is
required to enable error propagation in `~astropy.nddata.CCDData`.

If you want access to the underlying uncertainty, use its ``.array`` attribute:

    >>> ccd.uncertainty.array  # doctest: +ELLIPSIS
    array(...)

Arithmetic with Images
----------------------

Methods are provided to perform arithmetic operations with a
`~astropy.nddata.CCDData` image and a number, an ``astropy``
`~astropy.units.Quantity` (a number with units), or another
`~astropy.nddata.CCDData` image.

Using these methods propagates errors correctly (if the errors are
uncorrelated), takes care of any necessary unit conversions, and applies masks
appropriately. Note that the metadata of the result is *not* set if the
operation is between two `~astropy.nddata.CCDData` objects.

    >>> result = ccd.multiply(0.2 * u.adu)
    >>> uncertainty_ratio = result.uncertainty.array[0, 0]/ccd.uncertainty.array[0, 0]
    >>> round(uncertainty_ratio, 5)   # doctest: +FLOAT_CMP
    0.2
    >>> result.unit
    Unit("adu electron")

.. note::
    The affiliated package `ccdproc <https://ccdproc.readthedocs.io>`_ provides
    functions for many common data reduction operations. Those functions try to
    construct a sensible header for the result and provide a mechanism for
    logging the action of the function in the header.


The arithmetic operators ``*``, ``/``, ``+``, and ``-`` are *not* overridden.

.. note::
   If two images have different WCS values, the ``wcs`` on the first
   `~astropy.nddata.CCDData` object will be used for the resultant object.