SunPy coordinates¶
The SunPy coordinates submodule is an implementation of the common solar physics coordinate frames using the Astropy coordinates framework.
Getting Started¶
The easiest interface to the coordinates module is through the SkyCoord
class:
>>> import astropy.units as u
>>> from astropy.coordinates import SkyCoord
>>> from sunpy.coordinates import frames
>>> c = SkyCoord(-100*u.arcsec, 500*u.arcsec, frame=frames.Helioprojective)
>>> c = SkyCoord(x=-72241.0*u.km, y=361206.1*u.km, z=589951.4*u.km, frame=frames.Heliocentric)
>>> c = SkyCoord(70*u.deg, -30*u.deg, frame=frames.HeliographicStonyhurst)
>>> c
<SkyCoord (HeliographicStonyhurst: obstime=None): (lon, lat, radius) in (deg, deg, km)
(70., -30., 695700.)>
It is also possible to use strings to define the frame but in that case make sure to
explicitly import sunpy.coordinates
as it registers solar coordinate frames with
astropy coordinates.:
>>> import astropy.units as u
>>> from astropy.coordinates import SkyCoord
>>> import sunpy.coordinates
>>> c = SkyCoord(-100*u.arcsec, 500*u.arcsec, frame='helioprojective')
>>> c
<SkyCoord (Helioprojective: obstime=None, rsun=695700.0 km, observer=earth): (Tx, Ty) in arcsec
(-100., 500.)>
SunPy implements support for the following solar physics coordinate systems:
Helioprojective (Cartesian)
Helioprojective
Heliocentric
Heliocentric
Heliographic Stonyhurst
HeliographicStonyhurst
Heliographic Carrington
HeliographicCarrington
for a complete description of these frames see sunpy.coordinates.frames
, for
a more detailed description of the frames see Thompson (2006)
SkyCoord
and all other coordinates
objects
also support array coordinates. These work the same as single-value coordinates,
but they store multiple coordinates in a single object. When you’re going to
apply the same operation to many different coordinates, this is a better choice
than a list of SkyCoord
objects, because it will be
much faster than applying the operation to each
SkyCoord
in a for
loop.
>>> c = SkyCoord([-500, 400]*u.arcsec, [100, 200]*u.arcsec, frame=frames.Helioprojective)
>>> c
<SkyCoord (Helioprojective: obstime=None, rsun=695700.0 km, observer=earth): (Tx, Ty) in arcsec
[(-500., 100.), ( 400., 200.)]>
>>> c[0]
<SkyCoord (Helioprojective: obstime=None, rsun=695700.0 km, observer=earth): (Tx, Ty) in arcsec
(-500., 100.)>
Accessing Coordinates¶
Individual coordinates can be accessed via attributes on the SkyCoord object,
but the names of the components of the coordinates for each frame differ. For a full
description of all the properties of the frames see sunpy.coordinates.frames
.
Helioprojective
¶
For the helioprojective frame the coordinates are accessed as Tx
and Ty
representing theta x and y. These are the same coordinates that are often
referred to as ‘solar-x’ and ‘solar-y’.:
>>> c = SkyCoord(-500*u.arcsec, 100*u.arcsec, frame=frames.Helioprojective)
>>> c.Tx
<Longitude -500. arcsec>
>>> c.Ty
<Latitude 100. arcsec>
Heliocentric
¶
Heliocentric normally a Cartesian frame so the coordinates are accessed as x, y, z
:
>>> c = SkyCoord(-72241.0*u.km, 361206.1*u.km, 589951.4*u.km, frame=frames.Heliocentric)
>>> c.x
<Quantity -72241. km>
>>> c.y
<Quantity 361206.1 km>
>>> c.z
<Quantity 589951.4 km>
HeliographicStonyhurst
and HeliographicCarrington
¶
Both the heliographic frames use latitude, longitude and radius which are accessed as follows:
>>> c = SkyCoord(70*u.deg, -30*u.deg, frame=frames.HeliographicStonyhurst)
>>> c.lat
<Latitude -30. deg>
>>> c.lon
<Longitude 70. deg>
>>> c.radius
<Distance 695700. km>
Transforming Between Coordinate Frames¶
Both SkyCoord
and
BaseCoordinateFrame
instances have a
transform_to
method. This can be used to
transform the frame to any other frame, either implemented in SunPy or in
Astropy the astropy documentation for more
details. An example of transforming the center of the solar disk to Carrington
coordinates is:
>>> c = SkyCoord(0*u.arcsec, 0*u.arcsec, frame=frames.Helioprojective, obstime="2017-07-26")
>>> c
<SkyCoord (Helioprojective: obstime=2017-07-26T00:00:00.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate for 'earth'>): (Tx, Ty) in arcsec
(0., 0.)>
>>> c.transform_to(frames.HeliographicCarrington)
<SkyCoord (HeliographicCarrington: obstime=2017-07-26T00:00:00.000): (lon, lat, radius) in (deg, deg, AU)
(283.95274241, 5.31701821, 0.00465047)>
It is also possible to transform to any coordinate system implemented in Astropy. This can be used to find the position of the solar limb in AltAz equatorial coordinates:
>>> from astropy.coordinates import EarthLocation, AltAz
>>> time = '2017-07-11 15:00'
>>> greenbelt = EarthLocation(lat=39.0044*u.deg, lon=-76.8758*u.deg)
>>> greenbelt_frame = AltAz(obstime=time, location=greenbelt)
>>> west_limb = SkyCoord(900*u.arcsec, 0*u.arcsec, frame=frames.Helioprojective, obstime=time)
>>> west_limb.transform_to(greenbelt_frame)
<SkyCoord (AltAz: obstime=2017-07-11 15:00:00.000, location=(1126916.53031967, -4833386.58391627, 3992696.622115747) m, pressure=0.0 hPa, temperature=0.0 deg_C, relative_humidity=0, obswl=1.0 micron): (az, alt, distance) in (deg, deg, m)
(111.40839101, 57.16645715, 1.51860261e+11)>
Observer Location Information¶
Both Helioprojective
and
Heliocentric
frames are defined by the location of
the observer. For example in Helioprojective
the
observer is at the origin of the coordinate system. This information is encoded
in the Helioprojective
and
Heliocentric
frames as the observer
attribute,
which is itself an instance of the
HeliographicStonyhurst
frame. The default observer
location is set to the position of the Earth (using
get_body_heliographic_stonyhurst
) as long as the
obstime
attribute is specified. If the obstime
attribute is not set then
you will be unable to transform the frame unless an explicit observer
is
specified, as the time is required to calculate the location of the Earth. The
location of the observer is automatically populated from meta data when
coordinate frames are created using map.
It is possible to convert from a Helioprojective
frame with one observer location to another
Helioprojective
frame with a different observer
location, by converting through Heliographic
, this
does involve making an assumption of the radius of the Sun to calculate the
position on the solar sphere. The conversion can be performed as follows:
# Input coordinate
>>> hpc1 = SkyCoord(0*u.arcsec, 0*u.arcsec, observer="earth", obstime="2017-07-26", frame=frames.Helioprojective)
Define a new Helioprojective frame with a different observer.
>>> import sunpy.coordinates
>>> hpc_out = sunpy.coordinates.Helioprojective(observer="venus", obstime="2017-07-26")
Perform the transformation from one to the other.
>>> hpc2 = hpc1.transform_to(hpc_out)
An example with two maps, named aia
and stereo
:
>>> hpc1 = SkyCoord(0*u.arcsec, 0*u.arcsec, frame=aia.coordinate_frame)
>>> hpc2 = hpc1.transform_to(stereo.coordinate_frame)
Design of the Coordinates Module¶
This module works by defining a collection of Frames
(sunpy.coordinates.frames
), which exists on a transformation graph, where the
transformations between the coordinate frames are then defined and registered
with the transformation graph (sunpy.coordinates.transformations
). It is also
possible to transform SunPy frames to Astropy frames.
Positions within these Frames
are stored as a Representation
of a
coordinate, a representation being a description of a point in a Cartesian,
spherical or cylindrical system (astropy.coordinates.representation
). A frame
that contains a representation of one or many points is said to have been
‘realized’.
For a more in depth look at the design and concepts of the Astropy coordinates system see Overview of astropy.coordinates Concepts
Frames and SkyCoord¶
The SkyCoord
class is a high level wrapper around the
astropy.coordinates
package. It provides an easier way to create and transform
coordinates, by using string representations for frames rather than the classes
themselves and some other usability improvements, for more information see the
SkyCoord
documentation.
The main advantage provided by SkyCoord
is the support it
provides for caching Frame attributes. Frame attributes are extra data specified
with a frame, some examples in sunpy.coordinates
are obstime
or
observer
for observer location. Only the frames where this data is
meaningful have these attributes, i.e. only the Helioprojective frames have
observer
. However, when you transform into another frame and then back to a
projective frame using SkyCoord
it will remember the attributes previously
provided, and repopulate the final frame with them. If you were to do
transformations using the Frames alone this would not happen.
The most important implication for this in sunpy.coordinates
is the rsun
parameter in the projective frames. If you create a projective frame with a
rsun
attribute, if you convert back to a projective frame it will be set
correctly. It should also be noted that, if you create a Heliographic frame and
then transform to a projective frame with an rsun
attribute, it will not
match the radius
coordinate in the Heliographic frame. This is because you may
mean to be describing a point above the defined ‘surface’ of the Sun.
Coordinates and WCS¶
The sunpy.coordinates
package provides a mapping between FITS-WCS CTYPE
convention and the coordinate frames as defined in sunpy.coordinates
. This is
used via the astropy.wcs.utils.wcs_to_celestial_frame
function, with which the
SunPy frames are registered upon being imported. This list is used by packages
such as wcsaxes
to convert from astropy.wcs.WCS
objects to coordinate
frames.
The sunpy.map.GenricMap
class creates astropy.wcs.WCS
objects as
amap.wcs
, however, it adds some extra attributes to the WCS
object to be able to fully specify the coordinate frame. It adds
heliographic_observer
and rsun
.
If you want to obtain a un-realized coordinate frame corresponding to a
GenericMap
object you can do the following:
>>> import sunpy.map
>>> from sunpy.data.sample import AIA_171_IMAGE
>>> amap = sunpy.map.Map(AIA_171_IMAGE)
>>> amap.observer_coordinate
<SkyCoord (HeliographicStonyhurst: obstime=2011-06-07T06:33:02.770): (lon, lat, radius) in (deg, deg, m)
(-0.00406308, 0.04787238, 1.51846026e+11)>
which is equivalent to:
>>> from astropy.wcs.utils import wcs_to_celestial_frame
>>> wcs_to_celestial_frame(amap.wcs)
<Helioprojective Frame (obstime=2011-06-07T06:33:02.770, rsun=696000000.0 m, observer=<HeliographicStonyhurst Coordinate (obstime=2011-06-07T06:33:02.770): (lon, lat, radius) in (deg, deg, m)
(-0.00406308, 0.04787238, 1.51846026e+11)>)>
sunpy.coordinates Package¶
This subpackage contains:
A robust framework for working with coordinate systems
Functions to obtain the locations of solar-system bodies
Functions to calculate Sun-specific coordinate information (
sunpy.coordinates.sun
)
The diagram below shows all of Sun-based and Earth-based coordinate systems
available through sunpy.coordinates
, as well as the transformations between
them. Each frame is labeled with the name of its class and its alias (useful
for converting other coordinates to them using attribute-style access).
The frames colored in cyan are implemented in astropy.coordinates
, and there
are other astronomical frames that can be transformed to that are not shown
below. See the documentation for astropy.coordinates
for more information.
-
SunPy frames:
-
Astropy frames:
-
AffineTransform: ➝
-
FunctionTransform: ➝
-
FunctionTransformWithFiniteDifference: ➝
-
StaticMatrixTransform: ➝
-
DynamicMatrixTransform: ➝
Functions¶
|
Return a |
|
Return a |
|
Queries JPL HORIZONS and returns a |
|
Deprecated since version 1.0. |
|
Deprecated since version 1.0. |
|
Deprecated since version 1.0. |
|
Deprecated since version 1.0. |
|
Deprecated since version 1.0. |
|
For a given frame, this function returns the corresponding WCS object. |
This function registers the coordinates frames to their FITS-WCS coordinate type values in the |
Classes¶
|
A coordinate or frame in the Heliocentric system, which is observer-based. |
|
A coordinate or frame in the Carrington Heliographic (HGC) system. |
|
A coordinate or frame in the Stonyhurst Heliographic (HGS) system. |
|
A coordinate or frame in the Helioprojective Cartesian (HPC) system, which is observer-based. |
A frame which is relative to some position and another frame. |
Class Inheritance Diagram¶
sunpy.coordinates.frames Module¶
Common solar physics coordinate systems.
This submodule implements various solar physics coordinate frames for use with
the astropy.coordinates
module.
Classes¶
|
A coordinate or frame in the Stonyhurst Heliographic (HGS) system. |
|
A coordinate or frame in the Carrington Heliographic (HGC) system. |
|
A coordinate or frame in the Heliocentric system, which is observer-based. |
|
A coordinate or frame in the Helioprojective Cartesian (HPC) system, which is observer-based. |
Class Inheritance Diagram¶
sunpy.coordinates.transformations Module¶
Coordinate Transformation Functions
This module contains the functions for converting one
sunpy.coordinates.frames
object to another.
Warning
The functions in this submodule should never be called directly, transforming
between coordinate frames should be done using the .transform_to
methods
on BaseCoordinateFrame
or
SkyCoord
instances.
Functions¶
|
Convert from Heliographic Stonyhurst to Heliographic Carrington. |
|
Convert from Heliographic Carrington to Heliographic Stonyhurst. |
|
Convert from Heliocentric Cartesian to Helioprojective Cartesian. |
|
Convert from Helioprojective Cartesian to Heliocentric Cartesian. |
|
Convert from Heliocentric Cartesian to Heliographic Stonyhurst. |
|
Convert from Heliographic Stonyhurst to Heliocentric Cartesian. |
|
This converts from HPC to HPC, with different observer location parameters. |
|
Convert from HCRS to Heliographic Stonyhurst (HGS). |
|
Convert from Heliographic Stonyhurst to HCRS. |
|
Convert between two Heliographic Stonyhurst frames. |
|
Convert between two Heliographic Carrington frames. |
|
Convert between two Heliocentric frames. |
sunpy.coordinates.ephemeris Module¶
Ephemeris calculations using SunPy coordinate frames
Functions¶
|
Return a |
|
Return a |
|
Deprecated since version 1.0. |
|
Deprecated since version 1.0. |
|
Deprecated since version 1.0. |
|
Deprecated since version 1.0. |
|
Deprecated since version 1.0. |
|
Queries JPL HORIZONS and returns a |
sunpy.coordinates.sun Module¶
Sun-specific coordinate calculations
Functions¶
|
Return the angular radius of the Sun as viewed from Earth. |
|
Returns the apparent position of the Sun (right ascension and declination) on the celestial sphere using the equatorial coordinate system, referred to the true equinox of date (as default). |
Return the Carrington rotation number. |
|
|
Returns the Sun’s true geometric longitude, referred to the mean equinox of date. |
|
Returns the Sun’s apparent longitude, referred to the true equinox of date. |
|
Returns the Sun’s true geometric latitude, referred to the mean equinox of date. |
|
Returns the Sun’s apparent latitude, referred to the true equinox of date. |
Returns the mean obliquity of the ecliptic, using the IAU 2006 definition. |
|
|
Returns the Sun’s true geometric right ascension relative to Earth, referred to the mean equinox of date (as default). |
|
Returns the Sun’s true geometric declination relative to Earth, referred to the mean equinox of date (as default). |
Returns the true obliquity of the ecliptic, using the IAU 2006 definition. |
|
|
Returns the Sun’s apparent right ascension relative to Earth, referred to the true equinox of date (as default). |
|
Returns the Sun’s apparent declination relative to Earth, referred to the true equinox of date (as default). |
|
Print out a summary of solar ephemeris. |
|
Return the B0 angle for the Sun at a specified time, which is the heliographic latitude of the Sun-disk center as seen from Earth. |
|
Return the L0 angle for the Sun at a specified time, which is the Carrington longitude of the Sun-disk center as seen from Earth. |
|
Return the position (P) angle for the Sun at a specified time, which is the angle between geocentric north and solar north as seen from Earth, measured eastward from geocentric north. |
|
Return the distance between the Sun and the Earth at a specified time. |
|
Return the orientation angle for the Sun from a specified Earth location and time. |
sunpy.coordinates.offset_frame Module¶
Classes¶
A frame which is relative to some position and another frame. |
Class Inheritance Diagram¶
sunpy.coordinates.wcs_utils Module¶
Functions¶
This function registers the coordinates frames to their FITS-WCS coordinate type values in the |
|
|
For a given frame, this function returns the corresponding WCS object. |
Attribution¶
Some of this documentation was adapted from Astropy under the terms of the BSD License.
This package was initially developed by Pritish Chakraborty as part of GSOC 2014 and Stuart Mumford.