Helioprojective

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

Bases: sunpy.coordinates.frames.SunPyBaseCoordinateFrame

A coordinate or frame in the Helioprojective Cartesian (HPC) system, which is observer-based.

  • The origin is the location of the observer.

  • theta_x is the angle relative to the plane containing the Sun-observer line and the Sun’s rotation axis, with positive values in the direction of the Sun’s west limb.

  • theta_y is the angle relative to the Sun’s equatorial plane, with positive values in the direction of the Sun’s north pole.

  • distance is the Sun-observer distance.

This system is frequently used in a projective form without distance specified. For observations looking very close to the center of the Sun, where the small-angle approximation is appropriate, theta_x and theta_y can be approximated as Cartesian components.

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

Helioprojective(theta_x, theta_y, obstime=obstime, observer=observer)
Helioprojective(theta_x, theta_y, distance, obstime=obstime, observer=observer)
Parameters
  • data (BaseRepresentation or None) – A representation object or None to have no data (or use the coordinate component arguments, see below).

  • Tx (Angle or Quantity) – The theta_x coordinate for this object. Not needed if data is given.

  • Ty (Angle or Quantity) – The theta_y coordinate for this object. Not needed if data is given.

  • distance (Quantity) – The distance coordinate from the observer 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.

  • rsun (Quantity) – The physical (length) radius of the Sun. Used to calculate the position of the limb for calculating distance from the observer to the coordinate. Defaults to the solar radius.

  • 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
>>> import sunpy.coordinates
>>> import astropy.units as u
>>> sc = SkyCoord(0*u.deg, 0*u.deg, 5*u.km,
...               obstime="2010/01/01T00:00:00", observer="earth", frame="helioprojective")
>>> sc
<SkyCoord (Helioprojective: obstime=2010-01-01T00:00:00.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate for 'earth'>): (Tx, Ty, distance) in (arcsec, arcsec, km)
    (0., 0., 5.)>
>>> sc = SkyCoord(0*u.deg, 0*u.deg,
...               obstime="2010/01/01T00:00:00", observer="earth", frame="helioprojective")
>>> sc
<SkyCoord (Helioprojective: obstime=2010-01-01T00:00:00.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate for 'earth'>): (Tx, Ty) in arcsec
    (0., 0.)>
>>> sc = SkyCoord(CartesianRepresentation(1*u.AU, 1e5*u.km, -2e5*u.km),
...               obstime="2011/01/05T00:00:50", observer="earth", frame="helioprojective")
>>> sc
<SkyCoord (Helioprojective: obstime=2011-01-05T00:00:50.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate for 'earth'>): (Tx, Ty, distance) in (arcsec, arcsec, AU)
    (137.87948623, -275.75878762, 1.00000112)>

Attributes Summary

default_differential

default_representation

frame_attributes

frame_specific_representation_info

name

observer

rsun

Methods Summary

calculate_distance()

Deprecated since version 1.1.

make_3d()

This method calculates the third coordinate of the Helioprojective frame.

Attributes Documentation

default_differential
default_representation
frame_attributes = {'observer': <sunpy.coordinates.frameattributes.ObserverCoordinateAttribute object>, 'obstime': <sunpy.coordinates.frameattributes.TimeFrameAttributeSunPy object>, 'rsun': <astropy.coordinates.attributes.Attribute object>}
frame_specific_representation_info
name = 'helioprojective'
observer = None
rsun = <Quantity 695700. km>

Methods Documentation

calculate_distance() [edit on github]

Deprecated since version 1.1: This is scheduled for removal in 2.1.

This method calculates the third coordinate of the Helioprojective frame. It assumes that the coordinate point is on the surface of the Sun.

If a point in the frame is off limb then NaN will be returned.

Returns

new_frame (Helioprojective) – A new frame instance with all the attributes of the original but now with a third coordinate.

make_3d()[source] [edit on github]

This method calculates the third coordinate of the Helioprojective frame. It assumes that the coordinate point is on the surface of the Sun.

If a point in the frame is off limb then NaN will be returned.

Returns

new_frame (Helioprojective) – A new frame instance with all the attributes of the original but now with a third coordinate.