File: getting_started.rst

package info (click to toggle)
astropy-regions 0.11-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 2,096 kB
  • sloc: python: 7,703; makefile: 117
file content (238 lines) | stat: -rw-r--r-- 6,762 bytes parent folder | download | duplicates (3)
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
.. include:: references.txt

.. _getting_started:

Getting Started
===============

Introduction
------------

The Regions package provides classes to represent:

* Regions defined using pixel coordinates (e.g.,
  `~regions.CirclePixelRegion`)
* Regions defined using celestial coordinates, but still in an Euclidean
  geometry (e.g., `~regions.CircleSkyRegion`)

To transform between sky and pixel regions, a `world coordinate system
<https://docs.astropy.org/en/stable/wcs/wcsapi.html>`_ object (e.g.,
`astropy.wcs.WCS`) is needed.

Regions also provides a unified interface for reading, writing,
parsing, and serializing regions data in different formats, including
the `DS9 Region Format <http://ds9.si.edu/doc/ref/region.html>`_,
`CRTF (CASA Region Text Format)
<https://casadocs.readthedocs.io/en/stable/notebooks/image_analysis.html
#Region-File-Format>`_, and `FITS Region Binary Table
<https://fits.gsfc.nasa.gov/registry/region.html>`_ format.


.. _getting_started-coord:

Coordinates
-----------

Pixel Coordinates
~~~~~~~~~~~~~~~~~

:class:`~regions.PixCoord` objects are used to represent pixel
coordinates. Pixel coordinates are defined with a scalar or an array of
``x`` and ``y`` Cartesian coordinates:

.. code-block:: python

    >>> from regions import PixCoord
    >>> pixcoord = PixCoord(x=42, y=43)
    >>> pixcoord
    PixCoord(x=42, y=43)
    >>> pixcoord.x
    42
    >>> pixcoord.y
    43
    >>> pixcoord.xy
    (42, 43)

`~regions.PixCoord` objects can also represent arrays of pixel
coordinates. These work in the same way as single-value coordinates, but
they store multiple coordinates in a single object. Let's create a 1D
array of pixel coordinates:

.. code-block:: python

    >>> pixcoord = PixCoord(x=[0, 1], y=[2, 3])
    >>> pixcoord
    PixCoord(x=[0 1], y=[2 3])
    >>> pixcoord.x
    array([0, 1])
    >>> pixcoord.y
    array([2, 3])
    >>> pixcoord.xy
    (array([0, 1]), array([2, 3]))

Let's now create a 2D array of pixel coordinates:

.. code-block:: python

    >>> pixcoord = PixCoord(x=[[1, 2, 3], [4, 5, 6]],
    ...                     y=[[11, 12, 13], [14, 15, 16]])
    >>> pixcoord
    PixCoord(x=[[1 2 3]
     [4 5 6]], y=[[11 12 13]
     [14 15 16]])


Sky Coordinates
~~~~~~~~~~~~~~~

:class:`~astropy.coordinates.SkyCoord` objects are used to represent
sky coordinates. `~astropy.coordinates.SkyCoord` is a very powerful
class that provides a flexible interface for celestial coordinate
representation, manipulation, and transformation between systems. See
the extensive :ref:`astropy-coordinates` documentation for more details.

Let's create a single sky coordinate:

.. code-block:: python

    >>> from astropy.coordinates import SkyCoord
    >>> skycoord = SkyCoord(42, 43, unit='deg', frame='galactic')
    >>> skycoord
    <SkyCoord (Galactic): (l, b) in deg
        (42., 43.)>

Sky coordinates also support array coordinates. These work in the same
way as single-value coordinates, but they store multiple coordinates in
a single object:

.. code-block:: python

    >>> skycoord = SkyCoord(ra=[10, 11, 12], dec=[41, 42, 43], unit='deg')
    >>> skycoord
    <SkyCoord (ICRS): (ra, dec) in deg
    [(10., 41.), (11., 42.), (12., 43.)]>

To represent angles both on the sky and in an image,
`~astropy.coordinates.Angle` objects or `~astropy.units.Quantity`
objects with angular units can be used.


Pixel/Sky Coordinate Transformations
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

To transform between pixel and sky coordinates, a `world coordinate system
<https://docs.astropy.org/en/stable/wcs/wcsapi.html>`_ object (e.g.,
`astropy.wcs.WCS`) is needed.

Let's start by creating a WCS object:

.. code-block:: python

    >>> from regions import make_example_dataset
    >>> dataset = make_example_dataset(data='simulated')
    >>> wcs = dataset.wcs

Now let's use this WCS to convert between sky and pixel coordinates:

.. code-block:: python

    >>> skycoord = SkyCoord(42, 43, unit='deg', frame='galactic')
    >>> pixcoord = PixCoord.from_sky(skycoord=skycoord, wcs=wcs)
    >>> pixcoord   # doctest: +FLOAT_CMP
    PixCoord(x=146.2575703393558, y=131.5998051082584)
    >>> pixcoord.to_sky(wcs=wcs)
    <SkyCoord (Galactic): (l, b) in deg
        (42., 43.)>


Pixel Regions
-------------

Pixel regions are regions that are defined using pixel coordinates and
sizes in pixels. The regions package provides a set of "pixel-based"
regions classes, for example, `~regions.CirclePixelRegion`:

.. code-block:: python

    >>> from regions import PixCoord, CirclePixelRegion
    >>> center = PixCoord(x=42, y=43)
    >>> radius = 4.2
    >>> region = CirclePixelRegion(center, radius)

You can print the region to get some information about its properties:

.. code-block:: python

    >>> print(region)
    Region: CirclePixelRegion
    center: PixCoord(x=42, y=43)
    radius: 4.2

You can access its properties via attributes:

.. code-block:: python

   >>> region.center
   PixCoord(x=42, y=43)
   >>> region.radius
   4.2

See the :ref:`shapes` documentation for the complete list of pixel-based
regions and to learn more about :class:`~regions.Region` objects and
their capabilities.


Sky Regions
-----------

Sky regions are regions that are defined using celestial coordinates.
Please note they are **not** defined as regions on the celestial sphere,
but rather are meant to represent shapes on an image. They simply use
sky coordinates instead of pixel coordinates to define their position.
The remaining shape parameters are converted to pixels using the pixel
scale of the image.

Let's create a sky region:

.. code-block:: python

    >>> from astropy.coordinates import Angle, SkyCoord
    >>> from regions import CircleSkyRegion
    >>> center = SkyCoord(42, 43, unit='deg')
    >>> radius = Angle(3, 'deg')
    >>> region = CircleSkyRegion(center, radius)

Alternatively, one can define the radius using a
`~astropy.units.Quantity` object with angular units:

.. code-block:: python

    >>> import astropy.units as u
    >>> from regions import CircleSkyRegion
    >>> center = SkyCoord(42, 43, unit='deg')
    >>> radius = 3.0 * u.deg
    >>> region = CircleSkyRegion(center, radius)

You can print the region to get some information about its properties:

.. code-block:: python

    >>> print(region)
    Region: CircleSkyRegion
    center: <SkyCoord (ICRS): (ra, dec) in deg
        (42., 43.)>
    radius: 3.0 deg

You can access its properties via attributes:

.. code-block:: python

   >>> region.center
    <SkyCoord (ICRS): (ra, dec) in deg
        (42., 43.)>
   >>> region.radius
   <Quantity 3. deg>

See the :ref:`shapes` documentation for the complete list of pixel-based
regions and to learn more about :class:`~regions.Region` objects and
their capabilities.