element_conversion#

Functions for converting between sets of orbital elements.

This module provide a variety of functions and classes to convert between different representations of translational and rotational states (e.g. Cartesian ↔ Keplerian).

Note

Rotations between different reference frames are provided in the frame_conversion module.

Notes#

Unless specified otherwise, the Keplerian elements are ordered as:

Index

Keplerian Element

0

Semi-major axis (except if eccentricity = 1.0, then represents semi-latus rectum)

1

Eccentricity

2

Inclination

3

Argument of periapsis

4

Longitude of ascending node

5

True anomaly

Unless specified otherwise, the spherical elements are ordered as:

Index

Spherical State Element

0

Radial distance

1

Latitude

2

Longitude

3

Speed

4

Flight-path angle

5

Heading angle

Unless specified otherwise, the Modified equinoctial elements (MEE; see here, element set k) are ordered as follows. The element \(I\), which defines the location of the MEE singularity is not treated as a state element, but is provided/determined separately:

Index

Modified equinoctial Element

0

Semi-parameter

1

f-element

2

g-element

3

h-element

4

k-element

5

True longitude

Functions#

cartesian_to_keplerian(cartesian_elements, ...)

Convert Cartesian to Keplerian elements.

keplerian_to_cartesian(keplerian_elements, ...)

Convert Keplerian elements to Cartesian.

keplerian_to_cartesian_elementwise(...)

Convert Keplerian elements to Cartesian, with elementwise input.

mean_to_true_anomaly(eccentricity, mean_anomaly)

Convert mean to true anomaly.

true_to_mean_anomaly(eccentricity, true_anomaly)

Convert true to mean anomaly.

true_to_eccentric_anomaly(true_anomaly, ...)

Convert true to eccentric anomaly.

eccentric_to_true_anomaly(eccentric_anomaly, ...)

Convert eccentric to true anomaly.

eccentric_to_mean_anomaly(eccentric_anomaly, ...)

Convert eccentric to mean anomaly.

mean_to_eccentric_anomaly(eccentricity, ...)

Convert mean to eccentric anomaly.

elapsed_time_to_delta_mean_anomaly(...)

Convert elapsed time to the corresponding change in mean anomaly along a Keplerian orbit.

delta_mean_anomaly_to_elapsed_time(...)

Convert change in mean anomaly along a Keplerian orbit to the corresponding elapsed time.

mean_motion_to_semi_major_axis(mean_motion, ...)

Convert mean motion to corresponding semi-major axis (in a Keplerian orbit).

semi_major_axis_to_mean_motion(...)

Convert semi-major axis to corresponding mean motion (along a Keplerian orbit).

keplerian_to_mee_manual_singularity(...)

Convert Keplerian to Modified equinoctial elements.

keplerian_to_mee(keplerian_elements)

Convert Keplerian to Modified equinoctial elements.

flip_mee_singularity(keplerian_elements)

Function to determine 'optimal' location of the singularity-flipping modified equinoctial element.

mee_to_keplerian(...)

Convert Modified equinoctial to Keplerian elements.

cartesian_to_mee(cartesian_elements, ...)

Convert Cartesian to Modified equinoctial elements.

cartesian_to_mee_manual_singularity(...)

Convert Cartesian to Modified equinoctial elements.

mee_to_cartesian(...)

Convert Modified equinoctial to Cartesian elements.

quaternion_entries_to_rotation_matrix(...)

Converts an array of four quaternion elements to the equivalent rotation matrix.

rotation_matrix_to_quaternion_entries(...)

Converts a rotation matrix to the equivalent array of four quaternion elements.

cartesian_to_spherical(cartesian_elements)

Convert Cartesian to spherical elements.

spherical_to_cartesian(spherical_elements)

Convert spherical elements to Cartesian.

spherical_to_cartesian_elementwise(...)

Convert Spherical elements to Cartesian, with elementwise input.

cartesian_to_keplerian(cartesian_elements: numpy.ndarray[numpy.float64[6, 1]], gravitational_parameter: float) numpy.ndarray[numpy.float64[6, 1]]#

Convert Cartesian to Keplerian elements.

Note

See module level documentation for the standard ordering convention of Keplerian elements used.

Parameters:
  • cartesian_elements (numpy.ndarray) – Cartesian state that is to be converted to Keplerian elements

  • gravitational_parameter (float) – Gravitational parameter of central body used for conversion

Returns:

Keplerian elements, as computed from Cartesian element input.

Return type:

numpy.ndarray

keplerian_to_cartesian(keplerian_elements: numpy.ndarray[numpy.float64[6, 1]], gravitational_parameter: float) numpy.ndarray[numpy.float64[6, 1]]#

Convert Keplerian elements to Cartesian.

Note

See module level documentation for the standard ordering convention of Keplerian elements used.

Parameters:
  • keplerian_elements (numpy.ndarray) – Keplerian state that is to be converted to Cartesian elements

  • gravitational_parameter (float) – Gravitational parameter of central body used for conversion

Returns:

Cartesian elements, as computed from Keplerian element input.

Return type:

numpy.ndarray

keplerian_to_cartesian_elementwise(semi_major_axis: float, eccentricity: float, inclination: float, argument_of_periapsis: float, longitude_of_ascending_node: float, true_anomaly: float, gravitational_parameter: float) numpy.ndarray[numpy.float64[6, 1]]#

Convert Keplerian elements to Cartesian, with elementwise input.

Note

The final Keplerian element is always the true anomaly.

Parameters:
  • semi_major_axis (float) – Semi-major axis (except if eccentricity = 1.0, then represents semi-latus rectum)

  • eccentricity (float) – Eccentricity

  • inclination (float) – Inclination

  • argument_of_periapsis (float) – Argument of periapsis

  • longitude_of_ascending_node (float) – Longitude of ascending node

  • true_anomaly (float) – True anomaly

  • gravitational_parameter (float) – Gravitational parameter of central body used for conversion

Returns:

Cartesian elements, as computed from Keplerian element input.

Return type:

numpy.ndarray

mean_to_true_anomaly(eccentricity: float, mean_anomaly: float, use_default_initial_guess: bool = True, non_default_initial_guess: float = nan, root_finder: tudatpy.kernel.math.root_finders.RootFinderCore = None) float#

Convert mean to true anomaly.

Convert the mean anomaly of the orbit to its true anomaly. This conversion first converts mean to eccentric anomaly (hyperbolic eccentric anomaly, if eccentricity is larger than 1, elliptical eccentric anomaly if it is smaller than 1), and subsequently to true anomaly.

Parameters:
  • eccentricity (float) – Value of the orbital eccentricity

  • mean_anomaly (float) – Value of the mean anomaly

  • use_default_initial_guess (bool, default = True) – Boolean to determine whether the user-defined initial guess (for mean-to-eccentric anomaly conversion) is used, or an automatically generated one.

  • non_default_initial_guess (float, default = NaN) – User-defined initial guess for mean-to-eccentric anomaly conversion, to be used only if use_default_initial_guess is set to False.

  • root_finder (RootFinder, default = None) – User-defined root finder, overriding default root-finding algorithm for mean-to-eccentric anomaly conversion (default is used if this input is left empty)

Returns:

Value of the true anomaly

Return type:

float

true_to_mean_anomaly(eccentricity: float, true_anomaly: float) float#

Convert true to mean anomaly.

Convert the true anomaly of the orbit to its mean anomaly. This conversion first converts true to eccentric anomaly (hyperbolic eccentric anomaly, if eccentricity is larger than 1, elliptical eccentric anomaly if it is smaller than 1), and subsequently to mean anomaly.

Parameters:
  • eccentricity (float) – Value of the orbital eccentricity

  • true_anomaly (float) – Value of the true anomaly

Returns:

Value of the mean anomaly

Return type:

float

true_to_eccentric_anomaly(true_anomaly: float, eccentricity: float) float#

Convert true to eccentric anomaly.

Parameters:
  • eccentricity (float) – Value of the orbital eccentricity

  • true_anomaly (float) – Value of the true anomaly

Returns:

Hyperbolic eccentric anomaly, if eccentricity is larger than 1, elliptical eccentric anomaly if it is smaller than 1

Return type:

float

eccentric_to_true_anomaly(eccentric_anomaly: float, eccentricity: float) float#

Convert eccentric to true anomaly.

Parameters:
  • eccentric_anomaly (float) – Hyperbolic eccentric anomaly, if eccentricity is larger than 1, elliptical eccentric anomaly if it is smaller than 1

  • eccentricity (float) – Value of the orbital eccentricity

Returns:

Value of the true anomaly

Return type:

float

eccentric_to_mean_anomaly(eccentric_anomaly: float, eccentricity: float) float#

Convert eccentric to mean anomaly.

Parameters:
  • eccentric_anomaly (float) – Hyperbolic eccentric anomaly, if eccentricity is larger than 1, elliptical eccentric anomaly if it is smaller than 1

  • eccentricity (float) – Value of the orbital eccentricity

Returns:

Value of the mean anomaly

Return type:

float

mean_to_eccentric_anomaly(eccentricity: float, mean_anomaly: float, use_default_initial_guess: bool = True, non_default_initial_guess: float = nan, root_finder: tudatpy.kernel.math.root_finders.RootFinderCore = None) float#

Convert mean to eccentric anomaly.

Parameters:
  • eccentricity (float) – Value of the orbital eccentricity

  • mean_anomaly (float) – Value of the mean anomaly

  • use_default_initial_guess (bool, default = True) – Boolean to determine whether the user-defined initial guess is used for conversion, or an automatically generated one.

  • non_default_initial_guess (float, default = NaN) – User-defined initial guess for conversion, to be used only if use_default_initial_guess is set to False.

  • root_finder (RootFinder, default = None) – User-defined root finder, overriding default root-finding algorithm for conversion (default is used if this input is left empty)

Returns:

Value of the eccentric anomaly

Return type:

float

elapsed_time_to_delta_mean_anomaly(elapsed_time: float, gravitational_parameter: float, semi_major_axis: float) float#

Convert elapsed time to the corresponding change in mean anomaly along a Keplerian orbit.

Parameters:
  • elapsed_time (float) – Elapsed time (in seconds)

  • gravitational_parameter (float) – Gravitational parameter of central body

  • semi_major_axis (float) – Semi-major axis of orbit

Returns:

Total change in mean anomaly along the Kepler orbit, accumulated in the provided time.

Return type:

float

delta_mean_anomaly_to_elapsed_time(mean_anomaly_change: float, gravitational_parameter: float, semi_major_axis: float) float#

Convert change in mean anomaly along a Keplerian orbit to the corresponding elapsed time.

Parameters:
  • mean_anomaly_change (float) – Total change in mean anomaly along the Kepler orbit

  • gravitational_parameter (float) – Gravitational parameter of central body

  • semi_major_axis (float) – Semi-major axis of orbit

Returns:

Time required for the provided mean anomaly change to be accumulated

Return type:

float

mean_motion_to_semi_major_axis(mean_motion: float, gravitational_parameter: float) float#

Convert mean motion to corresponding semi-major axis (in a Keplerian orbit).

Parameters:
  • mean_motion (float) – Orbital mean motion

  • gravitational_parameter (float) – Gravitational parameter of central body

Returns:

Semi-major axis corresponding to mean motion

Return type:

float

semi_major_axis_to_mean_motion(semi_major_axis: float, gravitational_parameter: float) float#

Convert semi-major axis to corresponding mean motion (along a Keplerian orbit).

Parameters:
  • semi_major_axis (float) – Semi-major axis of orbit

  • gravitational_parameter (float) – Gravitational parameter of central body

Returns:

Semi-major axis corresponding to mean motion

Return type:

float

keplerian_to_mee_manual_singularity(keplerian_elements: numpy.ndarray[numpy.float64[6, 1]], singularity_at_zero_inclination: bool) numpy.ndarray[numpy.float64[6, 1]]#

Convert Keplerian to Modified equinoctial elements.

Convert Keplerian to Modified equinoctial elements (without intermediate step to Cartesian elements). The singularity-flipping element \(I\) is to be provided manually for this function

Note

See module level documentation for the standard ordering convention of Modified Equinoctial elements used.

Parameters:
  • keplerian_elements (numpy.ndarray) – Keplerian elements that are to be converted to Modified equinoctial elements

  • singularity_at_zero_inclination (bool) – Singularity at 0 degrees inclination if True, 180 degrees if False

Returns:

Modified equinoctial elements, as computed from Keplerian element input.

Return type:

numpy.ndarray

keplerian_to_mee(keplerian_elements: numpy.ndarray[numpy.float64[6, 1]]) numpy.ndarray[numpy.float64[6, 1]]#

Convert Keplerian to Modified equinoctial elements.

Convert Keplerian to Modified equinoctial elements (without intermediate step to Cartesian elements). The singularity-flipping element \(I\) is computed automatically by this function (using flip_mee_singularity())

Note

See module level documentation for the standard ordering convention of Modified Equinoctial elements used.

Parameters:

keplerian_elements (numpy.ndarray) – Keplerian elements that are to be converted to Modified equinoctial elements

Returns:

Modified equinoctial elements, as computed from Keplerian element input (with element \(I\) defined by flip_mee_singularity()).

Return type:

numpy.ndarray

flip_mee_singularity(keplerian_elements: numpy.ndarray[numpy.float64[6, 1]]) bool#

Function to determine ‘optimal’ location of the singularity-flipping modified equinoctial element.

Function to determine ‘optimal’ location of the singularity-flipping modified equinoctial element \(I\), if orbit inclination is less than 90 degrees, it puts the singularity at 180 degrees, if it is larger than 90 degrees, it puts it at 0 degrees.

Parameters:

keplerian_elements (numpy.ndarray) – Keplerian elements that are to be converted to Modified equinoctial elements

Returns:

Singularity at 0 degrees inclination if false, 180 degrees if true

Return type:

bool

mee_to_keplerian(modified_equinoctial_elements: numpy.ndarray[numpy.float64[6, 1]], singularity_at_zero_inclination: bool) numpy.ndarray[numpy.float64[6, 1]]#

Convert Modified equinoctial to Keplerian elements.

Modified equinoctial elements to Keplerian (without intermediate step to Cartesian elements).

Note

See module level documentation for the standard ordering convention of Modified Equinoctial elements used.

Parameters:
  • modified_equinoctial_elements (numpy.ndarray) – Modified equinoctial elements that are to be converted to Keplerian elements

  • singularity_at_zero_inclination (bool) – Singularity at 0 degrees inclination if false, 180 degrees if true

Returns:

Keplerian elements, as computed from Modified equinoctial element input.

Return type:

numpy.ndarray

cartesian_to_mee(cartesian_elements: numpy.ndarray[numpy.float64[6, 1]], gravitational_parameter: float) numpy.ndarray[numpy.float64[6, 1]]#

Convert Cartesian to Modified equinoctial elements.

Convert cartesian to Modified equinoctial elements. The singularity-flipping element \(I\) is computed automatically by this function (using flip_mee_singularity())

Note

See module level documentation for the standard ordering convention of Modified Equinoctial elements used.

Parameters:
  • cartesian_elements (numpy.ndarray) – Cartesian elements that are to be converted to Modified equinoctial elements

  • gravitational_parameter (float) – Gravitational parameter of central body

Returns:

Modified equinoctial elements, as computed from Cartesian element input.

Return type:

numpy.ndarray

cartesian_to_mee_manual_singularity(cartesian_elements: numpy.ndarray[numpy.float64[6, 1]], gravitational_parameter: float, singularity_at_zero_inclination: bool) numpy.ndarray[numpy.float64[6, 1]]#

Convert Cartesian to Modified equinoctial elements.

Convert cartesian to Modified equinoctial elements. The singularity-flipping element \(I\) is to be provided manually for this function

Note

See module level documentation for the standard ordering convention of Modified Equinoctial elements used.

Parameters:
  • cartesian_elements (numpy.ndarray) – Cartesian elements that are to be converted to Modified equinoctial elements

  • gravitational_parameter (float) – Gravitational parameter of central body

  • singularity_at_zero_inclination (bool) – Singularity at 0 degrees inclination if false, 180 degrees if true

Returns:

Modified equinoctial elements, as computed from Cartesian element input.

Return type:

numpy.ndarray

mee_to_cartesian(modified_equinoctial_elements: numpy.ndarray[numpy.float64[6, 1]], gravitational_parameter: float, singularity_at_zero_inclination: bool) numpy.ndarray[numpy.float64[6, 1]]#

Convert Modified equinoctial to Cartesian elements.

Note

See module level documentation for the standard ordering convention of Modified Equinoctial elements used.

Parameters:
  • modified_equinoctial_elements (numpy.ndarray) – Modified equinoctial elements that are to be converted to Cartesian elements

  • gravitational_parameter (float) – Gravitational parameter of central body

  • singularity_at_zero_inclination (bool) – Singularity at 0 degrees inclination if false, 180 degrees if true

Returns:

Cartesian elements, as computed from Modified equinoctial element input.

Return type:

numpy.ndarray

quaternion_entries_to_rotation_matrix(quaternion_entries: numpy.ndarray[numpy.float64[4, 1]]) numpy.ndarray[numpy.float64[3, 3]]#

Converts an array of four quaternion elements to the equivalent rotation matrix.

Function to convert an array of four quaternion elements to the equivalent rotation matrix. These quaternion elements are for instance used when propagating rotational dynamics in Tudat, and this function can be used to convert the numerical results to a usable rotation matrix. See our user guide for more details.

Parameters:

quaternion_entries (numpy.ndarray) – Quaternion elements, as per the convention used in the Eigen library

Returns:

Rotation matrix defining the equivalent rotation.

Return type:

numpy.ndarray

rotation_matrix_to_quaternion_entries(rotation_matrix: numpy.ndarray[numpy.float64[3, 3]]) numpy.ndarray[numpy.float64[4, 1]]#

Converts a rotation matrix to the equivalent array of four quaternion elements.

Inverse function of quaternion_entries_to_rotation_matrix().

Parameters:

rotation_matrix (numpy.ndarray) – Rotation matrix

Returns:

Equivalent quaternion elements, as per the convention used in the Eigen library

Return type:

numpy.ndarray

cartesian_to_spherical(cartesian_elements: numpy.ndarray[numpy.float64[6, 1]]) numpy.ndarray[numpy.float64[6, 1]]#

Convert Cartesian to spherical elements.

Note

See module level documentation for the standard ordering convention of spherical state elements used.

Parameters:

cartesian_elements (numpy.ndarray) – Cartesian state that is to be converted to spherical elements

Returns:

Spherical elements, as computed from Cartesian element input.

Return type:

numpy.ndarray

spherical_to_cartesian(spherical_elements: numpy.ndarray[numpy.float64[6, 1]]) numpy.ndarray[numpy.float64[6, 1]]#

Convert spherical elements to Cartesian.

Note

See module level documentation for the standard ordering convention of spherical state elements used.

Parameters:

spherical_elements (numpy.ndarray) – Spherical state that is to be converted to Cartesian elements

Returns:

Cartesian elements, as computed from spherical element input.

Return type:

numpy.ndarray

spherical_to_cartesian_elementwise(radial_distance: float, latitude: float, longitude: float, speed: float, flight_path_angle: float, heading_angle: float) numpy.ndarray[numpy.float64[6, 1]]#

Convert Spherical elements to Cartesian, with elementwise input.

Parameters:
  • radial_distance (float) – Distance from origin of central body

  • latitude (float) – Central body-fixed latitude

  • longitude (float) – Central body-fixed longitude

  • speed (float) – Central body-fixed speed (norm of velocity vector). Note that this is not the norm of the inertial velocity

  • flight_path_angle (float) – Flight-path angle (of central body-fixed velocity vector)

  • heading_angle (float) – Heading angle (of central body-fixed velocity vector)

Returns:

Cartesian elements, as computed from spherical element input.

Return type:

numpy.ndarray