Astronomical Coordinate Systems (astropy.coordinates)

Introduction

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

Getting Started

The best way to start using coordinates is to use the SkyCoord class. SkyCoord objects are instantiated by passing in positions (and optional velocities) with specified units and a coordinate frame. Sky positions are commonly passed in as Quantity objects and the frame is specified with the string name. As an example of creating a SkyCoord to represent an ICRS (Right ascension [RA], Declination [Dec]) sky position:

>>> 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')

The initializer for SkyCoord is very flexible and supports inputs provided in a number of convenient formats. The following ways of initializing a coordinate are all equivalent to the above:

>>> 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 rules to follow when creating a coordinate object:

  • Coordinate values can be provided either as unnamed positional arguments or via keyword arguments like ra and dec, or l and b (depending on the frame).

  • The coordinate frame keyword is optional because it defaults to ICRS.

  • Angle units must be specified for all components, either by passing in a Quantity object (e.g., 10.5*u.degree), by including them in the value (e.g., '+41d12m00s'), or via the unit keyword.

SkyCoord and all other coordinates objects also support array coordinates. These work in the same way as single-value coordinates, but they store multiple coordinates in a single object. When you are 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. Like the underlying ndarray instances that contain the data, SkyCoord objects can be sliced, reshaped, etc.:

>>> c = SkyCoord(ra=[10, 11, 12, 13]*u.degree, dec=[41, -5, 42, 0]*u.degree)
>>> c  
<SkyCoord (ICRS): (ra, dec) in deg
    [(10., 41.), (11., -5.), (12., 42.), (13.,  0.)]>
>>> c[1]  
<SkyCoord (ICRS): (ra, dec) in deg
    (11., -5.)>
>>> c.reshape(2, 2)  
<SkyCoord (ICRS): (ra, dec) in deg
    [[(10., 41.), (11., -5.)],
     [(12., 42.), (13.,  0.)]]>

Coordinate Access

Once you have a coordinate object you can access the components of that coordinate (e.g., RA, Dec) to get string representations of the full coordinate.

The component values are accessed using (typically lowercase) named attributes that depend on the coordinate frame (e.g., ICRS, Galactic, etc.). For the default, ICRS, the coordinate component names are ra and dec:

>>> 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 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

One convenient 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.17424181, -21.57288557)>

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.68459154, 41.26917146)>

>>> 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.34209135, 41.13232112)>

This form of transform_to also makes it possible 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 Determining and plotting the altitude/azimuth of a celestial object.

Some coordinate frames such as AltAz require Earth rotation information (UT1-UTC offset and/or polar motion) when transforming to/from other frames. These Earth rotation values are automatically downloaded from the International Earth Rotation and Reference Systems (IERS) service when required. See IERS data access (astropy.utils.iers) for details of this process.

Representation

So far we have been using a spherical coordinate representation in all of our 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_type for either SkyCoord objects or low-level frame coordinate objects:

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

>>> c.representation_type = 'cylindrical'
>>> c  
<SkyCoord (ICRS): (rho, phi, z) in (kpc, deg, kpc)
    (2.23606798, 63.43494882, 3.)>

For all of the details see Representations.

Distance

SkyCoord and the individual frame classes also support specifying a distance from the frame origin. The origin depends on the particular coordinate frame; this can be, for example, centered on the earth, centered on the solar system barycenter, etc. Two angles and a distance specify a unique point in 3D space, which also allows converting the coordinates to a Cartesian representation:

>>> c = SkyCoord(ra=10.68458*u.degree, dec=41.26917*u.degree, distance=770*u.kpc)
>>> c.cartesian.x  
<Quantity 568.71286542 kpc>
>>> c.cartesian.y  
<Quantity 107.3008974 kpc>
>>> c.cartesian.z  
<Quantity 507.88994292 kpc>

With distances assigned, SkyCoord convenience methods are more powerful, as they can make use of the 3D information. For example, to compute the physical, 3D separation between two points in space:

>>> 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.52286024 pc>

Convenience Methods

SkyCoord defines a number of convenience methods that support, for example, computing on-sky (i.e., angular) and 3D separations between two coordinates:

>>> 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.40453359 deg>

Or cross-matching catalog coordinates (detailed in Matching Catalogs):

>>> target_c = SkyCoord(ra=10*u.degree, dec=9*u.degree, frame='icrs')
>>> # read in coordinates from a catalog...
>>> catalog_c = ... 
>>> idx, sep, _ = target_c.match_to_catalog_sky(catalog_c) 

The astropy.coordinates sub-package 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("PSR J1012+5307")  
<SkyCoord (ICRS): (ra, dec) in deg
    (153.1393271, 53.117343)>

In some cases, the coordinates are embedded in the catalog name of the object. For such object names, from_name is able to parse the coordinates from the name if given the parse=True option. For slow connections, this may be much faster than a sesame query for the same object name. It’s worth noting, however, that the coordinates extracted in this way may differ from the database coordinates by a few deci-arcseconds, so only use this option if you do not need sub-arcsecond accuracy for your coordinates:

>>> SkyCoord.from_name("CRTS SSS100805 J194428-420209", parse=True)  
<SkyCoord (ICRS): (ra, dec) in deg
    (296.11666667, -42.03583333)>

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

>>> from astropy.coordinates import EarthLocation
>>> EarthLocation.of_site('Apache Point Observatory')  
<EarthLocation (-1463969.30185172, -5166673.34223433,  3434985.71204565) m>

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

For arbitrary Earth addresses (e.g., not observatory sites), use the of_address classmethod. Any address passed to this function uses Google maps to retrieve the latitude and longitude and can also (optionally) query Google maps to get the height of the location. As with Google maps, this works with fully specified addresses, location names, city names, etc.:

>>> EarthLocation.of_address('1002 Holy Grail Court, St. Louis, MO')
<EarthLocation (-26726.98216371, -4997009.8604809, 3950271.16507911) m>
>>> EarthLocation.of_address('1002 Holy Grail Court, St. Louis, MO',
...                          get_height=True)
<EarthLocation (-26727.6272786, -4997130.47437768, 3950367.15622108) m>
>>> EarthLocation.of_address('Danbury, CT')
<EarthLocation ( 1364606.64511651, -4593292.9428273,  4195415.93695139) m>

Note

from_name, of_site, and of_address are for convenience, and hence are by design relatively low precision. If you need more 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 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.

This functionality can be combined to do more complicated tasks like computing barycentric corrections to radial velocity observations (also a supported high-level SkyCoord method - see Radial Velocity Corrections):

>>> from astropy.time import Time
>>> obstime = Time('2017-2-14')
>>> target = SkyCoord.from_name('M31')  
>>> keck = EarthLocation.of_site('Keck')  
>>> target.radial_velocity_correction(obstime=obstime, location=keck).to('km/s')  
<Quantity -22.359784554780255 km / s>

Velocities (Proper Motions and Radial Velocities)

In addition to positional coordinates, coordinates supports storing and transforming velocities. These are available both via the lower-level coordinate frame classes, and (new in v3.0) via SkyCoord objects:

>>> sc = SkyCoord(1*u.deg, 2*u.deg, radial_velocity=20*u.km/u.s)
>>> sc  
<SkyCoord (ICRS): (ra, dec) in deg
    ( 1.,  2.)
 (radial_velocity) in km / s
    ( 20.,)>

For more details on velocity support (and limitations), see the Working with Velocities in Astropy Coordinates page.

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 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 ease-of-use of the SkyCoord class.

Examples:

See Determining and plotting the altitude/azimuth of a celestial object for an example of using the coordinates functionality to prepare for an observing run.

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 downloading 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??

Performance Tips

If you are using SkyCoord for many different coordinates, you will see much better performance if you create a single SkyCoord with arrays of coordinates as opposed to creating individual SkyCoord objects for each individual coordinate:

>>> coord = SkyCoord(ra_array, dec_array, unit='deg')  

In addition, looping over a SkyCoord object can be slow. If you need to transform the coordinates to a different frame, it is much faster to transform a single SkyCoord with arrays of values as opposed to looping over the SkyCoord and transforming them individually.

Finally, for more advanced users, note that you can use broadcasting to transform SkyCoord objects into frames with vector properties. For example:

>>> from astropy.coordinates import SkyCoord, EarthLocation
>>> from astropy import coordinates as coord
>>> from astropy.coordinates.tests.utils import randomly_sample_sphere
>>> from astropy.time import Time
>>> from astropy import units as u
>>> import numpy as np

>>> # 1000 random locations on the sky
>>> ra, dec, _ = randomly_sample_sphere(1000)
>>> coos = SkyCoord(ra, dec)

>>> # 300 times over the space of 10 hours
>>> times = Time.now() + np.linspace(-5, 5, 300)*u.hour

>>> # note the use of broadcasting so that 300 times are broadcast against 1000 positions
>>> lapalma = EarthLocation.from_geocentric(5327448.9957829, -1718665.73869569, 3051566.90295403, unit='m')
>>> aa_frame = coord.AltAz(obstime=times[:, np.newaxis], location=lapalma)

>>> # calculate alt-az of each object at each time.
>>> aa_coos = coos.transform_to(aa_frame)  

See Also

Some references that are 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.

Built-in Frame Classes

astropy.coordinates.builtin_frames Package

The diagram below shows all of the built in coordinate systems, 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 color of an edge in the graph (i.e. the transformations between two frames) is set by the type of transformation; the legend box defines the mapping from transform class name to color.

digraph AstropyCoordinateTransformGraph { FK5 [shape=oval label="FK5\n`fk5`"]; ICRS [shape=oval label="ICRS\n`icrs`"]; FK4NoETerms [shape=oval label="FK4NoETerms\n`fk4noeterms`"]; Galactic [shape=oval label="Galactic\n`galactic`"]; FK4 [shape=oval label="FK4\n`fk4`"]; Galactocentric [shape=oval label="Galactocentric\n`galactocentric`"]; LSR [shape=oval label="LSR\n`lsr`"]; CIRS [shape=oval label="CIRS\n`cirs`"]; GCRS [shape=oval label="GCRS\n`gcrs`"]; HCRS [shape=oval label="HCRS\n`hcrs`"]; BarycentricMeanEcliptic [shape=oval label="BarycentricMeanEcliptic\n`barycentricmeanecliptic`"]; HeliocentricMeanEcliptic [shape=oval label="HeliocentricMeanEcliptic\n`heliocentricmeanecliptic`"]; BarycentricTrueEcliptic [shape=oval label="BarycentricTrueEcliptic\n`barycentrictrueecliptic`"]; HeliocentricTrueEcliptic [shape=oval label="HeliocentricTrueEcliptic\n`heliocentrictrueecliptic`"]; HeliocentricEclipticIAU76 [shape=oval label="HeliocentricEclipticIAU76\n`heliocentriceclipticiau76`"]; CustomBarycentricEcliptic [shape=oval label="CustomBarycentricEcliptic\n`custombarycentricecliptic`"]; GalacticLSR [shape=oval label="GalacticLSR\n`galacticlsr`"]; Supergalactic [shape=oval label="Supergalactic\n`supergalactic`"]; AltAz [shape=oval label="AltAz\n`altaz`"]; ITRS [shape=oval label="ITRS\n`itrs`"]; PrecessedGeocentric [shape=oval label="PrecessedGeocentric\n`precessedgeocentric`"]; GeocentricMeanEcliptic [shape=oval label="GeocentricMeanEcliptic\n`geocentricmeanecliptic`"]; GeocentricTrueEcliptic [shape=oval label="GeocentricTrueEcliptic\n`geocentrictrueecliptic`"]; FK5 -> FK5[ color = "#1b9e77" ]; FK5 -> ICRS[ color = "#1b9e77" ]; FK5 -> FK4NoETerms[ color = "#1b9e77" ]; FK5 -> Galactic[ color = "#1b9e77" ]; FK4 -> FK4[ color = "#d95f02" ]; FK4 -> FK4NoETerms[ color = "#d95f02" ]; FK4NoETerms -> FK4NoETerms[ color = "#1b9e77" ]; FK4NoETerms -> FK4[ color = "#d95f02" ]; FK4NoETerms -> FK5[ color = "#1b9e77" ]; FK4NoETerms -> Galactic[ color = "#1b9e77" ]; ICRS -> Galactocentric[ color = "#555555" ]; ICRS -> LSR[ color = "#555555" ]; ICRS -> FK5[ color = "#1b9e77" ]; ICRS -> CIRS[ color = "#d95f02" ]; ICRS -> GCRS[ color = "#d95f02" ]; ICRS -> HCRS[ color = "#555555" ]; ICRS -> BarycentricMeanEcliptic[ color = "#1b9e77" ]; ICRS -> HeliocentricMeanEcliptic[ color = "#555555" ]; ICRS -> BarycentricTrueEcliptic[ color = "#1b9e77" ]; ICRS -> HeliocentricTrueEcliptic[ color = "#555555" ]; ICRS -> HeliocentricEclipticIAU76[ color = "#555555" ]; ICRS -> CustomBarycentricEcliptic[ color = "#1b9e77" ]; Galactocentric -> ICRS[ color = "#555555" ]; LSR -> ICRS[ color = "#555555" ]; Galactic -> GalacticLSR[ color = "#555555" ]; Galactic -> FK5[ color = "#1b9e77" ]; Galactic -> FK4NoETerms[ color = "#1b9e77" ]; Galactic -> Supergalactic[ color = "#7570b3" ]; GalacticLSR -> Galactic[ color = "#555555" ]; Supergalactic -> Galactic[ color = "#7570b3" ]; CIRS -> ICRS[ color = "#d95f02" ]; CIRS -> CIRS[ color = "#d95f02" ]; CIRS -> AltAz[ color = "#d95f02" ]; CIRS -> GCRS[ color = "#d95f02" ]; CIRS -> ITRS[ color = "#d95f02" ]; GCRS -> ICRS[ color = "#d95f02" ]; GCRS -> GCRS[ color = "#d95f02" ]; GCRS -> HCRS[ color = "#d95f02" ]; GCRS -> CIRS[ color = "#d95f02" ]; GCRS -> PrecessedGeocentric[ color = "#d95f02" ]; GCRS -> GeocentricMeanEcliptic[ color = "#d95f02" ]; GCRS -> GeocentricTrueEcliptic[ color = "#d95f02" ]; HCRS -> ICRS[ color = "#555555" ]; HCRS -> HCRS[ color = "#d95f02" ]; AltAz -> CIRS[ color = "#d95f02" ]; AltAz -> AltAz[ color = "#d95f02" ]; ITRS -> CIRS[ color = "#d95f02" ]; ITRS -> ITRS[ color = "#d95f02" ]; PrecessedGeocentric -> GCRS[ color = "#d95f02" ]; GeocentricMeanEcliptic -> GCRS[ color = "#d95f02" ]; BarycentricMeanEcliptic -> ICRS[ color = "#1b9e77" ]; HeliocentricMeanEcliptic -> ICRS[ color = "#555555" ]; GeocentricTrueEcliptic -> GCRS[ color = "#d95f02" ]; BarycentricTrueEcliptic -> ICRS[ color = "#1b9e77" ]; HeliocentricTrueEcliptic -> ICRS[ color = "#555555" ]; HeliocentricEclipticIAU76 -> ICRS[ color = "#555555" ]; CustomBarycentricEcliptic -> ICRS[ color = "#1b9e77" ]; overlap=false }

  • AffineTransform:

  • FunctionTransform:

  • FunctionTransformWithFiniteDifference:

  • StaticMatrixTransform:

  • DynamicMatrixTransform:

Classes

ICRS(*args[, copy, representation_type, …])

A coordinate or frame in the ICRS system.

FK5(*args[, copy, representation_type, …])

A coordinate or frame in the FK5 system.

FK4(*args[, copy, representation_type, …])

A coordinate or frame in the FK4 system.

FK4NoETerms(*args[, copy, …])

A coordinate or frame in the FK4 system, but with the E-terms of aberration removed.

Galactic(*args[, copy, representation_type, …])

A coordinate or frame in the Galactic coordinate system.

Galactocentric(*args, **kwargs)

A coordinate or frame in the Galactocentric system.

Supergalactic(*args[, copy, …])

Supergalactic Coordinates (see Lahav et al.

AltAz(*args, **kwargs)

A coordinate or frame in the Altitude-Azimuth system (Horizontal coordinates).

GCRS(*args[, copy, representation_type, …])

A coordinate or frame in the Geocentric Celestial Reference System (GCRS).

CIRS(*args[, copy, representation_type, …])

A coordinate or frame in the Celestial Intermediate Reference System (CIRS).

ITRS(*args[, copy, representation_type, …])

A coordinate or frame in the International Terrestrial Reference System (ITRS).

HCRS(*args[, copy, representation_type, …])

A coordinate or frame in a Heliocentric system, with axes aligned to ICRS.

PrecessedGeocentric(*args[, copy, …])

A coordinate frame defined in a similar manner as GCRS, but precessed to a requested (mean) equinox.

GeocentricMeanEcliptic(*args[, copy, …])

Geocentric mean ecliptic coordinates.

BarycentricMeanEcliptic(*args[, copy, …])

Barycentric mean ecliptic coordinates.

HeliocentricMeanEcliptic(*args[, copy, …])

Heliocentric mean ecliptic coordinates.

GeocentricTrueEcliptic(*args[, copy, …])

Geocentric true ecliptic coordinates.

BarycentricTrueEcliptic(*args[, copy, …])

Barycentric true ecliptic coordinates.

HeliocentricTrueEcliptic(*args[, copy, …])

Heliocentric true ecliptic coordinates.

SkyOffsetFrame(*args, **kwargs)

A frame which is relative to some specific position and oriented to match its frame.

GalacticLSR(*args[, copy, …])

A coordinate or frame in the Local Standard of Rest (LSR), axis-aligned to the Galactic frame.

LSR(*args[, copy, representation_type, …])

A coordinate or frame in the Local Standard of Rest (LSR).

BaseEclipticFrame(*args[, copy, …])

A base class for frames that have names and conventions like that of ecliptic frames.

BaseRADecFrame(*args[, copy, …])

A base class that defines default representation info for frames that represent longitude and latitude as Right Ascension and Declination following typical “equatorial” conventions.

HeliocentricEclipticIAU76(*args[, copy, …])

Heliocentric mean (IAU 1976) ecliptic coordinates.

CustomBarycentricEcliptic(*args[, copy, …])

Barycentric ecliptic coordinates with custom obliquity.

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.

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.

concatenate_representations(reps)

Combine multiple representation objects into a single instance by concatenating the data in each component.

get_body(body, time[, location, ephemeris])

Get a SkyCoord for a solar system body as observed from a location on Earth in the GCRS reference system.

get_body_barycentric(body, time[, ephemeris])

Calculate the barycentric position of a solar system body.

get_body_barycentric_posvel(body, time[, …])

Calculate the barycentric position and velocity of a solar system body.

get_constellation(coord[, short_name, …])

Determines the constellation(s) a given coordinate object contains.

get_icrs_coordinates(name[, parse])

Retrieve an ICRS object by using an online name resolving service to retrieve coordinates for the specified name.

get_moon(time[, location, ephemeris])

Get a SkyCoord for the Earth’s Moon as observed from a location on Earth in the GCRS reference system.

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.

make_transform_graph_docs(transform_graph)

Generates a string that can be used in other docstrings to include a transformation graph, showing the available transforms and coordinate systems.

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

AffineTransform(transform_func, fromsys, tosys)

A coordinate transformation specified as a function that yields a 3 x 3 cartesian transformation matrix and a tuple of displacement vectors.

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.

Attribute([default, secondary_attribute])

A non-mutable data descriptor to hold a frame attribute.

BarycentricMeanEcliptic(*args[, copy, …])

Barycentric mean ecliptic coordinates.

BarycentricTrueEcliptic(*args[, copy, …])

Barycentric true ecliptic coordinates.

BaseAffineTransform(fromsys, tosys[, …])

Base class for common functionality between the AffineTransform-type subclasses.

BaseCoordinateFrame(*args[, copy, …])

The base class for coordinate frames.

BaseDifferential(*args, **kwargs)

A base class representing differentials of representations.

BaseEclipticFrame(*args[, copy, …])

A base class for frames that have names and conventions like that of ecliptic frames.

BaseRADecFrame(*args[, copy, …])

A base class that defines default representation info for frames that represent longitude and latitude as Right Ascension and Declination following typical “equatorial” conventions.

BaseRepresentation(*args[, differentials])

Base for representing a point in a 3D coordinate system.

BaseRepresentationOrDifferential(*args, **kwargs)

3D coordinate representations and differentials.

BaseSphericalCosLatDifferential(*args, **kwargs)

Differentials from points on a spherical base representation.

BaseSphericalDifferential(*args, **kwargs)

BoundsError

Raised when an angle is outside of its user-specified bounds.

CIRS(*args[, copy, representation_type, …])

A coordinate or frame in the Celestial Intermediate Reference System (CIRS).

CartesianDifferential(d_x[, d_y, d_z, unit, …])

Differentials in of points in 3D cartesian coordinates.

CartesianRepresentation(x[, y, z, unit, …])

Representation of points in 3D cartesian coordinates.

CartesianRepresentationAttribute([default, …])

A frame attribute that is a CartesianRepresentation with specified units.

CartesianRepresentationFrameAttribute(*args, …)

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

CoordinateAttribute(frame[, default, …])

A frame attribute which is a coordinate object.

CoordinateTransform(fromsys, tosys[, …])

An object that transforms a coordinate from one system to another.

CustomBarycentricEcliptic(*args[, copy, …])

Barycentric ecliptic coordinates with custom obliquity.

CylindricalDifferential(d_rho, d_phi, d_z[, …])

Differential(s) of points in cylindrical coordinates.

CylindricalRepresentation(rho, phi, z[, …])

Representation of points in 3D cylindrical coordinates.

DifferentialAttribute([default, …])

A frame attribute which is a differential instance.

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[, copy, representation_type, …])

A coordinate or frame in the FK4 system.

FK4NoETerms(*args[, copy, …])

A coordinate or frame in the FK4 system, but with the E-terms of aberration removed.

FK5(*args[, copy, representation_type, …])

A coordinate or frame in the FK5 system.

FunctionTransform(func, fromsys, tosys[, …])

A coordinate transformation defined by a function that accepts a coordinate object and returns the transformed coordinate object.

FunctionTransformWithFiniteDifference(func, …)

A coordinate transformation that works like a FunctionTransform, but computes velocity shifts based on the finite-difference relative to one of the frame attributes.

GCRS(*args[, copy, representation_type, …])

A coordinate or frame in the Geocentric Celestial Reference System (GCRS).

Galactic(*args[, copy, representation_type, …])

A coordinate or frame in the Galactic coordinate system.

GalacticLSR(*args[, copy, …])

A coordinate or frame in the Local Standard of Rest (LSR), axis-aligned to the Galactic frame.

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.

GeocentricMeanEcliptic(*args[, copy, …])

Geocentric mean ecliptic coordinates.

GeocentricTrueEcliptic(*args[, copy, …])

Geocentric true ecliptic coordinates.

HCRS(*args[, copy, representation_type, …])

A coordinate or frame in a Heliocentric system, with axes aligned to ICRS.

HeliocentricEclipticIAU76(*args[, copy, …])

Heliocentric mean (IAU 1976) ecliptic coordinates.

HeliocentricMeanEcliptic(*args[, copy, …])

Heliocentric mean ecliptic coordinates.

HeliocentricTrueEcliptic(*args[, copy, …])

Heliocentric true ecliptic coordinates.

ICRS(*args[, copy, representation_type, …])

A coordinate or frame in the ICRS system.

ITRS(*args[, copy, representation_type, …])

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.

LSR(*args[, copy, representation_type, …])

A coordinate or frame in the Local Standard of Rest (LSR).

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.

PhysicsSphericalDifferential(d_phi, d_theta, d_r)

Differential(s) of 3D spherical coordinates using physics convention.

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[, copy, …])

A coordinate frame defined in a similar manner as GCRS, but precessed to a requested (mean) equinox.

QuantityAttribute(default[, …])

A frame attribute that is a quantity with specified units and shape (optionally).

QuantityFrameAttribute(*args, **kwargs)

RadialDifferential(*args, **kwargs)

Differential(s) of radial distances.

RadialRepresentation(distance[, …])

Representation of the distance of points from the origin.

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[, copy])

High-level object providing a flexible interface for celestial coordinate representation, manipulation, and transformation between systems.

SkyCoordInfo([bound])

Container for meta information like name, description, format.

SkyOffsetFrame(*args, **kwargs)

A frame which is relative to some specific position and oriented to match its frame.

SphericalCosLatDifferential(d_lon_coslat, …)

Differential(s) of points in 3D spherical coordinates.

SphericalDifferential(d_lon, d_lat, d_distance)

Differential(s) of points in 3D spherical coordinates.

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[, copy, …])

Supergalactic Coordinates (see Lahav et al.

TimeAttribute([default, secondary_attribute])

Frame attribute descriptor for quantities that are Time objects.

TimeFrameAttribute(*args, **kwargs)

TransformGraph()

A graph representing the paths between coordinate frames.

UnitSphericalCosLatDifferential(…[, copy])

Differential(s) of points on a unit sphere.

UnitSphericalDifferential(d_lon, d_lat[, copy])

Differential(s) of points on a unit sphere.

UnitSphericalRepresentation(lon, lat[, …])

Representation of points on a unit sphere.

UnknownSiteException(site, attribute[, …])

solar_system_ephemeris()

Default ephemerides for calculating positions of Solar-System bodies.

Class Inheritance Diagram

Inheritance diagram of astropy.coordinates.transformations.AffineTransform, astropy.coordinates.builtin_frames.altaz.AltAz, astropy.coordinates.angles.Angle, astropy.coordinates.attributes.Attribute, astropy.coordinates.builtin_frames.ecliptic.BarycentricMeanEcliptic, astropy.coordinates.builtin_frames.ecliptic.BarycentricTrueEcliptic, astropy.coordinates.transformations.BaseAffineTransform, astropy.coordinates.baseframe.BaseCoordinateFrame, astropy.coordinates.representation.BaseDifferential, astropy.coordinates.builtin_frames.ecliptic.BaseEclipticFrame, astropy.coordinates.builtin_frames.baseradec.BaseRADecFrame, astropy.coordinates.representation.BaseRepresentation, astropy.coordinates.representation.BaseRepresentationOrDifferential, astropy.coordinates.representation.BaseSphericalCosLatDifferential, astropy.coordinates.representation.BaseSphericalDifferential, astropy.coordinates.errors.BoundsError, astropy.coordinates.builtin_frames.cirs.CIRS, astropy.coordinates.representation.CartesianDifferential, astropy.coordinates.representation.CartesianRepresentation, astropy.coordinates.attributes.CartesianRepresentationAttribute, astropy.coordinates.attributes.CartesianRepresentationFrameAttribute, astropy.coordinates.transformations.CompositeTransform, astropy.coordinates.errors.ConvertError, astropy.coordinates.attributes.CoordinateAttribute, astropy.coordinates.transformations.CoordinateTransform, astropy.coordinates.builtin_frames.ecliptic.CustomBarycentricEcliptic, astropy.coordinates.representation.CylindricalDifferential, astropy.coordinates.representation.CylindricalRepresentation, astropy.coordinates.attributes.DifferentialAttribute, astropy.coordinates.distances.Distance, astropy.coordinates.transformations.DynamicMatrixTransform, astropy.coordinates.earth.EarthLocation, astropy.coordinates.attributes.EarthLocationAttribute, astropy.coordinates.builtin_frames.fk4.FK4, astropy.coordinates.builtin_frames.fk4.FK4NoETerms, astropy.coordinates.builtin_frames.fk5.FK5, astropy.coordinates.transformations.FunctionTransform, astropy.coordinates.transformations.FunctionTransformWithFiniteDifference, astropy.coordinates.builtin_frames.gcrs.GCRS, astropy.coordinates.builtin_frames.galactic.Galactic, astropy.coordinates.builtin_frames.lsr.GalacticLSR, astropy.coordinates.builtin_frames.galactocentric.Galactocentric, astropy.coordinates.baseframe.GenericFrame, astropy.coordinates.builtin_frames.ecliptic.GeocentricMeanEcliptic, astropy.coordinates.builtin_frames.ecliptic.GeocentricTrueEcliptic, astropy.coordinates.builtin_frames.hcrs.HCRS, astropy.coordinates.builtin_frames.ecliptic.HeliocentricEclipticIAU76, astropy.coordinates.builtin_frames.ecliptic.HeliocentricMeanEcliptic, 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.builtin_frames.lsr.LSR, astropy.coordinates.angles.Latitude, astropy.coordinates.angles.Longitude, astropy.coordinates.representation.PhysicsSphericalDifferential, astropy.coordinates.representation.PhysicsSphericalRepresentation, astropy.coordinates.builtin_frames.gcrs.PrecessedGeocentric, astropy.coordinates.attributes.QuantityAttribute, astropy.coordinates.attributes.QuantityFrameAttribute, astropy.coordinates.representation.RadialDifferential, astropy.coordinates.representation.RadialRepresentation, astropy.coordinates.errors.RangeError, astropy.coordinates.baseframe.RepresentationMapping, astropy.coordinates.sky_coordinate.SkyCoord, astropy.coordinates.sky_coordinate.SkyCoordInfo, astropy.coordinates.builtin_frames.skyoffset.SkyOffsetFrame, astropy.coordinates.representation.SphericalCosLatDifferential, astropy.coordinates.representation.SphericalDifferential, astropy.coordinates.representation.SphericalRepresentation, astropy.coordinates.transformations.StaticMatrixTransform, astropy.coordinates.builtin_frames.supergalactic.Supergalactic, astropy.coordinates.attributes.TimeAttribute, astropy.coordinates.attributes.TimeFrameAttribute, astropy.coordinates.transformations.TransformGraph, astropy.coordinates.representation.UnitSphericalCosLatDifferential, astropy.coordinates.representation.UnitSphericalDifferential, astropy.coordinates.representation.UnitSphericalRepresentation, astropy.coordinates.errors.UnknownSiteException, astropy.coordinates.solar_system.solar_system_ephemeris