Astronomical Coordinate Systems (astropy.coordinates)

Introduction

The coordinates package provides classes for representing a variety of celestial/spatial coordinates, as well as tools for converting between common coordinate systems in a uniform way.

Getting Started

The simplest way to use coordinates is to use the SkyCoord class. SkyCoord objects are instantiated with a flexible and natural approach that supports inputs provided in a number of convenient formats. The following ways of initializing a coordinate are all equivalent:

>>> from astropy import units as u
>>> from astropy.coordinates import SkyCoord

>>> c = SkyCoord(ra=10.625*u.degree, dec=41.2*u.degree, frame='icrs')
>>> c = SkyCoord(10.625, 41.2, frame='icrs', unit='deg')
>>> c = SkyCoord('00h42m30s', '+41d12m00s', frame='icrs')
>>> c = SkyCoord('00h42.5m', '+41d12m')
>>> c = SkyCoord('00 42 30 +41 12 00', unit=(u.hourangle, u.deg))
>>> c = SkyCoord('00:42.5 +41:12', unit=(u.hourangle, u.deg))
>>> c
<SkyCoord (ICRS): (ra, dec) in deg
    (10.625, 41.2)>

The examples above illustrate a few simple rules to follow when creating a coordinate object:

  • Coordinate values can be provided either as unnamed positional arguments or via keyword arguments like ra, dec, l, or b (depending on the frame).
  • Coordinate frame keyword is optional and defaults to ICRS.
  • Angle units must be specified, either in the values themselves (e.g. 10.5*u.degree or '+41d12m00s') or via the unit keyword.

SkyCoord and all other coordinates objects also support array coordinates. These work the same as single-value coordinates, but they store multiple coordinates in a single object. When you’re going to apply the same operation to many different coordinates (say, from a catalog), this is a better choice than a list of SkyCoord objects, because it will be much faster than applying the operation to each SkyCoord in a for loop.

>>> c = SkyCoord(ra=[10, 11]*u.degree, dec=[41, -5]*u.degree)
>>> c
<SkyCoord (ICRS): (ra, dec) in deg
    [(10.0, 41.0), (11.0, -5.0)]>
>>> c[1]
<SkyCoord (ICRS): (ra, dec) in deg
    (11.0, -5.0)>

Coordinate access

Once you have a coordinate object you can now access the components of that coordinate (e.g. RA, Dec) and get a specific string representation of the full coordinate.

The component values are accessed using aptly named attributes:

>>> c = SkyCoord(ra=10.68458*u.degree, dec=41.26917*u.degree)
>>> c.ra  
<Longitude 10.68458 deg>
>>> c.ra.hour  
0.7123053333333335
>>> c.ra.hms  
hms_tuple(h=0.0, m=42.0, s=44.299200000000525)
>>> c.dec  
<Latitude 41.26917 deg>
>>> c.dec.degree  
41.26917
>>> c.dec.radian  
0.7202828960652683

Coordinates can easily be converted to strings using the to_string() method:

>>> c = SkyCoord(ra=10.68458*u.degree, dec=41.26917*u.degree)
>>> c.to_string('decimal')
'10.6846 41.2692'
>>> c.to_string('dms')
'10d41m04.488s 41d16m09.012s'
>>> c.to_string('hmsdms')
'00h42m44.2992s +41d16m09.012s'

For additional information see the section on Working with Angles.

Transformation

The simplest way to transform to a new coordinate frame is by accessing the appropriately-named attribute. For instance to get the coordinate in the Galactic frame use:

>>> c_icrs = SkyCoord(ra=10.68458*u.degree, dec=41.26917*u.degree, frame='icrs')
>>> c_icrs.galactic  
<SkyCoord (Galactic): (l, b) in deg
    (121.174241811, -21.5728855724)>

For more control, you can use the transform_to method, which accepts a frame name, frame class, or frame instance:

>>> c_fk5 = c_icrs.transform_to('fk5')  # c_icrs.fk5 does the same thing
>>> c_fk5  
<SkyCoord (FK5: equinox=J2000.000): (ra, dec) in deg
    (10.6845915393, 41.2691714591)>

>>> from astropy.coordinates import FK5
>>> c_fk5.transform_to(FK5(equinox='J1975'))  # precess to a different equinox  
<SkyCoord (FK5: equinox=J1975.000): (ra, dec) in deg
    (10.3420913461, 41.1323211229)>

This form of transform_to also makes it straightforward to convert from celestial coordinates to AltAz coordinates, allowing the use of SkyCoord as a tool for planning observations. For a more complete example of this, see Example: Observation Planning.

Representation

So far we have been using a spherical coordinate representation in the all the examples, and this is the default for the built-in frames. Frequently it is convenient to initialize or work with a coordinate using a different representation such as cartesian or cylindrical. This can be done by setting the representation for either SkyCoord objects or low-level frame coordinate objects:

>>> c = SkyCoord(x=1, y=2, z=3, unit='kpc', representation='cartesian')
>>> c
<SkyCoord (ICRS): (x, y, z) in kpc
    (1.0, 2.0, 3.0)>
>>> c.x, c.y, c.z
(<Quantity 1.0 kpc>, <Quantity 2.0 kpc>, <Quantity 3.0 kpc>)

>>> c.representation = 'cylindrical'
>>> c  
<SkyCoord (ICRS): (rho, phi, z) in (kpc, deg, kpc)
    (2.2360679775, 63.4349488229, 3.0)>

For all the details see Representations.

Distance

Distance from the origin (which is system-dependent, but often the Earth center) can also be assigned to a SkyCoord. With two angles and a distance, a unique point in 3D space is available, which also allows conversion to the Cartesian representation of this location:

>>> from astropy.coordinates import Distance
>>> c = SkyCoord(ra=10.68458*u.degree, dec=41.26917*u.degree, distance=770*u.kpc)
>>> c.cartesian.x  
<Quantity 568.7128654235232 kpc>
>>> c.cartesian.y  
<Quantity 107.3008974042025 kpc>
>>> c.cartesian.z  
<Quantity 507.88994291875713 kpc>

With distances assigned, SkyCoord convenience methods are more powerful, as they can make use of the 3D information. For example:

>>> c1 = SkyCoord(ra=10*u.degree, dec=9*u.degree, distance=10*u.pc, frame='icrs')
>>> c2 = SkyCoord(ra=11*u.degree, dec=10*u.degree, distance=11.5*u.pc, frame='icrs')
>>> c1.separation_3d(c2)  
<Distance 1.5228602415117989 pc>

Convenience methods

SkyCoord defines a number of convenience methods as well, like on-sky separation between two coordinates and catalog matching (detailed in Matching Catalogs):

>>> c1 = SkyCoord(ra=10*u.degree, dec=9*u.degree, frame='icrs')
>>> c2 = SkyCoord(ra=11*u.degree, dec=10*u.degree, frame='fk5')
>>> c1.separation(c2)  # Differing frames handled correctly  
<Angle 1.4045335865905868 deg>

The astropy.coordinates subpackage also provides a quick way to get coordinates for named objects assuming you have an active internet connection. The from_name method of SkyCoord uses Sesame to retrieve coordinates for a particular named object:

>>> SkyCoord.from_name("M42")  
<SkyCoord (ICRS): (ra, dec) in deg
    (83.82208, -5.39111)>

For sites (primarily observatories) on the Earth, astropy.coordinates provides a similar quick way to get an EarthLocation:

>>> from astropy.coordinates import EarthLocation
>>> EarthLocation.of_site('Apache Point Observatory')  
<EarthLocation (-1463969.3018517173, -5166673.342234327, 3434985.7120456537) m>

To see the list of site names available, use astropy.coordinates.EarthLocation.get_site_names().

Note

from_name and of_site are for convenience, and hence are by design rather simple. If you need precise coordinates for an object you should find the appropriate reference and input the coordinates manually, or use more specialized functionality like that in the astroquery or astroplan affiliated packages.

Also note that these two methods retrieve data from the internet to determine the celestial or Earth coordinates. The online data may be updated, so if you need to guarantee that your scripts are reproducible in the long term, see the Usage tips/suggestions for methods that access remote resources section.

Overview of astropy.coordinates concepts

Note

The coordinates package from v0.4 onward builds from previous versions of the package, and more detailed information and justification of the design is available in APE (Astropy Proposal for Enhancement) 5.

Here we provide an overview of the package and associated framework. This background information is not necessary for simply using coordinates, particularly if you use the SkyCoord high- level class, but it is helpful for more advanced usage, particularly creating your own frame, transformations, or representations. Another useful piece of background information are some Important Definitions as they are used in coordinates.

coordinates is built on a three-tiered system of objects: representations, frames, and a high-level class. Representations classes are a particular way of storing a three-dimensional data point (or points), such as Cartesian coordinates or spherical polar coordinates. Frames are particular reference frames like FK5 or ICRS, which may store their data in different representations, but have well- defined transformations between each other. These transformations are all stored in the astropy.coordinates.frame_transform_graph, and new transformations can be created by users. Finally, the high-level class (SkyCoord) uses the frame classes, but provides a more accessible interface to these objects as well as various convenience methods and more string-parsing capabilities.

Separating these concepts makes it easier to extend the functionality of coordinates. It allows representations, frames, and transformations to be defined or extended separately, while still preserving the high-level capabilities and simplicity of the SkyCoord class.

Using astropy.coordinates

More detailed information on using the package is provided on separate pages, listed below.

In addition, another resource for the capabilities of this package is the astropy.coordinates.tests.test_api_ape5 testing file. It showcases most of the major capabilities of the package, and hence is a useful supplement to this document. You can see it by either looking at it directly if you downloaded a copy of the astropy source code, or typing the following in an IPython session:

In [1]: from astropy.coordinates.tests import test_api_ape5
In [2]: test_api_ape5??

See Also

Some references particularly useful in understanding subtleties of the coordinate systems implemented here include:

  • USNO Circular 179

    A useful guide to the IAU 2000/2003 work surrounding ICRS/IERS/CIRS and related problems in precision coordinate system work.

  • Standards Of Fundamental Astronomy

    The definitive implementation of IAU-defined algorithms. The “SOFA Tools for Earth Attitude” document is particularly valuable for understanding the latest IAU standards in detail.

  • IERS Conventions (2010)

    An exhaustive reference covering the ITRS, the IAU2000 celestial coordinates framework, and other related details of modern coordinate conventions.

  • Meeus, J. “Astronomical Algorithms”

    A valuable text describing details of a wide range of coordinate-related problems and concepts.

Reference/API

astropy.coordinates Package

This subpackage contains classes and functions for celestial coordinates of astronomical objects. It also contains a framework for conversions between coordinate systems.

The diagram below shows all of the coordinate systems built into the coordinates package, their aliases (useful for converting other coordinates to them using attribute-style access) and the pre-defined transformations between them. The user is free to override any of these transformations by defining new transformations between these systems, but the pre-defined transformations should be sufficient for typical usage.

The graph also indicates the priority for each transformation as a number next to the arrow. These priorities are used to decide the preferred order when two transformation paths have the same number of steps. These priorities are defined such that the path with a smaller total priority is favored.

digraph AstropyCoordinateTransformGraph { ICRS [shape=oval label="ICRS\n`icrs`"]; HeliocentricTrueEcliptic [shape=oval label="HeliocentricTrueEcliptic\n`heliocentrictrueecliptic`"]; GCRS [shape=oval label="GCRS\n`gcrs`"]; Galactocentric [shape=oval label="Galactocentric\n`galactocentric`"]; FK5 [shape=oval label="FK5\n`fk5`"]; BarycentricTrueEcliptic [shape=oval label="BarycentricTrueEcliptic\n`barycentrictrueecliptic`"]; CIRS [shape=oval label="CIRS\n`cirs`"]; PrecessedGeocentric [shape=oval label="PrecessedGeocentric\n`precessedgeocentric`"]; GeocentricTrueEcliptic [shape=oval label="GeocentricTrueEcliptic\n`geocentrictrueecliptic`"]; FK4NoETerms [shape=oval label="FK4NoETerms\n`fk4noeterms`"]; FK4 [shape=oval label="FK4\n`fk4`"]; Galactic [shape=oval label="Galactic\n`galactic`"]; Supergalactic [shape=oval label="Supergalactic\n`supergalactic`"]; ITRS [shape=oval label="ITRS\n`itrs`"]; AltAz [shape=oval label="AltAz\n`altaz`"]; ICRS -> HeliocentricTrueEcliptic[ label = "1.0" ]; ICRS -> GCRS[ label = "1.0" ]; ICRS -> Galactocentric[ label = "1.0" ]; ICRS -> FK5[ label = "1.0" ]; ICRS -> BarycentricTrueEcliptic[ label = "1.0" ]; ICRS -> CIRS[ label = "1.0" ]; HeliocentricTrueEcliptic -> ICRS[ label = "1.0" ]; GCRS -> ICRS[ label = "1.0" ]; GCRS -> CIRS[ label = "1.0" ]; GCRS -> PrecessedGeocentric[ label = "1.0" ]; GCRS -> GeocentricTrueEcliptic[ label = "1.0" ]; GCRS -> GCRS[ label = "1.0" ]; FK4NoETerms -> FK5[ label = "1.0" ]; FK4NoETerms -> FK4NoETerms[ label = "1.0" ]; FK4NoETerms -> FK4[ label = "1.0" ]; FK4NoETerms -> Galactic[ label = "1.0" ]; Galactocentric -> ICRS[ label = "1.0" ]; Galactic -> FK5[ label = "1.0" ]; Galactic -> FK4NoETerms[ label = "1.0" ]; Galactic -> Supergalactic[ label = "1.0" ]; FK5 -> FK5[ label = "1.0" ]; FK5 -> ICRS[ label = "1.0" ]; FK5 -> FK4NoETerms[ label = "1.0" ]; FK5 -> Galactic[ label = "1.0" ]; BarycentricTrueEcliptic -> ICRS[ label = "1.0" ]; Supergalactic -> Galactic[ label = "1.0" ]; CIRS -> ICRS[ label = "1.0" ]; CIRS -> CIRS[ label = "1.0" ]; CIRS -> ITRS[ label = "1.0" ]; CIRS -> AltAz[ label = "1.0" ]; CIRS -> GCRS[ label = "1.0" ]; PrecessedGeocentric -> GCRS[ label = "1.0" ]; GeocentricTrueEcliptic -> GCRS[ label = "1.0" ]; FK4 -> FK4[ label = "1.0" ]; FK4 -> FK4NoETerms[ label = "1.0" ]; ITRS -> CIRS[ label = "1.0" ]; ITRS -> ITRS[ label = "1.0" ]; AltAz -> CIRS[ label = "1.0" ]; AltAz -> AltAz[ label = "1.0" ]; overlap=false }

Note

The ecliptic coordinate systems (added in Astropy v1.1) have not been extensively tested for accuracy or consistency with other implementations of ecliptic coordinates. We welcome contributions to add such testing, but in the meantime, users who depend on consistency with other implementations may wish to check test inputs against good datasets before using Astropy’s ecliptic coordinates.

Functions

cartesian_to_spherical(x, y, z) Converts 3D rectangular cartesian coordinates to spherical polar coordinates.
concatenate(coords) Combine multiple coordinate objects into a single SkyCoord.
get_constellation(coord[, short_name, ...]) Determines the constellation(s) a given coordinate object contains.
get_icrs_coordinates(name) Retrieve an ICRS object by using an online name resolving service to retrieve coordinates for the specified name.
get_sun(time) Determines the location of the sun at a given time (or times, if the input is an array Time object), in geocentric coordinates.
match_coordinates_3d(matchcoord, catalogcoord) Finds the nearest 3-dimensional matches of a coordinate or coordinates in a set of catalog coordinates.
match_coordinates_sky(matchcoord, catalogcoord) Finds the nearest on-sky matches of a coordinate or coordinates in a set of catalog coordinates.
search_around_3d(coords1, coords2, distlimit) Searches for pairs of points that are at least as close as a specified distance in 3D space.
search_around_sky(coords1, coords2, seplimit) Searches for pairs of points that have an angular separation at least as close as a specified angle.
spherical_to_cartesian(r, lat, lon) Converts spherical polar coordinates to rectangular cartesian coordinates.

Classes

AltAz(*args, **kwargs) A coordinate or frame in the Altitude-Azimuth system (Horizontal coordinates).
Angle One or more angular value(s) with units equivalent to radians or degrees.
BarycentricTrueEcliptic(*args, **kwargs) Barycentric ecliptic coordinates.
BaseCoordinateFrame(*args, **kwargs) The base class for coordinate frames.
BaseRepresentation Base Representation object, for representing a point in a 3D coordinate system.
BoundsError Raised when an angle is outside of its user-specified bounds.
CIRS(*args, **kwargs) A coordinate or frame in the Celestial Intermediate Reference System (CIRS).
CartesianRepresentation(x[, y, z, copy]) Representation of points in 3D cartesian coordinates.
CompositeTransform(transforms, fromsys, tosys) A transformation constructed by combining together a series of single-step transformations.
ConvertError Raised if a coordinate system cannot be converted to another
CoordinateTransform(fromsys, tosys[, ...]) An object that transforms a coordinate from one system to another.
CylindricalRepresentation(rho, phi, z[, copy]) Representation of points in 3D cylindrical coordinates.
Distance A one-dimensional distance.
DynamicMatrixTransform(matrix_func, fromsys, ...) A coordinate transformation specified as a function that yields a 3 x 3 cartesian transformation matrix.
EarthLocation Location on the Earth.
EarthLocationAttribute([default, ...]) A frame attribute that can act as a EarthLocation.
FK4(*args, **kwargs) A coordinate or frame in the FK4 system.
FK4NoETerms(*args, **kwargs) A coordinate or frame in the FK4 system, but with the E-terms of aberration removed.
FK5(*args, **kwargs) A coordinate or frame in the FK5 system.
FrameAttribute([default, secondary_attribute]) A non-mutable data descriptor to hold a frame attribute.
FunctionTransform(func, fromsys, tosys[, ...]) A coordinate transformation defined by a function that accepts a coordinate object and returns the transformed coordinate object.
GCRS(*args, **kwargs) A coordinate or frame in the Geocentric Celestial Reference System (GCRS).
Galactic(*args, **kwargs) Galactic Coordinates.
Galactocentric(*args, **kwargs) A coordinate or frame in the Galactocentric system.
GenericFrame(frame_attrs) A frame object that can’t store data but can hold any arbitrary frame attributes.
GeocentricTrueEcliptic(*args, **kwargs) Geocentric ecliptic coordinates.
HeliocentricTrueEcliptic(*args, **kwargs) Heliocentric ecliptic coordinates.
ICRS(*args, **kwargs) A coordinate or frame in the ICRS system.
ITRS(*args, **kwargs) A coordinate or frame in the International Terrestrial Reference System (ITRS).
IllegalHourError(hour) Raised when an hour value is not in the range [0,24).
IllegalHourWarning(hour[, alternativeactionstr]) Raised when an hour value is 24.
IllegalMinuteError(minute) Raised when an minute value is not in the range [0,60].
IllegalMinuteWarning(minute[, ...]) Raised when a minute value is 60.
IllegalSecondError(second) Raised when an second value (time) is not in the range [0,60].
IllegalSecondWarning(second[, ...]) Raised when a second value is 60.
Latitude Latitude-like angle(s) which must be in the range -90 to +90 deg.
Longitude Longitude-like angle(s) which are wrapped within a contiguous 360 degree range.
PhysicsSphericalRepresentation(phi, theta, r) Representation of points in 3D spherical coordinates (using the physics convention of using phi and theta for azimuth and inclination from the pole).
PrecessedGeocentric(*args, **kwargs) A coordinate frame defined in a similar manner as GCRS, but precessed to a requested (mean) equinox.
QuantityFrameAttribute([default, ...]) A frame attribute that is a quantity with specified units and shape (optionally).
RangeError Raised when some part of an angle is out of its valid range.
RepresentationMapping This namedtuple is used with the frame_specific_representation_info attribute to tell frames what attribute names (and default units) to use for a particular representation.
SkyCoord(*args, **kwargs) High-level object providing a flexible interface for celestial coordinate representation, manipulation, and transformation between systems.
SphericalRepresentation(lon, lat, distance) Representation of points in 3D spherical coordinates.
StaticMatrixTransform(matrix, fromsys, tosys) A coordinate transformation defined as a 3 x 3 cartesian transformation matrix.
Supergalactic(*args, **kwargs) Supergalactic Coordinates (see Lahav et al.
TimeFrameAttribute([default, ...]) Frame attribute descriptor for quantities that are Time objects.
TransformGraph() A graph representing the paths between coordinate frames.
UnitSphericalRepresentation(lon, lat[, copy]) Representation of points on a unit sphere.
UnknownSiteException(site, attribute[, ...])

Class Inheritance Diagram

Inheritance diagram of astropy.coordinates.builtin_frames.altaz.AltAz, astropy.coordinates.angles.Angle, astropy.coordinates.builtin_frames.ecliptic.BarycentricTrueEcliptic, astropy.coordinates.baseframe.BaseCoordinateFrame, astropy.coordinates.representation.BaseRepresentation, astropy.coordinates.errors.BoundsError, astropy.coordinates.builtin_frames.cirs.CIRS, astropy.coordinates.representation.CartesianRepresentation, astropy.coordinates.transformations.CompositeTransform, astropy.coordinates.errors.ConvertError, astropy.coordinates.transformations.CoordinateTransform, astropy.coordinates.representation.CylindricalRepresentation, astropy.coordinates.distances.Distance, astropy.coordinates.transformations.DynamicMatrixTransform, astropy.coordinates.earth.EarthLocation, astropy.coordinates.baseframe.EarthLocationAttribute, astropy.coordinates.builtin_frames.fk4.FK4, astropy.coordinates.builtin_frames.fk4.FK4NoETerms, astropy.coordinates.builtin_frames.fk5.FK5, astropy.coordinates.baseframe.FrameAttribute, astropy.coordinates.transformations.FunctionTransform, astropy.coordinates.builtin_frames.gcrs.GCRS, astropy.coordinates.builtin_frames.galactic.Galactic, astropy.coordinates.builtin_frames.galactocentric.Galactocentric, astropy.coordinates.baseframe.GenericFrame, astropy.coordinates.builtin_frames.ecliptic.GeocentricTrueEcliptic, astropy.coordinates.builtin_frames.ecliptic.HeliocentricTrueEcliptic, astropy.coordinates.builtin_frames.icrs.ICRS, astropy.coordinates.builtin_frames.itrs.ITRS, astropy.coordinates.errors.IllegalHourError, astropy.coordinates.errors.IllegalHourWarning, astropy.coordinates.errors.IllegalMinuteError, astropy.coordinates.errors.IllegalMinuteWarning, astropy.coordinates.errors.IllegalSecondError, astropy.coordinates.errors.IllegalSecondWarning, astropy.coordinates.angles.Latitude, astropy.coordinates.angles.Longitude, astropy.coordinates.representation.PhysicsSphericalRepresentation, astropy.coordinates.builtin_frames.gcrs.PrecessedGeocentric, astropy.coordinates.baseframe.QuantityFrameAttribute, astropy.coordinates.errors.RangeError, astropy.coordinates.baseframe.RepresentationMapping, astropy.coordinates.sky_coordinate.SkyCoord, astropy.coordinates.representation.SphericalRepresentation, astropy.coordinates.transformations.StaticMatrixTransform, astropy.coordinates.builtin_frames.supergalactic.Supergalactic, astropy.coordinates.baseframe.TimeFrameAttribute, astropy.coordinates.transformations.TransformGraph, astropy.coordinates.representation.UnitSphericalRepresentation, astropy.coordinates.errors.UnknownSiteException