Helioprojective

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

Bases: astropy.coordinates.BaseCoordinateFrame

A coordinate or frame in the Helioprojective (Cartesian) system.

This is a projective coordinate system centered around the observer. It is a full spherical coordinate system with position given as longitude theta_x and latitude theta_y.

Parameters:
  • representation (BaseRepresentation or None.) – A representation object. If specified, other parameters must be in keyword form.
  • Tx (Angle or Quantity) – X-axis coordinate.
  • Ty (Angle or Quantity) – Y-axis coordinate.
  • distance (Quantity) – The radial distance from the observer to the coordinate point.
  • obstime (SunPy Time) – The date and time of the observation, used to convert to heliographic carrington coordinates.
  • observer (HeliographicStonyhurst) – The coordinate of the observer in the solar system.
  • 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.

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",
...               frame="helioprojective")
>>> sc
<SkyCoord (Helioprojective: obstime=2010-01-01 00:00:00, rsun=695508.0 km, observer=<HeliographicStonyhurst Coordinate (obstime=2010-01-01 00:00:00): (lon, lat, radius) in (deg, deg, AU)
    (0., -3.00724817, 0.98330294)>): (Tx, Ty, distance) in (arcsec, arcsec, km)
    (0., 0., 5.)>
>>> sc = SkyCoord(0*u.deg, 0*u.deg, obstime="2010/01/01T00:00:00", frame="helioprojective")
>>> sc
<SkyCoord (Helioprojective: obstime=2010-01-01 00:00:00, rsun=695508.0 km, observer=<HeliographicStonyhurst Coordinate (obstime=2010-01-01 00:00:00): (lon, lat, radius) in (deg, deg, AU)
    (0., -3.00724817, 0.98330294)>): (Tx, Ty) in arcsec
    (0., 0.)>

Attributes Summary

default_differential
default_representation
frame_attributes
frame_specific_representation_info
name
observer
obstime
rsun

Methods Summary

calculate_distance() This method calculates the third coordinate of the Helioprojective frame.

Attributes Documentation

default_differential
default_representation
frame_attributes = OrderedDict([('obstime', <sunpy.coordinates.frameattributes.TimeFrameAttributeSunPy object>), ('rsun', <astropy.coordinates.attributes.Attribute object>), ('observer', <sunpy.coordinates.frameattributes.ObserverCoordinateAttribute object>)])
frame_specific_representation_info
name = 'helioprojective'
observer = 'earth'
obstime = None
rsun = <Quantity 695508. km>

Methods Documentation

calculate_distance()[source] [edit on github]

This method calculates the third coordinate of the Helioprojective frame. It assumes that the coordinate point is on the disk of the Sun at the rsun radius.

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.