File: matchsep.rst

package info (click to toggle)
astropy 7.2.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 34,156 kB
  • sloc: python: 240,839; ansic: 55,852; lex: 8,621; sh: 3,318; xml: 2,399; makefile: 149
file content (349 lines) | stat: -rw-r--r-- 14,440 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
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
.. _astropy-coordinates-separations-matching:

Separations, Offsets, Catalog Matching, and Related Functionality
*****************************************************************

`astropy.coordinates` contains commonly-used tools for comparing or
matching coordinate objects. Of particular importance are those for
determining separations between coordinates and those for matching a
coordinate (or coordinates) to a catalog. These are mainly implemented
as methods on the coordinate objects.

In the examples below, we will assume that the following imports have already
been executed::

    >>> import astropy.units as u
    >>> from astropy.coordinates import SkyCoord

Separations
===========

The on-sky separation can be computed with
:meth:`~astropy.coordinates.BaseCoordinateFrame.separation`,
which computes the great-circle distance (*not* the small-angle approximation)::

    >>> c1 = SkyCoord('5h23m34.5s', '-69d45m22s', frame='icrs')
    >>> c2 = SkyCoord('0h52m44.8s', '-72d49m43s', frame='fk5')
    >>> sep = c1.separation(c2)
    >>> sep  # doctest: +FLOAT_CMP
    <Angle 20.74611448 deg>

The returned object is an `~astropy.coordinates.Angle` instance, so it
is possible to access the angle in any of several equivalent angular
units::

    >>> sep.radian  # doctest: +FLOAT_CMP
    np.float64(0.36208800460262563)
    >>> sep.hour  # doctest: +FLOAT_CMP
    np.float64(1.3830742984029318)
    >>> sep.arcminute  # doctest: +FLOAT_CMP
    np.float64(1244.7668685626384)
    >>> sep.arcsecond  # doctest: +FLOAT_CMP
    np.float64(74686.0121137583)

Also note that the two input coordinates were not in the same frame —
one is automatically converted to match the other, ensuring that even
though they are in different frames, the separation is determined
consistently.

In addition to the on-sky separation described above,
:meth:`~astropy.coordinates.BaseCoordinateFrame.separation_3d` will
determine the 3D distance between two coordinates that have ``distance``
defined::

    >>> c1 = SkyCoord('5h23m34.5s', '-69d45m22s', distance=70*u.kpc, frame='icrs')
    >>> c2 = SkyCoord('0h52m44.8s', '-72d49m43s', distance=80*u.kpc, frame='icrs')
    >>> sep = c1.separation_3d(c2)
    >>> sep  # doctest: +FLOAT_CMP
    <Distance 28.74398816 kpc>


Offsets
=======

Closely related to angular separations are offsets between coordinates. The key
distinction for offsets is generally the concept of a "from" and "to" coordinate
rather than the single scalar angular offset of a separation.
`~astropy.coordinates` contains conveniences to compute some of the common
offsets encountered in astronomy.

The first piece of such functionality is the
:meth:`~astropy.coordinates.BaseCoordinateFrame.position_angle` method.
This method computes the position angle between one
|SkyCoord| instance and another (passed as the argument) following the
astronomy convention (positive angles East of North)::

    >>> c1 = SkyCoord(1*u.deg, 1*u.deg, frame='icrs')
    >>> c2 = SkyCoord(2*u.deg, 2*u.deg, frame='icrs')
    >>> c1.position_angle(c2).to(u.deg)  # doctest: +FLOAT_CMP
    <Angle 44.97818294 deg>

The combination of :meth:`~astropy.coordinates.BaseCoordinateFrame.separation`
and :meth:`~astropy.coordinates.BaseCoordinateFrame.position_angle` thus give a
set of directional offsets.
To do the inverse operation — determining the new
"destination" coordinate given a separation and position angle — the
:meth:`~astropy.coordinates.SkyCoord.directional_offset_by` method is provided::

    >>> c1 = SkyCoord(1*u.deg, 1*u.deg, frame='icrs')
    >>> position_angle = 45 * u.deg
    >>> separation = 1.414 * u.deg
    >>> c1.directional_offset_by(position_angle, separation)  # doctest: +FLOAT_CMP
    <SkyCoord (ICRS): (ra, dec) in deg
        (2.0004075, 1.99964588)>

This technique is also useful for computing the midpoint (or indeed any point)
between two coordinates in a way that accounts for spherical geometry
(i.e., instead of averaging the RAs/Decs separately)::

    >>> coord1 = SkyCoord(0*u.deg, 0*u.deg, frame='icrs')
    >>> coord2 = SkyCoord(1*u.deg, 1*u.deg, frame='icrs')
    >>> pa = coord1.position_angle(coord2)
    >>> sep = coord1.separation(coord2)
    >>> coord1.directional_offset_by(pa, sep/2)  # doctest: +FLOAT_CMP
    <SkyCoord (ICRS): (ra, dec) in deg
        (0.49996192, 0.50001904)>

There is also a :meth:`~astropy.coordinates.SkyCoord.spherical_offsets_to`
method for computing angular offsets (e.g., small shifts like you might give a
telescope operator to move from a bright star to a fainter target)::

    >>> bright_star = SkyCoord('8h50m59.75s', '+11d39m22.15s', frame='icrs')
    >>> faint_galaxy = SkyCoord('8h50m47.92s', '+11d39m32.74s', frame='icrs')
    >>> dra, ddec = bright_star.spherical_offsets_to(faint_galaxy)
    >>> dra.to(u.arcsec)  # doctest: +FLOAT_CMP
    <Angle -173.78873354 arcsec>
    >>> ddec.to(u.arcsec)  # doctest: +FLOAT_CMP
    <Angle 10.60510342 arcsec>

The conceptual inverse of
:meth:`~astropy.coordinates.SkyCoord.spherical_offsets_to` is also available as
a method on any |SkyCoord| object:
:meth:`~astropy.coordinates.SkyCoord.spherical_offsets_by`, which accepts two
angular offsets (in longitude and latitude) and returns the coordinates at the
offset location::

    >>> target_star = SkyCoord(86.75309*u.deg, -31.5633*u.deg, frame='icrs')
    >>> target_star.spherical_offsets_by(1.3*u.arcmin, -0.7*u.arcmin)  # doctest: +FLOAT_CMP
    <SkyCoord (ICRS): (ra, dec) in deg
        (86.77852168, -31.57496415)>

.. _astropy-skyoffset-frames:

"Sky Offset" Frames
-------------------

To extend the concept of spherical offsets, `~astropy.coordinates` has
a frame class :class:`~astropy.coordinates.SkyOffsetFrame`
which creates distinct frames that are centered on a specific point.
These are known as "sky offset frames," as they are a convenient way to create
a frame centered on an arbitrary position on the sky suitable for computing
positional offsets (e.g., for astrometry)::

    >>> from astropy.coordinates import SkyOffsetFrame, ICRS
    >>> center = ICRS(10*u.deg, 45*u.deg)
    >>> center.transform_to(SkyOffsetFrame(origin=center)) # doctest: +FLOAT_CMP
    <SkyOffsetICRS Coordinate (rotation=0.0 deg, origin=<ICRS Coordinate: (ra, dec) in deg
        (10., 45.)>): (lon, lat) in deg
        (0., 0.)>
    >>> target = ICRS(11*u.deg, 46*u.deg)
    >>> target.transform_to(SkyOffsetFrame(origin=center))  # doctest: +FLOAT_CMP
    <SkyOffsetICRS Coordinate (rotation=0.0 deg, origin=<ICRS Coordinate: (ra, dec) in deg
        (10., 45.)>): (lon, lat) in deg
        (0.69474685, 1.00428706)>


Alternatively, the convenience method
:meth:`~astropy.coordinates.SkyCoord.skyoffset_frame` lets you create a sky
offset frame from an existing |SkyCoord|::

    >>> center = SkyCoord(10*u.deg, 45*u.deg)
    >>> aframe = center.skyoffset_frame()
    >>> target.transform_to(aframe)  # doctest: +FLOAT_CMP
    <SkyOffsetICRS Coordinate (rotation=0.0 deg, origin=<ICRS Coordinate: (ra, dec) in deg
        (10., 45.)>): (lon, lat) in deg
        (0.69474685, 1.00428706)>
    >>> other = SkyCoord(9*u.deg, 44*u.deg, frame='fk5')
    >>> other.transform_to(aframe)  # doctest: +FLOAT_CMP
    <SkyCoord (SkyOffsetICRS: rotation=0.0 deg, origin=<ICRS Coordinate: (ra, dec) in deg
        (10., 45.)>): (lon, lat) in deg
        (-0.71943945, -0.99556216)>

.. note::

    While sky offset frames *appear* to be all the same class, this not the
    case: the sky offset frame for each different type of frame for ``origin`` is
    actually a distinct class. E.g., ``SkyOffsetFrame(origin=ICRS(...))``
    yields an object of class ``SkyOffsetICRS``, *not* ``SkyOffsetFrame``.
    While this is not important for most uses of this class, it is important for
    things like type-checking, because something like
    ``SkyOffsetFrame(origin=ICRS(...)).__class__ is SkyOffsetFrame`` will
    *not* be ``True``, as it would be for most classes.

This same frame is also useful as a tool for defining frames that are relative
to a specific, known object useful for hierarchical physical systems like galaxy
groups. For example, objects around M31 are sometimes shown in a coordinate
frame aligned with standard ICRA RA/Dec, but on M31::

    >>> m31 = SkyCoord(10.6847083*u.deg, 41.26875*u.deg, frame='icrs')
    >>> ngc147 = SkyCoord(8.3005*u.deg, 48.5087389*u.deg, frame='icrs')
    >>> ngc147_inm31 = ngc147.transform_to(m31.skyoffset_frame())
    >>> xi, eta = ngc147_inm31.lon, ngc147_inm31.lat
    >>> xi  # doctest: +FLOAT_CMP
    <Longitude -1.59206948 deg>
    >>> eta  # doctest: +FLOAT_CMP
    <Latitude 7.26183757 deg>

.. note::

    Currently, distance information in the ``origin`` of a
    :class:`~astropy.coordinates.SkyOffsetFrame` is not
    used to compute any part of the transform. The ``origin`` is only used for
    on-sky rotation. This may change in the future, however.


.. _astropy-coordinates-matching:

Matching Catalogs
=================

`~astropy.coordinates` leverages the coordinate framework to make it
possible to find the closest coordinates in a catalog to a desired set
of other coordinates. For example, assuming ``ra1``/``dec1`` and
``ra2``/``dec2`` are NumPy arrays loaded from some file:

.. testsetup::
    >>> ra1 = [5.3517]
    >>> dec1 = [-5.2328]
    >>> distance1 = 1344
    >>> ra2 = [6.459]
    >>> dec2 = [-16.4258]
    >>> distance2 = 8.611

.. doctest-requires:: scipy

    >>> c = SkyCoord(ra=ra1*u.degree, dec=dec1*u.degree)
    >>> catalog = SkyCoord(ra=ra2*u.degree, dec=dec2*u.degree)
    >>> idx, d2d, d3d = c.match_to_catalog_sky(catalog)

The distances returned ``d3d`` are 3-dimensional distances.
Unless both source (``c``) and catalog (``catalog``) coordinates have
associated distances, this quantity assumes that all sources are at a distance
of 1 (dimensionless).

You can also find the nearest 3D matches, different from the on-sky
separation shown above only when the coordinates were initialized with
a ``distance``:

.. doctest-requires:: scipy

    >>> c = SkyCoord(ra=ra1*u.degree, dec=dec1*u.degree, distance=distance1*u.kpc)
    >>> catalog = SkyCoord(ra=ra2*u.degree, dec=dec2*u.degree, distance=distance2*u.kpc)
    >>> idx, d2d, d3d = c.match_to_catalog_3d(catalog)

Now ``idx`` are indices into ``catalog`` that are the closest objects to each
of the coordinates in ``c``, ``d2d`` are the on-sky distances between them, and
``d3d`` are the 3-dimensional distances. Because coordinate objects support
indexing, ``idx`` enables easy access to the matched set of coordinates in
the catalog:

.. doctest-requires:: scipy

    >>> d3d # doctest: +FLOAT_CMP
    <Quantity [1335.55538257] kpc>
    >>> matches = catalog[idx]
    >>> matches.separation_3d(c) # doctest: +FLOAT_CMP
    <Distance [1335.55538257] kpc>
    >>> dra, ddec = c.spherical_offsets_to(matches)

This functionality can also be accessed from the
:func:`~astropy.coordinates.match_coordinates_sky` and
:func:`~astropy.coordinates.match_coordinates_3d` functions. These
will work on either |SkyCoord| objects *or* the lower-level frame classes:

.. doctest-requires:: scipy

    >>> from astropy.coordinates import match_coordinates_sky
    >>> idx, d2d, d3d = match_coordinates_sky(c, catalog)
    >>> idx, d2d, d3d = match_coordinates_sky(c.frame, catalog.frame)

It is possible to impose a separation constraint (e.g., the maximum separation to be
considered a match) by creating a boolean mask with ``d2d`` or ``d3d``. For example:

.. doctest-requires:: scipy

    >>> max_sep = 1.0 * u.arcsec
    >>> idx, d2d, d3d = c.match_to_catalog_3d(catalog)
    >>> sep_constraint = d2d < max_sep
    >>> c_matches = c[sep_constraint]
    >>> catalog_matches = catalog[idx[sep_constraint]]

Now, ``c_matches`` and ``catalog_matches`` are the matched sources in ``c``
and ``catalog``, respectively, which are separated by less than 1 arcsecond.

.. _astropy-searching-coordinates:

Searching around Coordinates
============================

Closely related functionality can be used to search for *all* coordinates within
a certain distance (either 3D distance or on-sky) of another set of coordinates.
The ``search_around_*`` methods (and functions) provide this functionality,
with an interface very similar to ``match_coordinates_*``:

..  doctest-requires:: scipy

    >>> import numpy as np
    >>> idxc, idxcatalog, d2d, d3d = catalog.search_around_sky(c, 1*u.deg)
    >>> np.all(d2d < 1*u.deg)
    np.True_

.. doctest-requires:: scipy

    >>> idxc, idxcatalog, d2d, d3d = catalog.search_around_3d(c, 1*u.kpc)
    >>> np.all(d3d < 1*u.kpc)
    np.True_

The key difference for these methods is that there can be multiple (or no)
matches in ``catalog`` around any locations in ``c``. Hence, indices into both
``c`` and ``catalog`` are returned instead of just indices into ``catalog``.
These can then be indexed back into the two |SkyCoord| objects, or, for that
matter, any array with the same order:

..  doctest-requires:: scipy

    >>> np.all(c[idxc].separation(catalog[idxcatalog]) == d2d)
    np.True_
    >>> np.all(c[idxc].separation_3d(catalog[idxcatalog]) == d3d)
    np.True_
    >>> print(catalog_objectnames[idxcatalog]) #doctest: +SKIP
    ['NGC 1234' 'NGC 4567' ...]

Note, though, that this dual-indexing means that ``search_around_*`` does not
work well if one of the coordinates is a scalar, because the returned index
would not make sense for a scalar::

    >>> scalarc = SkyCoord(ra=1*u.deg, dec=2*u.deg, distance=distance1*u.kpc)
    >>> idxscalarc, idxcatalog, d2d, d3d = catalog.search_around_sky(scalarc, 1*u.deg) # doctest: +SKIP
    ValueError: One of the inputs to search_around_sky is a scalar.

As a result (and because the ``search_around_*`` algorithm is inefficient in
the scalar case), the best approach for this scenario is to instead
use the ``separation*`` methods:

..  doctest-requires:: scipy

    >>> d2d = scalarc.separation(catalog)
    >>> catalogmsk = d2d < 1*u.deg
    >>> d3d = scalarc.separation_3d(catalog)
    >>> catalog3dmsk = d3d < 1*u.kpc

The resulting ``catalogmsk`` or ``catalog3dmsk`` variables are boolean arrays
rather than arrays of indices, but in practice they usually can be used in
the same way as ``idxcatalog`` from the above examples. If you definitely do
need indices instead of boolean masks, you can do:

..  doctest-requires:: scipy

    >>> idxcatalog = np.where(catalogmsk)[0]
    >>> idxcatalog3d = np.where(catalog3dmsk)[0]