Celestial Coordinates

class coord.CelestialCoord(ra, dec=None)[source]

This class defines a position on the celestial sphere, normally given by two angles, ra and dec.

This class can be used to perform various calculations in spherical coordinates, such as finding the angular distance between two points in the sky, calculating the angles in spherical triangles, projecting from sky coordinates onto a Euclidean tangent plane, etc.

Initialization:

A CelestialCoord object is constructed from the right ascension and declination:

coord.CelestialCoord.__init__()

>>> c = CelestialCoord(ra=12*hours, dec=31*degrees)
>>> print(c)
coord.CelestialCoord(3.141592653589793 radians, 0.5410520681182421 radians)

Attributes:

A CelestialCoord has the following (read-only) attributes:

ra

The right ascension (an Angle instance)

dec

The declination (an Angle instance)

>>> print(c.ra / degrees, c.dec / degrees)
180.0 31.0

In addition there is a convenience access property that returns ra and dec in radians.

rad

A tuple (ra.rad, dec.rad)

>>> print(c.rad)
(3.141592653589793, 0.5410520681182421)

Spherical Geometry:

The basic spherical geometry operations are available to work with spherical triangles

For three coordinates cA, cB, cC making a spherical triangle, one can calculate the sides and angles via:

>>> cA = CelestialCoord(0 * degrees, 0 * degrees)
>>> cB = CelestialCoord(0 * degrees, 10 * degrees)
>>> cC = CelestialCoord(10 * degrees, 0 * degrees)
>>> a = cB.distanceTo(cC)
>>> b = cC.distanceTo(cA)
>>> c = cA.distanceTo(cB)
>>> print(a.deg, b.deg, c.deg)
14.106044260566366 10.0 10.0
>>> A = cA.angleBetween(cB, cC)
>>> B = cB.angleBetween(cC, cA)
>>> C = cC.angleBetween(cA, cB)
>>> print(A.deg, B.deg, C.deg)
90.0 45.43854858674231 45.43854858674231

Projections:

Local tangent plane projections of an area of the sky can be performed using the project method:

coord.CelestialCoord.project()

>>> center = CelestialCoord(ra=10*hours, dec=30*degrees)
>>> sky_coord = CelestialCoord(ra=10.5*hours, dec=31*degrees)
>>> print(sky_coord)
coord.CelestialCoord(2.748893571891069 radians, 0.5410520681182421 radians)
>>> u, v = center.project(sky_coord)
>>> print(u.deg, v.deg)
-6.452371275343261 1.21794987288635

and back:

coord.CelestialCoord.deproject()

>>> sky_coord = center.deproject(u,v)
>>> print(sky_coord)
coord.CelestialCoord(2.748893571891069 radians, 0.5410520681182421 radians)

where u and v are Angles and center and sky_coord are CelestialCoords.

__init__(ra, dec=None)[source]
Parameters
  • ra – The right ascension. Must be an Angle instance.

  • dec – The declination. Must be an Angle instance.

property ra

A read-only attribute, giving the Right Ascension as an Angle

property dec

A read-only attribute, giving the Declination as an Angle

property rad

A convenience property, giving a tuple (ra.rad, dec.rad)

get_xyz()[source]

Get the (x,y,z) coordinates on the unit sphere corresponding to this (RA, Dec).

\[\begin{split}x &= \cos(dec) \cos(ra) \\ y &= \cos(dec) \sin(ra) \\ z &= \sin(dec)\end{split}\]
Returns

a tuple (x,y,z)

static from_xyz(x, y, z)[source]

Construct a CelestialCoord from a given (x,y,z) position in three dimensions.

The 3D (x,y,z) position does not need to fall on the unit sphere. The RA, Dec will be inferred from the relations:

\[\begin{split}x &= r \cos(dec) \cos(ra) \\ y &= r \cos(dec) \sin(ra) \\ z &= r \sin(dec)\end{split}\]

where \(r\) is arbitrary.

Parameters
  • x – The x position in 3 dimensions. Corresponds to r cos(dec) cos(ra)

  • y – The y position in 3 dimensions. Corresponds to r cos(dec) sin(ra)

  • z – The z position in 3 dimensions. Corresponds to r sin(dec)

Returns

a CelestialCoord instance

static radec_to_xyz(ra, dec, r=1.0)[source]

Convert ra, dec (in radians) to 3D x,y,z coordinates on the unit sphere.

The connection between (ra,dec) and (x,y,z) are given by the following formulae: .. math:

x &= r \cos(dec) \cos(ra)  \\
y &= r \cos(dec) \sin(ra)  \\
z &= r \sin(dec)

For a single ra,dec pair, the following are essentially equivalent:

>>> ra = 12*hours/radians       # May be any angle measured
>>> dec = 31*degrees/radians    # in radians
>>> CelestialCoord.radec_to_xyz(ra, dec)
(-0.8571673007021123, 1.0497271911386187e-16, 0.5150380749100542)
>>> CelestialCoord(ra * radians, dec * radians).get_xyz()
(-0.8571673007021123, 1.0497271911386187e-16, 0.5150380749100542)

However, the advantage of this function is that the input values may be numpy arrays, in which case, the return values will also be numpy arrays.

Parameters
  • ra – The right ascension(s) in radians. May be a numpy array.

  • dec – The declination(s) in radians. May be a numpy array.

  • r – The distance(s) from Earth (default 1.). May be a numpy array.

Returns

x, y, z as a tuple.

static xyz_to_radec(x, y, z, return_r=False)[source]

Convert 3D x,y,z coordinates to ra, dec (in radians).

The connection between (ra,dec) and (x,y,z) are given by the following formulae: .. math:

x &= r \cos(dec) \cos(ra)  \\
y &= r \cos(dec) \sin(ra)  \\
z &= r \sin(dec)

For a single (x,y,z) position, the following are essentially equivalent:

>>> x = 0.839       # May be any any 3D location
>>> y = 0.123       # Not necessarily on unit sphere
>>> z = 0.530
>>> CelestialCoord.xyz_to_radec(x, y, z)
(0.14556615088111796, 0.558616191048523)
>>> c = CelestialCoord.from_xyz(x, y, z)
>>> c.ra.rad, c.dec.rad
(0.145566150881118, 0.558616191048523)

However, the advantage of this function is that the input values may be numpy arrays, in which case, the return values will also be numpy arrays.

Parameters
  • x – The x position(s) in 3 dimensions. May be a numpy array.

  • y – The y position(s) in 3 dimensions. May be a numpy array.

  • z – The z position(s) in 3 dimensions. May be a numpy array.

  • return_r – Whether to return r as well as ra, dec. (default: False)

Returns

ra, dec as a tuple. Or if return_r is True, (ra, dec, r).

normal()[source]

Return the coordinate in the “normal” convention of having 0 <= ra < 24 hours.

This convention is not enforced on construction, so this function exists to make it easy to convert if desired.

Functions such as from_galactic and from_xyz will return normal coordinates.

distanceTo(coord2)[source]

Returns the great circle distance between this coord and another one. The return value is an Angle object

Parameters

coord2 – The CelestialCoord to calculate the distance to.

Returns

the great circle distance to coord2.

greatCirclePoint(coord2, theta)[source]

Returns a point on the great circle connecting self and coord2.

Two points, c1 and c2, on the unit sphere define a great circle (so long as the two points are not either coincident or antipodal). We can define points on this great circle by their angle from c1, such that the angle for c2 has 0 < theta2 < pi. I.e. theta increases from 0 as the points move from c1 towards c2.

This function then returns the coordinate on this great circle (where c1 is self and c2 is coord2) that corresponds to the given angle theta.

Parameters
  • coord2 – Another CelestialCoord defining the great circle to use.

  • theta – The Angle along the great circle corresponding to the desired point.

Returns

the corresponding CelestialCoord

angleBetween(coord2, coord3)[source]

Find the open angle at the location of the current coord between coord2 and coord3.

The current coordinate along with the two other coordinates form a spherical triangle on the sky. This function calculates the angle between the two sides at the location of the current coordinate.

Note that this returns a signed angle. The angle is positive if the sweep direction from coord2 to coord3 is counter-clockwise (as observed from Earth). It is negative if the direction is clockwise.

Parameters
  • coord2 – A second CelestialCoord

  • coord3 – A third CelestialCoord

Returns

the angle between the great circles joining the other two coordinates to the current coordinate.

area(coord2, coord3)[source]

Find the area of a spherical triangle in steradians.

The current coordinate along with the two other coordinates form a spherical triangle on the sky. This function calculates the area of that spherical triangle, which is measured in steradians (i.e. surface area of the triangle on the unit sphere).

Parameters
  • coord2 – A second CelestialCoord

  • coord3 – A third CelestialCoord

Returns

the area in steradians of the given spherical triangle.

project(coord2, projection=None)[source]

Use the currect coord as the center point of a tangent plane projection to project the coord2 coordinate onto that plane.

This function return a tuple (u,v) in the Euclidean coordinate system defined by a tangent plane projection around the current coordinate, with +v pointing north and +u pointing west. (i.e. to the right on the sky if +v is up.)

There are currently four options for the projection, which you can specify with the optional projection keyword argument:

gnomonic

[default] uses a gnomonic projection (i.e. a projection from the center of the sphere, which has the property that all great circles become straight lines. For more information, see http://mathworld.wolfram.com/GnomonicProjection.html This is the usual TAN projection used by most FITS images.

stereographic

uses a stereographic proejection, which preserves angles, but not area. For more information, see http://mathworld.wolfram.com/StereographicProjection.html

lambert

uses a Lambert azimuthal projection, which preserves area, but not angles. For more information, see http://mathworld.wolfram.com/LambertAzimuthalEqual-AreaProjection.html

postel

uses a Postel equidistant proejection, which preserves distances from the projection point, but not area or angles. For more information, see http://mathworld.wolfram.com/AzimuthalEquidistantProjection.html

The distance or angle errors increase with distance from the projection point of course.

Parameters
  • coord2 – The coordinate to project onto the tangent plane.

  • projection – The name of the projection to be used. [default: gnomonic, see above for other options]

Returns

(u,v) as Angle instances

project_rad(ra, dec, projection=None)[source]

This is basically identical to the project() function except that the input ra, dec are given in radians rather than packaged as a CelestialCoord object and the returned u,v are given in radians.

The main advantage to this is that it will work if ra and dec are NumPy arrays, in which case the output u, v will also be NumPy arrays.

Parameters
  • ra – The right ascension in radians to project onto the tangent plane.

  • dec – The declination in radians to project onto the tangent plane.

  • projection – The name of the projection to be used. [default: gnomonic, see project docstring for other options]

Returns

(u,v) in radians

deproject(u, v, projection=None)[source]

Do the reverse process from the project() function.

i.e. This takes in a position (u,v) and returns the corresponding celestial coordinate, using the current coordinate as the center point of the tangent plane projection.

Parameters
  • u – The u position on the tangent plane to deproject (must be an Angle instance)

  • v – The v position on the tangent plane to deproject (must be an Angle instance)

  • projection – The name of the projection to be used. [default: gnomonic, see project docstring for other options]

Returns

the corresponding CelestialCoord for that position.

deproject_rad(u, v, projection=None)[source]

This is basically identical to the deproject() function except that the output ra, dec are returned as a tuple (ra, dec) in radians rather than packaged as a CelestialCoord object and u and v are in radians rather than Angle instances.

The main advantage to this is that it will work if u and v are NumPy arrays, in which case the output ra, dec will also be NumPy arrays.

Parameters
  • u – The u position in radians on the tangent plane to deproject

  • v – The v position in radians on the tangent plane to deproject

  • projection – The name of the projection to be used. [default: gnomonic, see project docstring for other options]

Returns

the corresponding RA, Dec in radians

jac_deproject(u, v, projection=None)[source]

Return the jacobian of the deprojection.

i.e. if the input position is (u,v) then the return matrix is

\[\begin{split}J \equiv \begin{bmatrix} J00 & J01 \\ J10 & J11 \end{bmatrix} = \begin{bmatrix} d\textrm{ra}/du \cos(\textrm{dec}) & d\textrm{ra}/dv \cos(\textrm{dec}) \\ d\textrm{dec}/du & d\textrm{dec}/dv \end{bmatrix}\end{split}\]
Parameters
  • u – The u position (as an Angle instance) on the tangent plane

  • v – The v position (as an Angle instance) on the tangent plane

  • projection – The name of the projection to be used. [default: gnomonic, see project docstring for other options]

Returns

the Jacobian as a 2x2 numpy array [[J00, J01], [J10, J11]]

jac_deproject_rad(u, v, projection=None)[source]

Equivalent to jac_deproject, but the inputs are in radians and may be numpy arrays.

Parameters
  • u – The u position (in radians) on the tangent plane

  • v – The v position (in radians) on the tangent plane

  • projection – The name of the projection to be used. [default: gnomonic, see project docstring for other options]

Returns

the Jacobian as a 2x2 numpy array [[J00, J01], [J10, J11]]

precess(from_epoch, to_epoch)[source]
This function precesses equatorial ra and dec from one epoch to another.

It is adapted from a set of fortran subroutines based on (a) pages 30-34 of the Explanatory Supplement to the AE, (b) Lieske, et al. (1977) A&A 58, 1-16, and (c) Lieske (1979) A&A 73, 282-284.

Parameters
  • from_epoch – The epoch of the current coordinate

  • to_epoch – The epoch of the returned precessed coordinate

Returns

a CelestialCoord object corresponding to the precessed position.

galactic(epoch=2000.0)[source]

Get the longitude and latitude in galactic coordinates corresponding to this position.

Parameters

epoch – The epoch of the current coordinate. [default: 2000.]

Returns

the longitude and latitude as a tuple (el, b), given as Angle instances.

static from_galactic(el, b, epoch=2000.0)[source]

Create a CelestialCoord from the given galactic coordinates

Parameters
  • el – The longitude in galactic coordinates (an Angle instance)

  • b – The latitude in galactic coordinates (an Angle instance)

  • epoch – The epoch of the returned coordinate. [default: 2000.]

Returns

the CelestialCoord corresponding to these galactic coordinates.

ecliptic(epoch=2000.0, date=None)[source]

Get the longitude and latitude in ecliptic coordinates corresponding to this position.

The epoch parameter is used to get an accurate value for the (time-varying) obliquity of the ecliptic. The formulae for this are quite straightforward. It requires just a single parameter for the transformation, the obliquity of the ecliptic (the Earth’s axial tilt).

Parameters
  • epoch – The epoch to be used for estimating the obliquity of the ecliptic, if date is None. But if date is given, then use that to determine the epoch. [default: 2000.]

  • date – If a date is given as a python datetime object, then return the position in ecliptic coordinates with respect to the sun position at that date. If None, then return the true ecliptic coordiantes. [default: None]

Returns

the longitude and latitude as a tuple (lambda, beta), given as Angle instances.

static from_ecliptic(lam, beta, epoch=2000.0, date=None)[source]

Create a CelestialCoord from the given ecliptic coordinates

Parameters
  • lam – The longitude in ecliptic coordinates (an Angle instance)

  • beta – The latitude in ecliptic coordinates (an Angle instance)

  • epoch – The epoch to be used for estimating the obliquity of the ecliptic, if date is None. But if date is given, then use that to determine the epoch. [default: 2000.]

  • date – If a date is given as a python datetime object, then return the position in ecliptic coordinates with respect to the sun position at that date. If None, then return the true ecliptic coordiantes. [default: None]

Returns

the CelestialCoord corresponding to these ecliptic coordinates.

coord._CelestialCoord(ra, dec)[source]

Equivalent to CeletialCoord(ra,dec), but without the normal sanity checks.

Parameters
  • ra – The right ascension. Must be an Angle instance.

  • dec – The declination. Must be an Angle instance.