HeliographicStonyhurst#

class sunpy.coordinates.frames.HeliographicStonyhurst(*args, **kwargs)[source]#

Bases: BaseHeliographic

A coordinate or frame in the Stonyhurst Heliographic (HGS) system.

  • The origin is the center of the Sun.

  • The Z-axis (+90 degrees latitude) is aligned with the Sun’s north pole.

  • The X-axis (0 degrees longitude and 0 degrees latitude) is aligned with the projection of the Sun-Earth line onto the Sun’s equatorial plane.

This system is also know as the Heliocentric Earth Equatorial (HEEQ) system when represented using Cartesian components.

A new instance can be created using the following signatures (note that if supplied, obstime and representation_type must be keyword arguments):

HeliographicStonyhurst(lon, lat, obstime=obstime)
HeliographicStonyhurst(lon, lat, radius, obstime=obstime)
HeliographicStonyhurst(x, y, z, representation_type='cartesian', obstime=obstime)
Parameters:
  • data (BaseRepresentation or None) – A representation object or None to have no data (or use the coordinate component arguments, see below).

  • lon (Angle or Quantity, optional) – The longitude coordinate for this object (lat must also be given and data must be None). Not needed if data is given.

  • lat (Angle or Quantity, optional) – The latitude coordinate for this object (lon must also be given and data must be None). Not needed if data is given.

  • radius (Quantity, optional) – The radial distance coordinate from Sun center for this object. Defaults to the radius of the Sun. Not needed if data is given.

  • rsun (Quantity) – The radius of the Sun in length units. Used to convert a 2D coordinate (i.e., no radius component) to a 3D coordinate by assuming that the coordinate is on the surface of the Sun. Defaults to the photospheric radius as defined in sunpy.sun.constants.

  • obstime (tuple, list, str, pandas.Timestamp, pandas.Series, pandas.DatetimeIndex, datetime.datetime, datetime.date, numpy.datetime64, numpy.ndarray, astropy.time.Time) – The time of the observation. This is used to determine the position of solar-system bodies (e.g., the Sun and the Earth) as needed to define the origin and orientation of the frame.

  • representation_type (BaseRepresentation, str, optional) – A representation class or string name of a representation class. This may change the valid coordinate component arguments from the defaults (see above). For example, passing representation_type='cartesian' will make the frame expect Cartesian coordinate component arguments (typically, x, y, and z).

  • copy (bool, optional) – If True (default), make copies of the input coordinate arrays.

Examples

>>> from astropy.coordinates import SkyCoord
>>> import sunpy.coordinates
>>> import astropy.units as u
>>> sc = SkyCoord(1*u.deg, 1*u.deg, 2*u.km,
...               frame="heliographic_stonyhurst",
...               obstime="2010/01/01T00:00:45")
>>> sc
<SkyCoord (HeliographicStonyhurst: obstime=2010-01-01T00:00:45.000, rsun=695700.0 km): (lon, lat, radius) in (deg, deg, km)
    (1., 1., 2.)>
>>> sc.frame
<HeliographicStonyhurst Coordinate (obstime=2010-01-01T00:00:45.000, rsun=695700.0 km): (lon, lat, radius) in (deg, deg, km)
    (1., 1., 2.)>
>>> sc = SkyCoord(HeliographicStonyhurst(-10*u.deg, 2*u.deg))
>>> sc
<SkyCoord (HeliographicStonyhurst: obstime=None, rsun=695700.0 km): (lon, lat) in deg
    (-10., 2.)>
>>> sc = SkyCoord(CartesianRepresentation(0*u.km, 45*u.km, 2*u.km),
...               obstime="2011/01/05T00:00:50",
...               frame="heliographic_stonyhurst")
>>> sc
<SkyCoord (HeliographicStonyhurst: obstime=2011-01-05T00:00:50.000, rsun=695700.0 km): (lon, lat, radius) in (deg, deg, km)
(90., 2.54480438, 45.04442252)>

Attributes Summary

default_differential

Default representation for differential data (e.g., velocity)

default_representation

Default representation for position data

frame_attributes

frame_specific_representation_info

Mapping for frame-specific component names

name

Attributes Documentation

default_differential#

Default representation for differential data (e.g., velocity)

default_representation#

Default representation for position data

frame_attributes = {'obstime': <sunpy.coordinates.frameattributes.TimeFrameAttributeSunPy object>, 'rsun': <astropy.coordinates.attributes.QuantityAttribute object>}#
frame_specific_representation_info#

Mapping for frame-specific component names

name = 'heliographic_stonyhurst'#