Skip to content

Commit

Permalink
Allow setting radius in Helioprojective.assume_spherical_screen
Browse files Browse the repository at this point in the history
  • Loading branch information
svank committed Mar 26, 2024
1 parent 4c17a50 commit 150840b
Show file tree
Hide file tree
Showing 3 changed files with 53 additions and 7 deletions.
1 change: 1 addition & 0 deletions changelog/7532.feature.rst
Original file line number Diff line number Diff line change
@@ -0,0 +1 @@
Allow the screen radius to be set when using :meth:`~sunpy.coordinates.Helioprojective.assume_spherical_screen`.
20 changes: 13 additions & 7 deletions sunpy/coordinates/frames.py
Original file line number Diff line number Diff line change
Expand Up @@ -679,14 +679,14 @@ def is_visible(self, *, tolerance: u.m = 1*u.m):

@classmethod
@contextmanager
def assume_spherical_screen(cls, center, only_off_disk=False):
def assume_spherical_screen(cls, center, only_off_disk=False, radius='center_to_sun'):
"""
Context manager to interpret 2D coordinates as being on the inside of a spherical screen.
The radius of the screen is the distance between the specified ``center`` and Sun center.
This ``center`` does not have to be the same as the observer location for the coordinate
frame. If they are the same, then this context manager is equivalent to assuming that the
helioprojective "zeta" component is zero.
The radius of the screen is, by default, the distance between the specified ``center`` and
Sun center. This ``center`` does not have to be the same as the observer location for the
coordinate frame. If they are the same, then this context manager is equivalent to assuming
that the helioprojective "zeta" component is zero.
This replaces the default assumption where 2D coordinates are mapped onto the surface of the
Sun.
Expand All @@ -698,6 +698,9 @@ def assume_spherical_screen(cls, center, only_off_disk=False):
only_off_disk : `bool`, optional
If `True`, apply this assumption only to off-disk coordinates, with on-disk coordinates
still mapped onto the surface of the Sun. Defaults to `False`.
radius : `~astropy.units.Quantity`
The radius of the spherical screen. The default sets the radius to the distance from the
screen center to the Sun.
Examples
--------
Expand Down Expand Up @@ -733,10 +736,13 @@ def assume_spherical_screen(cls, center, only_off_disk=False):
try:
old_spherical_screen = cls._spherical_screen # nominally None

center_hgs = center.transform_to(HeliographicStonyhurst(obstime=center.obstime))
if radius == 'center_to_sun':
center_hgs = center.transform_to(HeliographicStonyhurst(obstime=center.obstime))
radius = center_hgs.radius

Check warning on line 741 in sunpy/coordinates/frames.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/frames.py#L739-L741

Added lines #L739 - L741 were not covered by tests

cls._spherical_screen = {
'center': center,
'radius': center_hgs.radius,
'radius': radius,
'only_off_disk': only_off_disk
}
yield
Expand Down
39 changes: 39 additions & 0 deletions sunpy/coordinates/tests/test_transformations.py
Original file line number Diff line number Diff line change
Expand Up @@ -120,6 +120,45 @@ def test_hpc_hpc_null():
assert hpc_out.observer == hpc_new.observer


def test_hpc_hpc_spherical_screen():
D0 = 20*u.R_sun
L0 = 67.5*u.deg
Tx0 = 45*u.deg
observer_in = HeliographicStonyhurst(lat=0*u.deg, lon=0*u.deg, radius=D0)

Check warning on line 127 in sunpy/coordinates/tests/test_transformations.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/tests/test_transformations.py#L124-L127

Added lines #L124 - L127 were not covered by tests
# Once our coordinate is placed on the screen, observer_out will be looking
# directly along the line containing itself, the coordinate, and the Sun
observer_out = HeliographicStonyhurst(lat=0*u.deg, lon=L0, radius=2*D0)

Check warning on line 130 in sunpy/coordinates/tests/test_transformations.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/tests/test_transformations.py#L130

Added line #L130 was not covered by tests

sc_in = SkyCoord(Tx0, 0*u.deg, observer=observer_in,

Check warning on line 132 in sunpy/coordinates/tests/test_transformations.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/tests/test_transformations.py#L132

Added line #L132 was not covered by tests
frame='helioprojective')

with Helioprojective.assume_spherical_screen(observer_in):
sc_3d = sc_in.make_3d()
sc_out = sc_in.transform_to(Helioprojective(observer=observer_out))

Check warning on line 137 in sunpy/coordinates/tests/test_transformations.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/tests/test_transformations.py#L135-L137

Added lines #L135 - L137 were not covered by tests

assert quantity_allclose(sc_3d.distance, D0)

Check warning on line 139 in sunpy/coordinates/tests/test_transformations.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/tests/test_transformations.py#L139

Added line #L139 was not covered by tests

assert quantity_allclose(sc_out.Tx, 0*u.deg, atol=1e-6*u.deg)
assert quantity_allclose(sc_out.Ty, 0*u.deg, atol=1e-6*u.deg)

Check warning on line 142 in sunpy/coordinates/tests/test_transformations.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/tests/test_transformations.py#L141-L142

Added lines #L141 - L142 were not covered by tests
# Law of Cosines to compute the coordinate's distance from the Sun, and then
# r_expected is the distance from observer_out to the coordinate
radius_expected = 2 * D0 - np.sqrt(2 * D0**2 - 2 * D0 * D0 * np.cos(45*u.deg))
assert quantity_allclose(sc_out.distance, radius_expected)

Check warning on line 146 in sunpy/coordinates/tests/test_transformations.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/tests/test_transformations.py#L145-L146

Added lines #L145 - L146 were not covered by tests

# Now test with a very large screen, letting us approximate the two
# observers as being the same (aside from a different zero point for Tx)
r_s = 1e9 * u.lightyear
with Helioprojective.assume_spherical_screen(observer_in, radius=r_s):
sc_3d = sc_in.make_3d()
sc_out = sc_in.transform_to(Helioprojective(observer=observer_out))

Check warning on line 153 in sunpy/coordinates/tests/test_transformations.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/tests/test_transformations.py#L150-L153

Added lines #L150 - L153 were not covered by tests

assert quantity_allclose(sc_3d.distance, r_s)

Check warning on line 155 in sunpy/coordinates/tests/test_transformations.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/tests/test_transformations.py#L155

Added line #L155 was not covered by tests

assert quantity_allclose(sc_out.Tx, Tx0 + L0)
assert quantity_allclose(sc_out.Ty, 0*u.deg)
assert quantity_allclose(sc_out.distance, r_s)

Check warning on line 159 in sunpy/coordinates/tests/test_transformations.py

View check run for this annotation

Codecov / codecov/patch

sunpy/coordinates/tests/test_transformations.py#L157-L159

Added lines #L157 - L159 were not covered by tests


def test_hcrs_hgs():
# Get the current Earth location in HCRS
adate = parse_time('2015/05/01 01:13:00')
Expand Down

0 comments on commit 150840b

Please sign in to comment.