Heliocentric

class sunpy.coordinates.Heliocentric(*args, **kwargs)[source] [edit on github]

Bases: sunpy.coordinates.frames.SunPyBaseCoordinateFrame

A coordinate or frame in the Heliocentric system, which is observer-based.

  • The origin is the center of the Sun.

  • The Z-axis is aligned with the Sun-observer line.

  • The Y-axis is aligned with the component of the vector to the Sun’s north pole that is perpendicular to the Z-axis.

This frame defaults to a Cartesian component representation, which is known as Heliocentric Cartesian (HCC). This frame can also be represented using cylindrical components, where where rho is the impact parameter and psi is the position angle. psi is measured relative to the west limb, rather than solar north, so is shifted by 90 degrees compared to the convention of the Heliocentric Radial (HCR) system.

A new instance can be created using the following signatures (note that obstime and representation_type must be supplied as keywords):

Heliocentric(x, y, z, obstime)
Heliocentric(rho, psi, z, obstime, representation_type='cylindrical')
Parameters
  • data (BaseRepresentation or None) – A representation object or None to have no data (or use the coordinate component arguments, see below).

  • x (Quantity, optional) – X-axis coordinate for this object. Not needed if data is given.

  • y (Quantity, optional) – Y-axis coordinate for this object. Not needed if data is given.

  • z (Quantity, optional) – Z-axis coordinate for this object. Not needed if data is given.

  • observer (HeliographicStonyhurst, str) – The location of the observer. If a string is provided, it must be a solar system body that can be parsed by get_body_heliographic_stonyhurst at the time obstime. Defaults to Earth center.

  • obstime (tuple, 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, CartesianRepresentation
>>> import sunpy.coordinates
>>> import astropy.units as u
>>> sc = SkyCoord(CartesianRepresentation(10*u.km, 1*u.km, 2*u.km),
...               obstime="2011/01/05T00:00:50", frame="heliocentric")
>>> sc
<SkyCoord (Heliocentric: obstime=2011-01-05T00:00:50.000, observer=<HeliographicStonyhurst Coordinate for 'earth'>): (x, y, z) in km
    (10., 1., 2.)>
>>> sc = SkyCoord([1,2]*u.km, [3,4]*u.m, [5,6]*u.cm, frame="heliocentric", obstime="2011/01/01T00:00:54")
>>> sc
<SkyCoord (Heliocentric: obstime=2011-01-01T00:00:54.000, observer=<HeliographicStonyhurst Coordinate for 'earth'>): (x, y, z) in (km, m, cm)
    [(1., 3., 5.), (2., 4., 6.)]>
>>> sc = SkyCoord(CylindricalRepresentation(10*u.km, 60*u.deg, 10*u.km),
...               obstime="2011/01/05T00:00:50",
...               frame="heliocentric")
>>> sc
<SkyCoord (Heliocentric: obstime=2011-01-05T00:00:50.000, observer=<HeliographicStonyhurst Coordinate for 'earth'>): (x, y, z) in km
    (5., 8.66025404, 10.)>

Attributes Summary

default_differential

default_representation

frame_attributes

frame_specific_representation_info

name

observer

Attributes Documentation

default_differential
default_representation
frame_attributes = {'observer': <sunpy.coordinates.frameattributes.ObserverCoordinateAttribute object>, 'obstime': <sunpy.coordinates.frameattributes.TimeFrameAttributeSunPy object>}
frame_specific_representation_info
name = 'heliocentric'
observer = 'earth'