Precession

class sofia_redux.scan.coordinate_systems.epoch.precession.Precession(from_epoch, to_epoch)[source]

Bases: ABC

Initialize a precession object.

The precession object is used to precess equatorial coordinates from one epoch to another using the procedure found in:

Lederle, T. and Schwan, H., “Procedure for computing the apparent

places of fundamental stars (APFS) from 1984 onwards”, Astronomy and Astrophysics, vol. 134, no. 1, pp. 1–6, 1984.

Parameters:
from_epochEpoch

The epoch to precess from.

to_epochEpoch

The epoch to precess to.

Attributes Summary

JULIAN_CENTURY

YEAR

YEAR_TO_CENTURY

singular_epoch

Return whether the precession consists of single time epochs.

Methods Summary

calculate_matrix()

Precess the coordinates from the from_epoch to the to_epoch.

copy()

Return a copy of the precession.

precess(equatorial)

Precess the coordinates in-place.

r2(phi)

Calculate the R2 matrix.

r3(phi)

Calculate the R3 matrix.

Attributes Documentation

JULIAN_CENTURY = <Quantity 36525. d>
YEAR = <Quantity 365.24219879 d>
YEAR_TO_CENTURY = 0.009999786414510608
singular_epoch

Return whether the precession consists of single time epochs.

Returns:
is_singularbool

Methods Documentation

calculate_matrix()[source]

Precess the coordinates from the from_epoch to the to_epoch.

The precession matrix is a 3x3 matrix given by:

m = x1.x2.x3

where “.” indicates dot matrix multiplication and

x1 = | cos(z) sin(z) 0 |
-sin(z) cos(z) 0 |
0 0 1 |
x2 = | cos(theta) 0 sin(theta) |
0 1 0 |
sin(theta) 0 cos(theta) |
x3 = | cos(eta) sin(eta) 0 |
-sin(eta) cos(eta) 0 |
0 0 1 |
z = ((2305.6997 + (1.39744 + 0.000060 * tau) * tau +

(1.09543 + 0.000390 * tau + 0.018326 * t) * t) * t)

theta = ((2003.8746 - (0.85405 + 0.000370 * tau) * tau -

(0.42707 + 0.000370 * tau + 0.041803 * t) * t) * t)

eta = ((2305.6997 + (1.39744 + 0.000060 * tau) * tau +

(0.30201 - 0.000270 * tau + 0.017996 * t) * t) * t)

tau = (epoch1_year - 2000) * year_to_century

t = (epoch2_year - epoch1_year) * year_to_century

year_to_century = 365.24219879 / 36525

Note that z, theta and eta are in arcseconds.

Returns:
None
copy()[source]

Return a copy of the precession.

Returns:
Precession
precess(equatorial)[source]

Precess the coordinates in-place.

Parameters:
equatorialEquatorialCoordinates

The equatorial coordinates to precess.

Returns:
None
static r2(phi)[source]

Calculate the R2 matrix.

The R2 matrix is a 3x3 array of the form:

[[cos(phi), 0, sin(phi)],

[0, 1, 0 ], [sin(phi), 0, cos(phi)]]

Parameters:
phiastropy.units.Quantity

The phi angle. Either a single value or an array of shape (shape,) may be provided.

Returns:
R2numpy.ndarray (float)

The R2 matrix. If a single phi` was provided, the array will be of shape (3, 3). Otherwise, it will be of shape (phi.shape, 3, 3).

static r3(phi)[source]

Calculate the R3 matrix.

The R3 matrix is a 3x3 array of the form:

[[cos(phi), sin(phi), 0],

[-sin(phi), cos(phi), 0], [0, 0, 1]]

Parameters:
phiastropy.units.Quantity

The phi angle.

Returns:
R2numpy.ndarray (float)

The R3 matrix. If a single phi` was provided, the array will be of shape (3, 3). Otherwise, it will be of shape (phi.shape, 3, 3).