astro

  • Astronomical and nutation routines

  • Computes the basic astronomical mean longitudes and other fundamental orbital parameters

  • Computes the solar and lunar positions in Earth-Centered Earth-Fixed (ECEF) coordinates

Calling Sequence

import pyTMD.astro
S,H,P,N,PP = pyTMD.astro.mean_longitudes(MJD)

Source code

pyTMD.astro.mean_longitudes(MJD: ndarray, **kwargs)[source]

Computes the basic astronomical mean longitudes: \(S\), \(H\), \(P\), \(N\) and \(P_s\) [33, 51]

Note \(N\) is not \(N'\), i.e. \(N\) is decreasing with time.

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

method: str, default ‘Cartwright’

method for calculating mean longitudes

  • 'Cartwright': use coefficients from David Cartwright

  • 'Meeus': use coefficients from Meeus Astronomical Algorithms

  • 'ASTRO5': use Meeus Astronomical coefficients from ASTRO5

  • 'IERS': convert from IERS Delaunay arguments

Returns:
S: np.ndarray

mean longitude of moon (degrees)

H: np.ndarray

mean longitude of sun (degrees)

P: np.ndarray

mean longitude of lunar perigee (degrees)

N: np.ndarray

mean longitude of ascending lunar node (degrees)

Ps: np.ndarray

longitude of solar perigee (degrees)

pyTMD.astro.doodson_arguments(MJD: ndarray, equinox: bool = False, apply_correction: bool = True)[source]

Computes astronomical phase angles for the six Doodson Arguments: \(\tau\), \(S\), \(H\), \(P\), \(N'\), and \(P_s\) [16, 33]

Follows IERS conventions for the Doodson arguments [42]

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

equinox: bool, default False

use equinox method for calculating mean lunar time

apply_correction: bool, default True

Apply correction for mean lunar longitude

Returns:
TAU: np.ndarray

mean lunar time (radians)

S: np.ndarray

mean longitude of the moon (radians)

H: np.ndarray

mean longitude of the sun (radians)

P: np.ndarray

mean longitude of lunar perigee (radians)

Np: np.ndarray

negative mean longitude of the ascending node (radians)

Ps: np.ndarray

mean longitude of solar perigee (radians)

pyTMD.astro.delaunay_arguments(MJD: ndarray)[source]

Computes astronomical phase angles for the five primary Delaunay Arguments of Nutation: \(l\), \(l'\), \(F\), \(D\), and \(N\) [7, 33, 42]

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

Returns:
l: np.ndarray

mean anomaly of moon (radians)

lp: np.ndarray

mean anomaly of the sun (radians)

F: np.ndarray

mean argument of the moon (radians)

D: np.ndarray

mean elongation of the moon from the sun (radians)

N: np.ndarray

mean longitude of ascending lunar node (radians)

pyTMD.astro.schureman_arguments(P: ndarray, N: ndarray)[source]

Computes additional phase angles \(I\), \(\xi\), \(\nu\), \(R\), \(R_a\), \(\nu'\), and \(\nu''\) from Schureman [50]

See the explanation of symbols in appendix of Schureman [50]

Parameters:
P: np.ndarray

mean longitude of lunar perigee (radians)

N: np.ndarray

mean longitude of ascending lunar node (radians)

Returns:
I: np.ndarray

obliquity of lunar orbit with respect to Earth’s equator (radians)

xi: np.ndarray

longitude in the moon’s orbit of lunar intersection (radians)

nu: np.ndarray

right ascension of lunar intersection (radians)

Qa: np.ndarray

factor in amplitude for m1 constituent (radians)

Qu: np.ndarray

term in argument for m1 constituent (radians)

Ra: np.ndarray

factor in amplitude for l2 constituent (radians)

Ru: np.ndarray

term in argument for l2 constituent (radians)

nu_p: np.ndarray

term in argument for k1 constituent (radians)

nu_s: np.ndarray

term in argument for k2 constituent (radians)

pyTMD.astro.mean_obliquity(MJD: ndarray)[source]

Mean obliquity of the ecliptic [6, 7]

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

Returns:
epsilon: np.ndarray

mean obliquity of the ecliptic (radians)

pyTMD.astro.equation_of_time(MJD: ndarray)[source]

Approximate calcuation of the difference between apparent and mean solar times [33, 55]

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

Returns:
E: np.ndarray

equation of time (radians)

pyTMD.astro.solar_ecef(MJD: ndarray, **kwargs)[source]

Wrapper function for calculating the positional coordinates of the sun in an Earth-centric, Earth-Fixed (ECEF) frame [33, 34, 40]

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

ephemerides: str, default ‘approximate’

Method for calculating solar ephemerides

  • 'approximate': low-resolution ephemerides

  • 'JPL': computed ephemerides from JPL kernels

**kwargs: dict

Keyword options for ephemeris calculation

Returns:
X, Y, Z: np.ndarray

ECEF coordinates of the sun (meters)

pyTMD.astro.solar_approximate(MJD, **kwargs)[source]

Computes approximate positional coordinates of the sun in an Earth-centric, Earth-Fixed (ECEF) frame [33, 34]

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

Returns:
X, Y, Z: np.ndarray

ECEF coordinates of the sun (meters)

pyTMD.astro.solar_ephemerides(MJD: ndarray, **kwargs)[source]

Computes positional coordinates of the sun in an Earth-centric, Earth-Fixed (ECEF) frame using JPL ephemerides [33, 40]

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

kernel: str or pathlib.Path

Path to JPL ephemerides kernel file

include_aberration: bool, default False

Correct for aberration effects

Returns:
X, Y, Z: np.ndarray

ECEF coordinates of the sun (meters)

pyTMD.astro.lunar_ecef(MJD: ndarray, **kwargs)[source]

Wrapper function for calculating the positional coordinates of the moon in an Earth-centric, Earth-Fixed (ECEF) frame [33, 34, 40]

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

ephemerides: str, default ‘approximate’

Method for calculating lunar ephemerides

  • 'approximate': low-resolution ephemerides

  • 'JPL': computed ephemerides from JPL kernels

**kwargs: dict

Keyword options for ephemeris calculation

Returns:
X, Y, Z: np.ndarray

ECEF coordinates of the moon (meters)

pyTMD.astro.lunar_approximate(MJD, **kwargs)[source]

Computes approximate positional coordinates of the moon in an Earth-centric, Earth-Fixed (ECEF) frame [33, 34]

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

Returns:
X, Y, Z: np.ndarray

ECEF coordinates of the moon (meters)

pyTMD.astro.lunar_ephemerides(MJD: ndarray, **kwargs)[source]

Computes positional coordinates of the moon in an Earth-centric, Earth-Fixed (ECEF) frame using JPL ephemerides [33, 40]

Parameters:
MJD: np.ndarray

Modified Julian Day (MJD) of input date

kernel: str or pathlib.Path

Path to JPL ephemerides kernel file

include_aberration: bool, default False

Correct for aberration effects

Returns:
X, Y, Z: np.ndarray

ECEF coordinates of the moon (meters)

pyTMD.astro.gast(T: float | ndarray)[source]

Greenwich Apparent Sidereal Time (GAST) [6, 7, 42]

Parameters:
T: np.ndarray

Centuries since 2000-01-01T12:00:00

pyTMD.astro.itrs(T: float | ndarray, include_polar_motion: bool = True)[source]

International Terrestrial Reference System (ITRS) [6, 7, 42]

An Earth-centered Earth-fixed (ECEF) coordinate system combining the Earth’s true equator and equinox of date, the Earth’s rotation with respect to the stars, and the polar wobble of the crust with respect to the pole of rotation

Parameters:
T: np.ndarray

Centuries since 2000-01-01T12:00:00

include_polar_motion: bool, default True

Include polar motion in the rotation matrix

pyTMD.astro._eqeq_complement(T: float | ndarray)[source]

Compute complementary terms of the equation of the equinoxes [6, 7, 42]

These include the combined effects of precession and nutation [26, 42, 55]

Parameters:
T: np.ndarray

Centuries since 2000-01-01T12:00:00

pyTMD.astro._icrs_rotation_matrix(T: float | ndarray, include_polar_motion: bool = True)[source]

Rotation matrix for tranforming from the International Celestial Reference System (ICRS) to the International Terrestrial Reference System (ITRS) [6, 7, 42]

Parameters:
T: np.ndarray

Centuries since 2000-01-01T12:00:00

include_polar_motion: bool, default True

Include polar motion in the rotation matrix

pyTMD.astro._frame_bias_matrix()[source]

Frame bias rotation matrix for converting from a dynamical reference system to the International Celestial Reference System (ICRS) [42, 55]

pyTMD.astro._nutation_angles(T: float | ndarray)[source]

Calculate nutation rotation angles using tables from IERS Conventions [42]

Parameters:
T: np.ndarray

Centuries since 2000-01-01T12:00:00

Returns:
dpsi: np.ndarray

Nutation in longitude

deps: np.ndarray

Obliquity of the ecliptic

pyTMD.astro._nutation_matrix(mean_obliquity: float | ndarray, true_obliquity: float | ndarray, psi: float | ndarray)[source]

Nutation rotation matrix [27, 42]

Parameters:
mean_obliquity: np.ndarray

Mean obliquity of the ecliptic

true_obliquity: np.ndarray

True obliquity of the ecliptic

psi: np.ndarray

Nutation in longitude

pyTMD.astro._polar_motion_matrix(T: float | ndarray)[source]

Polar motion (Earth Orientation Parameters) rotation matrix [42, 55]

Parameters:
T: np.ndarray

Centuries since 2000-01-01T12:00:00

pyTMD.astro._precession_matrix(T: float | ndarray)[source]

Precession rotation matrix [6, 7, 28]

Parameters:
T: np.ndarray

Centuries since 2000-01-01T12:00:00

pyTMD.astro._correct_aberration(position, velocity)[source]

Correct a relative position for aberration effects [27]

Parameters:
position: np.ndarray

Position vector in astronomical units

velocity: np.ndarray

Velocity vector in astronomical units per day

pyTMD.astro._parse_table_5_2e()[source]

Parse table with expressions for Greenwich Sidereal Time provided in Chapter 5 of Petit and Luzum [42]

pyTMD.astro._parse_table_5_3a()[source]

Parse table with IAU 2000A lunisolar and planetary components of nutation in longitude provided in Chapter 5 of Petit and Luzum [42]

pyTMD.astro._parse_table_5_3b()[source]

Parse table with IAU 2000A lunisolar and planetary components of nutation in obliquity provided in Chapter 5 of Petit and Luzum [42]