orbit_conversion

This file copies the orbit_conversion_utilities as well as the ecliptic rotation matrix from sorcha and adapts them to jax, implementing the autograd jacobians needed for this process

Attributes

OBLIQUITY_ECLIPTIC

ECL_TO_EQ_ROTATION_MATRIX

EQ_TO_ECL_ROTATION_MATRIX

jac_cometary_xyz

jac_keplerian_xyz

Functions

create_ecl_to_eq_rotation_matrix(ecl)

Creates a rotation matrix for transforming ecliptical coordinates

stumpff(x)

Computes the Stumpff function c_k(x) for k = 0, 1, 2, 3

root_function(s, mu, alpha, r0, r0dot, t)

Root function used in the Halley minimizer

halley_safe(x1, x2, mu, alpha, r0, r0dot, t[, xacc, maxit])

Applies the Halley root finding algorithm on the universal Kepler equation

universal_cartesian(mu, q, e, incl, longnode, argperi, ...)

Converts from a series of orbital elements into state vectors

principal_value(theta)

Computes the principal value of an angle

atan2_checkzero(x, y)

eccanom(e, trueanom, mu, alpha, p)

paranom(e, trueanom, mu, alpha, p)

hypanom(e, trueanom, mu, alpha, p)

universal_cometary(mu, x, y, z, vx, vy, vz, epochMJD_TDB)

Converts from a state vectors into cometary orbital elements

universal_keplerian(mu, x, y, z, vx, vy, vz, epochMJD_TDB)

Converts from a state vectors into Keplerian orbital elements

covariance_ecl_to_eq(covariance)

Converts a covariance matrix from ecliptic to equatorial coordinates.

covariance_eq_to_ecl(covariance)

Converts a covariance matrix from equatorial to ecliptic coordinates.

covariance_cometary_xyz(mu, x, y, z, vx, vy, vz, ...)

covariance_keplerian_xyz(mu, x, y, z, vx, vy, vz, ...)

covariance_xyz_cometary(mu, q, e, incl, longnode, ...)

covariance_xyz_keplerian(mu, a, e, incl, longnode, ...)

parse_covariance_row_to_CART(row, gm_total, gm_sun)

Parses a row of orbit data, unpacking the flattened covariance matrix

Module Contents

OBLIQUITY_ECLIPTIC[source]
create_ecl_to_eq_rotation_matrix(ecl)[source]

Creates a rotation matrix for transforming ecliptical coordinates to equatorial coordinates. A rotation matrix based on the solar system’s ecliptic obliquity is already provided as ECL_TO_EQ_ROTATION_MATRIX.

Parameters:

ecl (float) – The ecliptical obliquity.

Returns:

rotmat – rotation matrix for transofmring ecliptical coordinates to equatorial coordinates. Array has shape (3,3).

Return type:

numpy array/matrix of floats

ECL_TO_EQ_ROTATION_MATRIX[source]
EQ_TO_ECL_ROTATION_MATRIX[source]
stumpff(x)[source]

Computes the Stumpff function c_k(x) for k = 0, 1, 2, 3

Parameters:

x (float) – Argument of the Stumpff function

Returns:

  • c_0(x) (float)

  • c_1(x) (float)

  • c_2(x) (float)

  • c_3(x) (float)

root_function(s, mu, alpha, r0, r0dot, t)[source]

Root function used in the Halley minimizer Computes the zeroth, first, second, and third derivatives of the universal Kepler equation f

Parameters:
  • s (float) – Eccentric anomaly

  • mu (float) – Standard gravitational parameter GM

  • alpha (float) – Total energy

  • r0 (float) – Initial position

  • r0dot (float) – Initial velocity

  • t (float) – Time

Returns:

  • f (float) – universal Kepler equation)

  • fp (float) – (first derivative of f

  • fpp (float) – second derivative of f

  • fppp (float) – third derivative of f

halley_safe(x1, x2, mu, alpha, r0, r0dot, t, xacc=1e-14, maxit=100)[source]

Applies the Halley root finding algorithm on the universal Kepler equation

Parameters:
  • x1 (float) – Previous guess used in minimization

  • x2 (float) – Current guess for minimization

  • mu (float) – Standard gravitational parameter GM

  • alpha (float) – Total energy

  • r0 (float) – Initial position

  • r0dot (float) – Initial velocity

  • t (float) – Time

  • xacc (float) – Accuracy in x before algorithm declares convergence

  • maxit (int) – Maximum number of iterations

Returns:

  • boolean – True if minimization converged, False otherwise

  • float – Solution

  • float – First derivative of solution

universal_cartesian(mu, q, e, incl, longnode, argperi, tp, epochMJD_TDB)[source]

Converts from a series of orbital elements into state vectors using the universal variable formulation

The output vector will be oriented in the same system as the positional angles (i, Omega, omega)

Note that mu, q, tp and epochMJD_TDB must have compatible units As an example, if q is in au and tp/epoch are in days, mu must be in (au^3)/days^2

Parameters:
  • mu (float) – Standard gravitational parameter GM (see note above about units)

  • q (float) – Perihelion (see note above about units)

  • e (float) – Eccentricity

  • incl (float) – Inclination (radians)

  • longnode (float) – Longitude of ascending node (radians)

  • argperi (float) – Argument of perihelion (radians)

  • tp (float) – Time of perihelion passage in TDB scale (see note above about units)

  • epochMJD_TDB (float) – Epoch (in TDB) when the elements are defined (see note above about units)

Returns:

  • float – x coordinate

  • float – y coordinate

  • float – z coordinate

  • float – x velocity

  • float – y velocity

  • float – z velocity

principal_value(theta)[source]

Computes the principal value of an angle

Parameters:

theta (float) – Angle

Returns:

Principal value of angle

Return type:

float

atan2_checkzero(x, y)[source]
eccanom(e, trueanom, mu, alpha, p)[source]
paranom(e, trueanom, mu, alpha, p)[source]
hypanom(e, trueanom, mu, alpha, p)[source]
universal_cometary(mu, x, y, z, vx, vy, vz, epochMJD_TDB)[source]

Converts from a state vectors into cometary orbital elements using the universal variable formulation

The input vector will determine the orientation of the positional angles (i, Omega, omega)

Note that mu and the state vectors must have compatible units As an example, if x is in au and vx are in au/days, mu must be in (au^3)/days^2

Parameters:
  • mu (float) – Standard gravitational parameter GM (see note above about units)

  • x (float) – x coordinate

  • y (float) – y coordinate

  • z (float) – z coordinate

  • vx (float) – x velocity

  • vy (float) – y velocity

  • vz (float) – z velocity

  • (float) (epochMJD_TDB) – Epoch (in TDB) when the elements are defined (see note above about units)

Returns:

  • float – Perihelion (see note above about units)

  • float – Eccentricity

  • float – Inclination (radians)

  • float – Longitude of ascending node (radians)

  • float – Argument of perihelion (radians)

  • float – Time of perihelion passage in TDB scale (see note above about units)

universal_keplerian(mu, x, y, z, vx, vy, vz, epochMJD_TDB)[source]

Converts from a state vectors into Keplerian orbital elements using the universal variable formulation

The input vector will determine the orientation of the positional angles (i, Omega, omega)

Note that mu and the state vectors must have compatible units As an example, if x is in au and vx are in au/days, mu must be in (au^3)/days^2

Parameters:
  • mu (float) – Standard gravitational parameter GM (see note above about units)

  • x (float) – x coordinate

  • y (float) – y coordinate

  • z (float) – z coordinate

  • vx (float) – x velocity

  • vy (float) – y velocity

  • vz (float) – z velocity

  • (float) (epochMJD_TDB) – Epoch (in TDB) when the elements are defined (see note above about units)

Returns:

  • float – Semi-major axis (see note above about units)

  • float – Eccentricity

  • float – Inclination (radians)

  • float – Longitude of ascending node (radians)

  • float – Argument of perihelion (radians)

  • float – Mean anomaly (radians)

jac_cometary_xyz[source]
jac_keplerian_xyz[source]
covariance_ecl_to_eq(covariance)[source]

Converts a covariance matrix from ecliptic to equatorial coordinates.

Parameters:

covariance (numpy array) – The covariance matrix to convert.

Returns:

The converted covariance matrix.

Return type:

numpy array

covariance_eq_to_ecl(covariance)[source]

Converts a covariance matrix from equatorial to ecliptic coordinates.

Parameters:

covariance (numpy array) – The covariance matrix to convert.

Returns:

The converted covariance matrix.

Return type:

numpy array

covariance_cometary_xyz(mu, x, y, z, vx, vy, vz, epochMJD_TDB, covariance)[source]
covariance_keplerian_xyz(mu, x, y, z, vx, vy, vz, epochMJD_TDB, covariance)[source]
covariance_xyz_cometary(mu, q, e, incl, longnode, argperi, tp, epochMJD_TDB, covariance)[source]
covariance_xyz_keplerian(mu, a, e, incl, longnode, argperi, M, epochMJD_TDB, covariance)[source]
parse_covariance_row_to_CART(row, gm_total, gm_sun)[source]

Parses a row of orbit data, unpacking the flattened covariance matrix and converting it to an equatorial cartesian format regardless of the input format.

Note that there is not a meaningful distinction between cartesian and barycentric cartesian coordinates here.

The input format of the row is read from the “FORMAT” column and acceptable input formats are CART, COM, KEP, BCART, BCART_EQ, BCOM, and BKEP.

Parameters:
  • row (numpy structured array) – The row of data to parse the covariance matrix from.

  • gm_total (float) – The gravitational parameter for the total system.

  • gm_sun (float) – The gravitational parameter for the Sun.