Skip to content

Commit

Permalink
Remove SkyCoord separation methods
Browse files Browse the repository at this point in the history
By making very small changes to `BaseCoordinateFrame.separation()` and
`BaseCoordinateFrame.separation_3d()` methods it was possible to make
them handle `SkyCoord` instances just as well as `BaseCoordinateFrame`
instances. This allowed the `SkyCoord` methods to be removed because a
`SkyCoord` instance exposes the methods of it underlying frame.
  • Loading branch information
eerovaher committed Nov 29, 2023
1 parent 6b36e53 commit 6c6a50b
Show file tree
Hide file tree
Showing 5 changed files with 44 additions and 130 deletions.
19 changes: 14 additions & 5 deletions astropy/coordinates/baseframe.py
Original file line number Diff line number Diff line change
Expand Up @@ -1777,9 +1777,12 @@ def separation(self, other):
that ``self.separation(other)`` and ``other.separation(self)`` may
not give the same answer in this case.
For more on how to use this (and related) functionality, see the
examples in :doc:`astropy:/coordinates/matchsep`.
Parameters
----------
other : `~astropy.coordinates.BaseCoordinateFrame`
other : `~astropy.coordinates.BaseCoordinateFrame` or `~astropy.coordinates.SkyCoord`
The coordinate to get the separation to.
Returns
Expand All @@ -1798,8 +1801,11 @@ def separation(self, other):
from .angles import Angle, angular_separation

self_unit_sph = self.represent_as(r.UnitSphericalRepresentation)
other_transformed = other.transform_to(self)
other_unit_sph = other_transformed.represent_as(r.UnitSphericalRepresentation)
other_unit_sph = (

Check warning on line 1804 in astropy/coordinates/baseframe.py

View check run for this annotation

Codecov / codecov/patch

astropy/coordinates/baseframe.py#L1804

Added line #L1804 was not covered by tests
getattr(other, "frame", other)
.transform_to(self)
.represent_as(r.UnitSphericalRepresentation)
)

# Get the separation as a Quantity, convert to Angle in degrees
sep = angular_separation(
Expand All @@ -1812,9 +1818,12 @@ def separation_3d(self, other):
Computes three dimensional separation between this coordinate
and another.
For more on how to use this (and related) functionality, see the
examples in :doc:`astropy:/coordinates/matchsep`.
Parameters
----------
other : `~astropy.coordinates.BaseCoordinateFrame`
other : `~astropy.coordinates.BaseCoordinateFrame` or `~astropy.coordinates.SkyCoord`
The coordinate system to get the distance to.
Returns
Expand All @@ -1835,7 +1844,7 @@ def separation_3d(self, other):
)

# do this first just in case the conversion somehow creates a distance
other_in_self_system = other.transform_to(self)
other_in_self_system = getattr(other, "frame", other).transform_to(self)

Check warning on line 1847 in astropy/coordinates/baseframe.py

View check run for this annotation

Codecov / codecov/patch

astropy/coordinates/baseframe.py#L1847

Added line #L1847 was not covered by tests

if issubclass(other_in_self_system.__class__, r.UnitSphericalRepresentation):
raise ValueError(
Expand Down
121 changes: 13 additions & 108 deletions astropy/coordinates/sky_coordinate.py
Original file line number Diff line number Diff line change
@@ -1,7 +1,10 @@
from __future__ import annotations

import copy
import operator
import re
import warnings
from typing import TYPE_CHECKING

import erfa
import numpy as np
Expand Down Expand Up @@ -31,6 +34,9 @@
_parse_coordinate_data,
)

if TYPE_CHECKING:
from typing import Callable

__all__ = ["SkyCoord", "SkyCoordInfo"]


Expand Down Expand Up @@ -288,6 +294,10 @@ class name that allow for creating a |SkyCoord| object and transforming
# info property.
info = SkyCoordInfo()

# Methods implemented by the underlying frame
separation: Callable[[BaseCoordinateFrame | SkyCoord], Angle]
separation_3d: Callable[[BaseCoordinateFrame | SkyCoord], Distance]

def __init__(self, *args, copy=True, **kwargs):
# these are frame attributes set on this SkyCoord but *not* a part of
# the frame object this SkyCoord contains
Expand Down Expand Up @@ -1133,111 +1143,6 @@ def is_equivalent_frame(self, other):
)

# High-level convenience methods
def separation(self, other):
"""
Computes on-sky separation between this coordinate and another.
.. note::
If the ``other`` coordinate object is in a different frame, it is
first transformed to the frame of this object. This can lead to
unintuitive behavior if not accounted for. Particularly of note is
that ``self.separation(other)`` and ``other.separation(self)`` may
not give the same answer in this case.
For more on how to use this (and related) functionality, see the
examples in :doc:`astropy:/coordinates/matchsep`.
Parameters
----------
other : `~astropy.coordinates.SkyCoord` or `~astropy.coordinates.BaseCoordinateFrame`
The coordinate to get the separation to.
Returns
-------
sep : `~astropy.coordinates.Angle`
The on-sky separation between this and the ``other`` coordinate.
Notes
-----
The separation is calculated using the Vincenty formula, which
is stable at all locations, including poles and antipodes [1]_.
.. [1] https://en.wikipedia.org/wiki/Great-circle_distance
"""
from .angles import Angle, angular_separation

if not self.is_equivalent_frame(other):
try:
kwargs = (
{"merge_attributes": False} if isinstance(other, SkyCoord) else {}
)
other = other.transform_to(self, **kwargs)
except TypeError:
raise TypeError(
"Can only get separation to another SkyCoord "
"or a coordinate frame with data"
)

lon1 = self.spherical.lon
lat1 = self.spherical.lat
lon2 = other.spherical.lon
lat2 = other.spherical.lat

# Get the separation as a Quantity, convert to Angle in degrees
sep = angular_separation(lon1, lat1, lon2, lat2)
return Angle(sep, unit=u.degree)

def separation_3d(self, other):
"""
Computes three dimensional separation between this coordinate
and another.
For more on how to use this (and related) functionality, see the
examples in :doc:`astropy:/coordinates/matchsep`.
Parameters
----------
other : `~astropy.coordinates.SkyCoord` or `~astropy.coordinates.BaseCoordinateFrame`
The coordinate to get the separation to.
Returns
-------
sep : `~astropy.coordinates.Distance`
The real-space distance between these two coordinates.
Raises
------
ValueError
If this or the other coordinate do not have distances.
"""
if not self.is_equivalent_frame(other):
try:
kwargs = (
{"merge_attributes": False} if isinstance(other, SkyCoord) else {}
)
other = other.transform_to(self, **kwargs)
except TypeError:
raise TypeError(
"Can only get separation to another SkyCoord "
"or a coordinate frame with data"
)

if issubclass(self.data.__class__, UnitSphericalRepresentation):
raise ValueError(
"This object does not have a distance; cannot compute 3d separation."
)
if issubclass(other.data.__class__, UnitSphericalRepresentation):
raise ValueError(
"The other object does not have a distance; "
"cannot compute 3d separation."
)

c1 = self.cartesian.without_differentials()
c2 = other.cartesian.without_differentials()
return Distance((c1 - c2).norm())

def spherical_offsets_to(self, tocoord):
r"""
Computes angular offsets to go *from* this coordinate *to* another.
Expand Down Expand Up @@ -1274,7 +1179,7 @@ def spherical_offsets_to(self, tocoord):
See Also
--------
separation :
:meth:`~astropy.coordinates.BaseCoordinateFrame.separation` :
for the *total* angular offset (not broken out into components).
position_angle :
for the direction of the offset.
Expand Down Expand Up @@ -1522,7 +1427,7 @@ def search_around_sky(self, searcharoundcoords, seplimit):
This is intended for use on `~astropy.coordinates.SkyCoord` objects
with coordinate arrays, rather than a scalar coordinate. For a scalar
coordinate, it is better to use
`~astropy.coordinates.SkyCoord.separation`.
:meth:`~astropy.coordinates.BaseCoordinateFrame.separation`.
For more on how to use this (and related) functionality, see the
examples in :doc:`astropy:/coordinates/matchsep`.
Expand Down Expand Up @@ -1582,7 +1487,7 @@ def search_around_3d(self, searcharoundcoords, distlimit):
This is intended for use on `~astropy.coordinates.SkyCoord` objects
with coordinate arrays, rather than a scalar coordinate. For a scalar
coordinate, it is better to use
`~astropy.coordinates.SkyCoord.separation_3d`.
:meth:`~astropy.coordinates.BaseCoordinateFrame.separation_3d`.
For more on how to use this (and related) functionality, see the
examples in :doc:`astropy:/coordinates/matchsep`.
Expand Down
16 changes: 9 additions & 7 deletions docs/coordinates/common_errors.rst
Original file line number Diff line number Diff line change
Expand Up @@ -9,8 +9,9 @@ Object Separation
-----------------

When calculating the separation between objects, it is important to bear in mind that
:meth:`astropy.coordinates.SkyCoord.separation` gives a different answer depending
upon the order in which is used. For example::
:meth:`~astropy.coordinates.BaseCoordinateFrame.separation` gives a different
answer depending upon the order in which is used.
For example::

>>> import numpy as np
>>> from astropy import units as u
Expand All @@ -26,11 +27,12 @@ upon the order in which is used. For example::

Why do these give such different answers?

The reason is that :meth:`astropy.coordinates.SkyCoord.separation` gives the separation as measured
in the frame of the |SkyCoord| object. So ``star.separation(moon)`` gives the angular separation
in the ICRS frame. This is the separation as it would appear from the Solar System Barycenter. For a
geocentric observer, ``moon.separation(star)`` gives the correct answer, since ``moon`` is in a
geocentric frame.
The reason is that :meth:`~astropy.coordinates.BaseCoordinateFrame.separation`
gives the separation as measured in the frame of the |SkyCoord| object.
So ``star.separation(moon)`` gives the angular separation in the ICRS frame.
This is the separation as it would appear from the Solar System Barycenter.
For a geocentric observer, ``moon.separation(star)`` gives the correct answer,
since ``moon`` is in a geocentric frame.

AltAz calculations for Earth-based objects
------------------------------------------
Expand Down
12 changes: 5 additions & 7 deletions docs/coordinates/matchsep.rst
Original file line number Diff line number Diff line change
Expand Up @@ -18,9 +18,8 @@ been executed::
Separations
===========

The on-sky separation can be computed with the
:meth:`astropy.coordinates.BaseCoordinateFrame.separation` or
:meth:`astropy.coordinates.SkyCoord.separation` methods,
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')
Expand Down Expand Up @@ -48,8 +47,7 @@ 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` or
:meth:`astropy.coordinates.SkyCoord.separation_3d` methods will
:meth:`~astropy.coordinates.BaseCoordinateFrame.separation_3d` will
determine the 3D distance between two coordinates that have ``distance``
defined::

Expand Down Expand Up @@ -80,8 +78,8 @@ astronomy convention (positive angles East of North)::
>>> c1.position_angle(c2).to(u.deg) # doctest: +FLOAT_CMP
<Angle 44.97818294 deg>

The combination of :meth:`~astropy.coordinates.SkyCoord.separation` and
:meth:`~astropy.coordinates.SkyCoord.position_angle` thus give a set of
The combination of :meth:`~astropy.coordinates.BaseCoordinateFrame.separation`
and :meth:`~astropy.coordinates.SkyCoord.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::
Expand Down
6 changes: 3 additions & 3 deletions docs/coordinates/skycoord.rst
Original file line number Diff line number Diff line change
Expand Up @@ -1134,7 +1134,7 @@ Matching Within Tolerance
-------------------------

To test if coordinates are within a certain angular distance of one other, use the
`~astropy.coordinates.SkyCoord.separation` method::
:meth:`~astropy.coordinates.BaseCoordinateFrame.separation` method::

>>> sc1.icrs.fk4.separation(sc1).to(u.arcsec) # doctest: +SKIP
<Angle 7.98873629e-13 arcsec>
Expand Down Expand Up @@ -1251,8 +1251,8 @@ the available docstrings below:
- `~astropy.coordinates.SkyCoord.match_to_catalog_sky`,
- `~astropy.coordinates.SkyCoord.match_to_catalog_3d`,
- `~astropy.coordinates.SkyCoord.position_angle`,
- `~astropy.coordinates.SkyCoord.separation`,
- `~astropy.coordinates.SkyCoord.separation_3d`
- `~astropy.coordinates.BaseCoordinateFrame.separation`,
- `~astropy.coordinates.BaseCoordinateFrame.separation_3d`
- `~astropy.coordinates.SkyCoord.apply_space_motion`

Additional information and examples can be found in the section on
Expand Down

0 comments on commit 6c6a50b

Please sign in to comment.