Skip to content

Commit

Permalink
Merge pull request #15659 from eerovaher/deduplicate-separation
Browse files Browse the repository at this point in the history
Ensure `BaseCoordinateFrame` separation methods work with `SkyCoord` input
  • Loading branch information
mhvk committed Dec 4, 2023
2 parents bf35e88 + 4f9bab5 commit 7d8bfaf
Show file tree
Hide file tree
Showing 3 changed files with 27 additions and 5 deletions.
13 changes: 8 additions & 5 deletions astropy/coordinates/baseframe.py
Original file line number Diff line number Diff line change
Expand Up @@ -1779,7 +1779,7 @@ def separation(self, other):
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 +1798,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 = (
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 @@ -1814,7 +1817,7 @@ def separation_3d(self, other):
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 +1838,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)

if issubclass(other_in_self_system.__class__, r.UnitSphericalRepresentation):
raise ValueError(
Expand Down
15 changes: 15 additions & 0 deletions astropy/coordinates/tests/test_frames.py
Original file line number Diff line number Diff line change
Expand Up @@ -736,6 +736,21 @@ def test_sep():
i7.separation_3d(i3)


@pytest.mark.parametrize(
"method,expectation",
[
pytest.param("separation", 0.69815121 * u.deg, id="separation"),
pytest.param("separation_3d", 0.12184962 * u.pc, id="separation_3d"),
],
)
def test_seps_with_skycoord(method, expectation):
coords = (1 * u.deg, 2 * u.deg, 10 * u.pc)
assert_allclose(
getattr(FK5(*coords), method)(SkyCoord(*coords, frame=FK5, equinox="B1950")),
expectation,
)


def test_time_inputs():
"""
Test validation and conversion of inputs for equinox and obstime attributes.
Expand Down
4 changes: 4 additions & 0 deletions docs/changes/coordinates/15659.bugfix.rst
Original file line number Diff line number Diff line change
@@ -0,0 +1,4 @@
Previously passing a ``SkyCoord`` instance to the ``BaseCoordinateFrame``
``separation()`` or ``separation_3d()`` methods could produce wrong results,
depending on what additional frame attributes were defined on the ``SkyCoord``,
but now ``SkyCoord`` input can be used safely.

0 comments on commit 7d8bfaf

Please sign in to comment.