keplerian_sim

Trait OrbitTrait

Source
pub trait OrbitTrait {
Show 33 methods // Required methods fn get_semi_major_axis(&self) -> f64; fn get_semi_minor_axis(&self) -> f64; fn get_semi_latus_rectum(&self) -> f64; fn get_linear_eccentricity(&self) -> f64; fn get_apoapsis(&self) -> f64; fn set_apoapsis(&mut self, apoapsis: f64) -> Result<(), ApoapsisSetterError>; fn set_apoapsis_force(&mut self, apoapsis: f64); fn get_transformation_matrix(&self) -> Matrix3x2<f64>; fn get_eccentric_anomaly(&self, mean_anomaly: f64) -> f64; fn get_true_anomaly_at_eccentric_anomaly( &self, eccentric_anomaly: f64, ) -> f64; fn get_mean_anomaly_at_time(&self, t: f64) -> f64; fn get_altitude_at_angle(&self, angle: f64) -> f64; fn get_eccentricity(&self) -> f64; fn set_eccentricity(&mut self, eccentricity: f64); fn get_periapsis(&self) -> f64; fn set_periapsis(&mut self, periapsis: f64); fn get_inclination(&self) -> f64; fn set_inclination(&mut self, inclination: f64); fn get_arg_pe(&self) -> f64; fn set_arg_pe(&mut self, arg_pe: f64); fn get_long_asc_node(&self) -> f64; fn set_long_asc_node(&mut self, long_asc_node: f64); fn get_mean_anomaly_at_epoch(&self) -> f64; fn set_mean_anomaly_at_epoch(&mut self, mean_anomaly: f64); // Provided methods fn get_true_anomaly(&self, mean_anomaly: f64) -> f64 { ... } fn get_eccentric_anomaly_at_time(&self, t: f64) -> f64 { ... } fn get_true_anomaly_at_time(&self, t: f64) -> f64 { ... } fn get_position_at_angle(&self, angle: f64) -> (f64, f64, f64) { ... } fn get_flat_position_at_angle(&self, angle: f64) -> (f64, f64) { ... } fn get_altitude_at_time(&self, t: f64) -> f64 { ... } fn get_position_at_time(&self, t: f64) -> (f64, f64, f64) { ... } fn get_flat_position_at_time(&self, t: f64) -> (f64, f64) { ... } fn tilt_flat_position(&self, x: f64, y: f64) -> (f64, f64, f64) { ... }
}
Expand description

A trait that defines the methods that a Keplerian orbit must implement.

This trait is implemented by both Orbit and CompactOrbit.

§Examples

use keplerian_sim::{Orbit, OrbitTrait, CompactOrbit};
 
fn accepts_orbit(orbit: &impl OrbitTrait) {
    println!("That's an orbit!");
}
 
fn main() {
    let orbit = Orbit::new_default();
    accepts_orbit(&orbit);
 
    let compact = CompactOrbit::new_default();
    accepts_orbit(&compact);
}

This example will fail to compile:

let not_orbit = (0.0, 1.0);
accepts_orbit(&not_orbit);

Required Methods§

Source

fn get_semi_major_axis(&self) -> f64

Gets the semi-major axis of the orbit.

In an elliptic orbit, the semi-major axis is the average of the apoapsis and periapsis.
This function uses a generalization which uses eccentricity instead.

This function returns infinity for parabolic orbits, and negative values for hyperbolic orbits.

Learn more: https://en.wikipedia.org/wiki/Semi-major_and_semi-minor_axes

§Example
use keplerian_sim::{Orbit, OrbitTrait};
 
let mut orbit = Orbit::new_default();
orbit.set_periapsis(50.0);
orbit.set_apoapsis_force(100.0);
let sma = orbit.get_semi_major_axis();
let expected = 75.0;
assert!((sma - expected).abs() < 1e-6);
Source

fn get_semi_minor_axis(&self) -> f64

Gets the semi-minor axis of the orbit.

In an elliptic orbit, the semi-minor axis is half of the maximum “width” of the orbit.

Learn more: https://en.wikipedia.org/wiki/Semi-major_and_semi-minor_axes

Source

fn get_semi_latus_rectum(&self) -> f64

Source

fn get_linear_eccentricity(&self) -> f64

Gets the linear eccentricity of the orbit, in meters.

In an elliptic orbit, the linear eccentricity is the distance between its center and either of its two foci (focuses).

§Example
use keplerian_sim::{Orbit, OrbitTrait};
 
let mut orbit = Orbit::new_default();
orbit.set_periapsis(50.0);
orbit.set_apoapsis_force(100.0);
 
// Let's say the periapsis is at x = -50.
// The apoapsis would be at x = 100.
// The midpoint would be at x = 25.
// The parent body - one of its foci - is always at the origin (x = 0).
// This means the linear eccentricity is 25.
 
let linear_eccentricity = orbit.get_linear_eccentricity();
let expected = 25.0;
 
assert!((linear_eccentricity - expected).abs() < 1e-6);
Source

fn get_apoapsis(&self) -> f64

Gets the apoapsis of the orbit.
Returns infinity for parabolic orbits.
Returns negative values for hyperbolic orbits.

§Examples
use keplerian_sim::{Orbit, OrbitTrait};
 
let mut orbit = Orbit::new_default();
orbit.set_eccentricity(0.5); // Elliptic
assert!(orbit.get_apoapsis() > 0.0);
 
orbit.set_eccentricity(1.0); // Parabolic
assert!(orbit.get_apoapsis().is_infinite());
 
orbit.set_eccentricity(2.0); // Hyperbolic
assert!(orbit.get_apoapsis() < 0.0);
Source

fn set_apoapsis(&mut self, apoapsis: f64) -> Result<(), ApoapsisSetterError>

Sets the apoapsis of the orbit.
Errors when the apoapsis is less than the periapsis, or less than zero.
If you want a setter that does not error, use set_apoapsis_force, which will try its best to interpret what you might have meant, but may have undesirable behavior.

§Examples
use keplerian_sim::{Orbit, OrbitTrait};
 
let mut orbit = Orbit::new_default();
orbit.set_periapsis(50.0);
 
assert!(
    orbit.set_apoapsis(100.0)
        .is_ok()
);
 
let result = orbit.set_apoapsis(25.0);
assert!(result.is_err());
assert!(
    result.unwrap_err() ==
    keplerian_sim::ApoapsisSetterError::ApoapsisLessThanPeriapsis
);
 
let result = orbit.set_apoapsis(-25.0);
assert!(result.is_err());
assert!(
    result.unwrap_err() ==
    keplerian_sim::ApoapsisSetterError::ApoapsisNegative
);
Source

fn set_apoapsis_force(&mut self, apoapsis: f64)

Sets the apoapsis of the orbit, with a best-effort attempt at interpreting possibly-invalid values.
This function will not error, but may have undesirable behavior:

  • If the given apoapsis is less than the periapsis but more than zero, the orbit will be flipped and the periapsis will be set to the given apoapsis.
  • If the given apoapsis is less than zero, the orbit will be hyperbolic instead.

If these behaviors are undesirable, consider creating a custom wrapper around set_eccentricity instead.

Source

fn get_transformation_matrix(&self) -> Matrix3x2<f64>

Gets the transformation matrix needed to tilt a 2D vector into the tilted orbital plane.

§Example
use keplerian_sim::{Orbit, OrbitTrait};
 
let orbit = Orbit::new_default();
let matrix = orbit.get_transformation_matrix();
 
assert_eq!(matrix, keplerian_sim::Matrix3x2 {
    e11: 1.0, e12: 0.0,
    e21: 0.0, e22: 1.0,
    e31: 0.0, e32: 0.0,
});
Source

fn get_eccentric_anomaly(&self, mean_anomaly: f64) -> f64

Gets the eccentric anomaly at a given mean anomaly in the orbit.

The method to get the eccentric anomaly often uses numerical methods like Newton’s method, and so it is not very performant.
It is recommended to cache this value if you can.

When the orbit is open (has an eccentricity of at least 1), the hyperbolic eccentric anomaly would be returned instead.

The eccentric anomaly is an angular parameter that defines the position of a body that is moving along an elliptic Kepler orbit.

- Wikipedia

Source

fn get_true_anomaly_at_eccentric_anomaly(&self, eccentric_anomaly: f64) -> f64

Gets the true anomaly at a given eccentric anomaly in the orbit.

This function is faster than the function which takes mean anomaly as input, as the eccentric anomaly is hard to calculate.

Source

fn get_mean_anomaly_at_time(&self, t: f64) -> f64

Gets the mean anomaly at a given time in the orbit.

The mean anomaly is the fraction of an elliptical orbit’s period that has elapsed since the orbiting body passed periapsis, expressed as an angle which can be used in calculating the position of that body in the classical two-body problem.

- Wikipedia

Source

fn get_altitude_at_angle(&self, angle: f64) -> f64

Gets the altitude of the body from its parent at a given angle (true anomaly) in the orbit.

§Example
use keplerian_sim::{Orbit, OrbitTrait};
 
let mut orbit = Orbit::new_default();
orbit.set_periapsis(100.0);
orbit.set_eccentricity(0.0);
 
let altitude = orbit.get_altitude_at_angle(0.0);
 
assert_eq!(altitude, 100.0);
Source

fn get_eccentricity(&self) -> f64

Gets the eccentricity of the orbit.

The eccentricity of an orbit is a measure of how much it deviates from a perfect circle.

An eccentricity of 0 means the orbit is a perfect circle.
Between 0 and 1, the orbit is elliptic, and has an oval shape.
An orbit with an eccentricity of 1 is said to be parabolic.
If it’s greater than 1, the orbit is hyperbolic.

For hyperbolic trajectories, the higher the eccentricity, the straighter the path.

Wikipedia on conic section eccentricity: https://en.wikipedia.org/wiki/Eccentricity_(mathematics)
(Keplerian orbits are conic sections, so the concepts still apply)

Source

fn set_eccentricity(&mut self, eccentricity: f64)

Sets the eccentricity of the orbit.

The eccentricity of an orbit is a measure of how much it deviates from a perfect circle.

An eccentricity of 0 means the orbit is a perfect circle.
Between 0 and 1, the orbit is elliptic, and has an oval shape.
An orbit with an eccentricity of 1 is said to be parabolic.
If it’s greater than 1, the orbit is hyperbolic.

For hyperbolic trajectories, the higher the eccentricity, the straighter the path.

Wikipedia on conic section eccentricity: https://en.wikipedia.org/wiki/Eccentricity_(mathematics)
(Keplerian orbits are conic sections, so the concepts still apply)

Source

fn get_periapsis(&self) -> f64

Gets the periapsis of the orbit.

The periapsis of an orbit is the distance at the closest point to the parent body.

More simply, this is the “minimum altitude” of an orbit.

Wikipedia: https://en.wikipedia.org/wiki/Apsis

Source

fn set_periapsis(&mut self, periapsis: f64)

Sets the periapsis of the orbit.

The periapsis of an orbit is the distance at the closest point to the parent body.

More simply, this is the “minimum altitude” of an orbit.

Wikipedia: https://en.wikipedia.org/wiki/Apsis

Source

fn get_inclination(&self) -> f64

Gets the inclination of the orbit in radians.

The inclination of an orbit is the angle between the plane of the orbit and the reference plane.

In simple terms, it tells you how “tilted” the orbit is.

Wikipedia: https://en.wikipedia.org/wiki/Orbital_inclination

Source

fn set_inclination(&mut self, inclination: f64)

Sets the inclination of the orbit in radians.

The inclination of an orbit is the angle between the plane of the orbit and the reference plane.

In simple terms, it tells you how “tilted” the orbit is.

Wikipedia: https://en.wikipedia.org/wiki/Orbital_inclination

Source

fn get_arg_pe(&self) -> f64

Gets the argument of periapsis of the orbit in radians.

Wikipedia:
The argument of periapsis is the angle from the body’s ascending node to its periapsis, measured in the direction of motion.
https://en.wikipedia.org/wiki/Argument_of_periapsis

In simple terms, it tells you how, and in which direction, the orbit “tilts”.

Source

fn set_arg_pe(&mut self, arg_pe: f64)

Sets the argument of periapsis of the orbit in radians.

Wikipedia:
The argument of periapsis is the angle from the body’s ascending node to its periapsis, measured in the direction of motion.
https://en.wikipedia.org/wiki/Argument_of_periapsis

In simple terms, it tells you how, and in which direction, the orbit “tilts”.

Source

fn get_long_asc_node(&self) -> f64

Gets the longitude of ascending node of the orbit in radians.

Wikipedia:
The longitude of ascending node is the angle from a specified reference direction, called the origin of longitude, to the direction of the ascending node, as measured in a specified reference plane.
https://en.wikipedia.org/wiki/Longitude_of_the_ascending_node

In simple terms, it tells you how, and in which direction, the orbit “tilts”.

Source

fn set_long_asc_node(&mut self, long_asc_node: f64)

Sets the longitude of ascending node of the orbit in radians.

Wikipedia:
The longitude of ascending node is the angle from a specified reference direction, called the origin of longitude, to the direction of the ascending node, as measured in a specified reference plane.
https://en.wikipedia.org/wiki/Longitude_of_the_ascending_node

In simple terms, it tells you how, and in which direction, the orbit “tilts”.

Source

fn get_mean_anomaly_at_epoch(&self) -> f64

Gets the mean anomaly of the orbit at a certain epoch.

For elliptic orbits, it’s measured in radians and so are bounded between 0 and tau; anything out of range will get wrapped around.
For hyperbolic orbits, it’s unbounded.

Wikipedia:
The mean anomaly at epoch, M_0, is defined as the instantaneous mean anomaly at a given epoch, t_0.
https://en.wikipedia.org/wiki/Mean_anomaly#Mean_anomaly_at_epoch

In simple terms, this modifies the “offset” of the orbit progression.

Source

fn set_mean_anomaly_at_epoch(&mut self, mean_anomaly: f64)

Sets the mean anomaly of the orbit at a certain epoch.

For elliptic orbits, it’s measured in radians and so are bounded between 0 and tau; anything out of range will get wrapped around.
For hyperbolic orbits, it’s unbounded.

Wikipedia:
The mean anomaly at epoch, M_0, is defined as the instantaneous mean anomaly at a given epoch, t_0.
https://en.wikipedia.org/wiki/Mean_anomaly#Mean_anomaly_at_epoch

In simple terms, this modifies the “offset” of the orbit progression.

Provided Methods§

Source

fn get_true_anomaly(&self, mean_anomaly: f64) -> f64

Gets the true anomaly at a given mean anomaly in the orbit.

The true anomaly is derived from the eccentric anomaly, which uses numerical methods and so is not very performant.
It is recommended to cache this value if you can.

The true anomaly is the angle between the direction of periapsis and the current position of the body, as seen from the main focus of the ellipse.

- Wikipedia

Source

fn get_eccentric_anomaly_at_time(&self, t: f64) -> f64

Gets the eccentric anomaly at a given time in the orbit.

The method to get the eccentric anomaly often uses numerical methods like Newton’s method, and so it is not very performant.
It is recommended to cache this value if you can.

When the orbit is open (has an eccentricity of at least 1), the hyperbolic eccentric anomaly would be returned instead.

The eccentric anomaly is an angular parameter that defines the position of a body that is moving along an elliptic Kepler orbit.

- Wikipedia

Source

fn get_true_anomaly_at_time(&self, t: f64) -> f64

Gets the true anomaly at a given time in the orbit.

The true anomaly is derived from the eccentric anomaly, which uses numerical methods and so is not very performant.
It is recommended to cache this value if you can.

The true anomaly is the angle between the direction of periapsis and the current position of the body, as seen from the main focus of the ellipse.

- Wikipedia

Source

fn get_position_at_angle(&self, angle: f64) -> (f64, f64, f64)

Gets the 3D position at a given angle (true anomaly) in the orbit.

The angle is expressed in radians, and ranges from 0 to tau.
Anything out of range will get wrapped around.

§Example
use keplerian_sim::{Orbit, OrbitTrait};
 
let mut orbit = Orbit::new_default();
orbit.set_periapsis(100.0);
orbit.set_eccentricity(0.0);
 
let pos = orbit.get_position_at_angle(0.0);
 
assert_eq!(pos, (100.0, 0.0, 0.0));
Source

fn get_flat_position_at_angle(&self, angle: f64) -> (f64, f64)

Gets the 2D position at a given angle (true anomaly) in the orbit.

This ignores “orbital tilting” parameters, namely the inclination and the longitude of ascending node.

The angle is expressed in radians, and ranges from 0 to tau.
Anything out of range will get wrapped around.

§Example
use keplerian_sim::{Orbit, OrbitTrait};
 
let mut orbit = Orbit::new_default();
orbit.set_periapsis(100.0);
orbit.set_eccentricity(0.0);
 
let pos = orbit.get_flat_position_at_angle(0.0);
 
assert_eq!(pos, (100.0, 0.0));
Source

fn get_altitude_at_time(&self, t: f64) -> f64

Gets the altitude of the body from its parent at a given time in the orbit.

This involves calculating the true anomaly at a given time, and so is not very performant.
It is recommended to cache this value when possible.

For closed orbits (with an eccentricity less than 1), the t (time) value ranges from 0 to 1.
Anything out of range will get wrapped around.

For open orbits (with an eccentricity of at least 1), the t (time) value is unbounded.
Note that due to floating-point imprecision, values of extreme magnitude may not be accurate.

Source

fn get_position_at_time(&self, t: f64) -> (f64, f64, f64)

Gets the 3D position at a given time in the orbit.

This involves calculating the true anomaly at a given time, and so is not very performant.
It is recommended to cache this value when possible.

For closed orbits (with an eccentricity less than 1), the t (time) value ranges from 0 to 1.
Anything out of range will get wrapped around.

For open orbits (with an eccentricity of at least 1), the t (time) value is unbounded.
Note that due to floating-point imprecision, values of extreme magnitude may not be accurate.

Source

fn get_flat_position_at_time(&self, t: f64) -> (f64, f64)

Gets the 2D position at a given time in the orbit.

This involves calculating the true anomaly at a given time, and so is not very performant. It is recommended to cache this value when possible.

This ignores “orbital tilting” parameters, namely the inclination and longitude of ascending node.

For closed orbits (with an eccentricity less than 1), the t (time) value ranges from 0 to 1.
Anything out of range will get wrapped around.

For open orbits (with an eccentricity of at least 1), the t (time) value is unbounded.
Note that due to floating-point imprecision, values of extreme magnitude may not be accurate.

Source

fn tilt_flat_position(&self, x: f64, y: f64) -> (f64, f64, f64)

Tilts a 2D position into 3D, using the orbital parameters.

This uses the “orbital tilting” parameters, namely the inclination and longitude of ascending node, to tilt that position into the same plane that the orbit resides in.

This function performs 10x faster in the cached version of the Orbit struct, as it doesn’t need to recalculate the transformation matrix needed to transform 2D vector.

Implementors§