HelioprojectiveRadial#

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

Bases: SunPyBaseCoordinateFrame

A coordinate or frame in the Helioprojective Radial system.

This is an observer-based spherical coordinate system, with:

  • psi is the position angle of the coordinate, measured eastward from solar north

  • delta is the declination angle, which is the impact angle (the angle between the observer-Sun line and the observer-coordinate line) minus 90 degrees

  • r is the observer-coordinate distance

Note

The declination angle, rather than the impact angle, is used as a component in order to match the FITS WCS definition. The impact angle can be readily retrieved using the theta property.

Parameters:
  • data (BaseRepresentation or None) – A representation object or None to have no data (or use the coordinate component arguments, see below).

  • psi (Angle or Quantity) – The position angle. Not needed if data is given.

  • delta (Angle or Quantity) – The declination angle. Not needed if data is given.

  • r (Angle or Quantity) – The observer-coordinate distance. 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.

  • 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, 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.

See also

Helioprojective

Examples

>>> from astropy.coordinates import SkyCoord
>>> import sunpy.coordinates
>>> import astropy.units as u
>>> sc = SkyCoord(0*u.deg, -90*u.deg, 5*u.km,
...               obstime="2010/01/01T00:00:00", observer="earth", frame="helioprojectiveradial")
>>> sc
<SkyCoord (HelioprojectiveRadial: obstime=2010-01-01T00:00:00.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate for 'earth'>): (psi, delta, distance) in (deg, deg, km)
    (0., -90., 5.)>
>>> sc.theta
<Angle 0. arcsec>
>>> sc = SkyCoord(30*u.deg, -89.9*u.deg,
...               obstime="2010/01/01T00:00:00", observer="earth", frame="helioprojectiveradial")
>>> sc
<SkyCoord (HelioprojectiveRadial: obstime=2010-01-01T00:00:00.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate for 'earth'>): (psi, delta) in deg
    (30., -89.9)>
>>> sc.theta
<Angle 360. arcsec>
>>> sc = SkyCoord(CartesianRepresentation(1e5*u.km, -2e5*u.km, -1*u.AU),
...               obstime="2011/01/05T00:00:50", observer="earth", frame="helioprojectiveradial")
>>> sc
<SkyCoord (HelioprojectiveRadial: obstime=2011-01-05T00:00:50.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate for 'earth'>): (psi, delta, distance) in (deg, deg, km)
    (296.56505118, -89.91435897, 1.49598038e+08)>
>>> sc.theta
<Angle 308.30772022 arcsec>

Reprojecting a Helioprojective Map to Helioprojective Radial

Reprojecting a Helioprojective Map to Helioprojective Radial

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

observer

A frame attribute

rsun

A frame attribute

theta

Returns the impact angle, which is the declination angle plus 90 degrees.

Methods Summary

make_3d()

Returns a 3D version of this coordinate.

Attributes Documentation

default_differential#

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

default_representation#

Default representation for position data

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

Mapping for frame-specific component names

name = 'helioprojectiveradial'#
observer#

A frame attribute

No default value

rsun#

A frame attribute

Default: 695700.0 km

theta#

Returns the impact angle, which is the declination angle plus 90 degrees.

Methods Documentation

make_3d()[source]#

Returns a 3D version of this coordinate.

If the coordinate is 2D, the default assumption is that the coordinate is on the surface of the Sun, and the distance component is calculated accordingly. Under this assumption, if the 2D coordinate is outside the disk, the distance component will be NaN.

The assumption can be changed using one of the screens in sunpy.coordinates.screens.

Returns:

HelioprojectiveRadial – The 3D version of this coordinate.