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#
In general, unless specified otherwise, the Keplerian elements are ordered as:
Index |
Keplerian Element |
|
Semi-major axis (except if eccentricity = |
|
Eccentricity |
|
Inclination |
|
Argument of periapsis |
|
Longitude of ascending node |
|
True anomaly |
Functions#
|
Convert Cartesian to Keplerian elements. |
|
Convert Keplerian elements to Cartesian. |
Convert Keplerian elements to Cartesian, with elementwise input. |
|
Converts an array of four quaternion elements to the equivalent rotation matrix. |
|
Converts a rotation matrix to the equivalent array of four quaternion elements. |
- 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
- 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
Keplerian elements, as computed from Cartesian element input.
- Return type
- 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
Keplerian elements, as computed from Cartesian element input.
- Return type
- 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
- 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