File: ndslicing.rst

package info (click to toggle)
astropy 7.0.1-3
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 35,328 kB
  • sloc: python: 233,437; ansic: 55,264; javascript: 17,680; lex: 8,621; sh: 3,317; xml: 2,287; makefile: 191
file content (162 lines) | stat: -rw-r--r-- 5,058 bytes parent folder | download | duplicates (2)
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
.. _nddata_slicing:

Slicing and Indexing NDData
***************************

Introduction
============

This page only deals with peculiarities that apply to
`~astropy.nddata.NDData`-like classes. For a tutorial about slicing/indexing see the
`python documentation <https://docs.python.org/3/tutorial/introduction.html#lists>`_
and `numpy documentation <https://numpy.org/doc/stable/reference/routines.indexing.html>`_.

.. warning::
    `~astropy.nddata.NDData` and `~astropy.nddata.NDDataRef` enforce almost no
    restrictions on the properties, so it might happen that some **valid but
    unusual** combinations of properties always result in an IndexError or
    incorrect results. In this case, see :ref:`nddata_subclassing` on how to
    customize slicing for a particular property.


Slicing NDDataRef
=================

Unlike `~astropy.nddata.NDData` the class `~astropy.nddata.NDDataRef`
implements slicing or indexing. The result will be wrapped inside the same
class as the sliced object.

Getting one element::

    >>> import numpy as np
    >>> from astropy.nddata import NDDataRef

    >>> data = np.array([1, 2, 3, 4])
    >>> ndd = NDDataRef(data)
    >>> ndd[1]
    NDDataRef(2)

Getting a sliced portion of the original::

    >>> ndd[1:3]  # Get element 1 (inclusive) to 3 (exclusive)
    NDDataRef([2, 3])

This will return a reference (and as such **not a copy**) of the original
properties, so changing a slice will affect the original::

    >>> ndd_sliced = ndd[1:3]
    >>> ndd_sliced.data[0] = 5
    >>> ndd_sliced
    NDDataRef([5, 3])
    >>> ndd
    NDDataRef([1, 5, 3, 4])

But only the one element that was indexed is affected (for example,
``ndd_sliced = ndd[1]``). The element is a scalar and changes will not
propagate to the original.

Slicing NDDataRef Including Attributes
======================================

In the case that a ``mask``, or ``uncertainty`` is present, this
attribute will be sliced too::

    >>> from astropy.nddata import StdDevUncertainty
    >>> data = np.array([1, 2, 3, 4])
    >>> mask = data > 2
    >>> uncertainty = StdDevUncertainty(np.sqrt(data))
    >>> ndd = NDDataRef(data, mask=mask, uncertainty=uncertainty)
    >>> ndd_sliced = ndd[1:3]

    >>> ndd_sliced.data
    array([2, 3])

    >>> ndd_sliced.mask
    array([False,  True]...)

    >>> ndd_sliced.uncertainty  # doctest: +FLOAT_CMP
    StdDevUncertainty([1.41421356, 1.73205081])

``unit`` and ``meta``, however, will be unaffected.

If any of the attributes are set but do not implement slicing, an info will be
printed and the property will be kept as is::

    >>> data = np.array([1, 2, 3, 4])
    >>> mask = False
    >>> uncertainty = StdDevUncertainty(0)
    >>> ndd = NDDataRef(data, mask=mask, uncertainty=uncertainty)
    >>> ndd_sliced = ndd[1:3]
    INFO: uncertainty cannot be sliced. [astropy.nddata.mixins.ndslicing]
    INFO: mask cannot be sliced. [astropy.nddata.mixins.ndslicing]

    >>> ndd_sliced.mask
    False


Slicing NDData with World Coordinates
-------------------------------------

If ``wcs`` is set, it must be either implement
`~astropy.wcs.wcsapi.BaseLowLevelWCS` or `~astropy.wcs.wcsapi.BaseHighLevelWCS`.
This means that only integer or range slices without a step are supported. So
slices like ``[::10]`` or array or boolean based slices will not work.

If you want to slice an ``NDData`` object called ``ndd`` without the WCS you can remove the
WCS from the ``NDData`` object by running:

    >>> ndd.wcs = None


Removing Masked Data
--------------------

.. warning::
    If ``wcs`` is set this will **NOT** be possible. But you can work around
    this by setting the wcs attribute to `None` with ``ndd.wcs = None`` before slicing.

By convention, the ``mask`` attribute indicates if a point is valid or invalid.
So we are able to get all valid data points by slicing with the mask.

Examples
^^^^^^^^

..
  EXAMPLE START
  Removing Masked Data in NDDataRef

To get all of the valid data points by slicing with the mask::

    >>> data = np.array([[1,2,3],[4,5,6],[7,8,9]])
    >>> mask = np.array([[0,1,0],[1,1,1],[0,0,1]], dtype=bool)
    >>> uncertainty = StdDevUncertainty(np.sqrt(data))
    >>> ndd = NDDataRef(data, mask=mask, uncertainty=uncertainty)
    >>> # don't forget that ~ or you'll get the invalid points
    >>> ndd_sliced = ndd[~ndd.mask]
    >>> ndd_sliced
    NDDataRef([1, 3, 7, 8])

    >>> ndd_sliced.mask
    array([False, False, False, False]...)

    >>> ndd_sliced.uncertainty  # doctest: +FLOAT_CMP
    StdDevUncertainty([1.        , 1.73205081, 2.64575131, 2.82842712])

Or all invalid points::

    >>> ndd_sliced = ndd[ndd.mask] # without the ~ now!
    >>> ndd_sliced
    NDDataRef([—, —, —, —, —])

    >>> ndd_sliced.mask
    array([ True,  True,  True,  True,  True]...)

    >>> ndd_sliced.uncertainty  # doctest: +FLOAT_CMP
    StdDevUncertainty([1.41421356, 2.        , 2.23606798, 2.44948974, 3.        ])

.. note::
    The result of this kind of indexing (boolean indexing) will always be
    one-dimensional!

..
  EXAMPLE END