solar_rotate_coordinate#

sunpy.physics.differential_rotation.solar_rotate_coordinate(coordinate, observer=None, time=None, **diff_rot_kwargs)[source]#

Given a coordinate on the Sun, calculate where that coordinate maps to as seen by a new observer at some later or earlier time, given that the input coordinate rotates according to the solar rotation profile.

The amount of solar rotation is based on the amount of time between the observation time of the input coordinate and the observation time of the new observer. The new observer is specified in one of two ways, either using the “observer” or “time” keywords.

If the “observer” keyword is set, it is used to specify the location of the new observer in space and time. The difference between the coordinate time and the new observer time is used to calculate the amount of solar rotation applied, and the location of the new observer in space is used to calculate where the rotated coordinate is as seen from the new observer.

If the “time” keyword is set, it is used to specify the number of seconds to rotate the coordinate by. Note that using the “time” keyword assumes that the new observer is on the Earth. This may be a reasonable assumption depending on the application.

Either the “observer” or “time” keyword must be specified, but both cannot be specified at the same time.

Parameters:
  • coordinate (SkyCoord) – Any valid coordinate which is transformable to Heliographic Stonyhurst.

  • observer (BaseCoordinateFrame, SkyCoord, None) – The location of the new observer in space and time (the observer must have an interpretable obstime property).

  • time (Time, TimeDelta, Quantity, None)

  • **diff_rot_kwargs (dict) – Keyword arguments are passed on as keyword arguments to diff_rot. Note that the keyword “frame_time” is automatically set to the value “sidereal”.

Returns:

coordinate (SkyCoord) – The locations of the input coordinates after the application of solar rotation as seen from the point-of-view of the new observer.

Notes

The translational motion of the Sun over the time interval will be ignored. See transform_with_sun_center().

Examples

>>> import astropy.units as u
>>> from astropy.coordinates import SkyCoord
>>> from sunpy.coordinates import Helioprojective, get_body_heliographic_stonyhurst
>>> from sunpy.physics.differential_rotation import solar_rotate_coordinate
>>> from sunpy.time import parse_time
>>> start_time = parse_time('2010-09-10 12:34:56')
>>> c = SkyCoord(-570*u.arcsec, 120*u.arcsec, obstime=start_time,
...              observer="earth", frame=Helioprojective)
>>> solar_rotate_coordinate(c, time=start_time + 25*u.hr)  
<SkyCoord (Helioprojective: obstime=2010-09-11T13:34:56.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate (obstime=2010-09-11T13:34:56.000, rsun=695700.0 km): (lon, lat, radius) in (deg, deg, AU)
    (-5.68434189e-14, 7.24318962, 1.00669016)>): (Tx, Ty, distance) in (arcsec, arcsec, AU)
    (-378.27830452, 105.70767875, 1.00245134)>
>>> new_observer = get_body_heliographic_stonyhurst("earth", start_time + 6*u.day)
>>> solar_rotate_coordinate(c, observer=new_observer)
<SkyCoord (Helioprojective: obstime=2010-09-16T12:34:56.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate (obstime=2010-09-16T12:34:56.000, rsun=695700.0 km): (lon, lat, radius) in (deg, deg, AU)
    (2.65061438e-14, 7.18706547, 1.00534174)>): (Tx, Ty, distance) in (arcsec, arcsec, AU)
    (620.42567049, 126.13662663, 1.00185786)>