Skip to main content

CompactOrbit2D

Struct CompactOrbit2D 

Source
pub struct CompactOrbit2D {
    pub eccentricity: f64,
    pub periapsis: f64,
    pub arg_pe: f64,
    pub mean_anomaly: f64,
    pub mu: f64,
}
Expand description

A minimal struct representing a 2D Keplerian orbit.

This struct minimizes memory footprint by not caching variables. Because of this, calculations can be slightly slower than the cached version of the 2D orbit. For this reason, you may consider using the Orbit2D struct instead.

Note that in 2D, the difference in performance between the cached and compact orbits is a lot less noticeable than in 3D.

§Example

use keplerian_sim::{CompactOrbit2D, OrbitTrait2D};

let orbit = CompactOrbit2D::new(
    // Initialize using eccentricity, periapsis,
    // argument of periapsis, mean anomaly at epoch,
    // and gravitational parameter

    // Eccentricity
    0.0,

    // Periapsis
    1.0,

    // Argument of periapsis
    0.0,

    // Mean anomaly at epoch
    0.0,

    // Gravitational parameter of the parent body
    1.0,
);

let orbit = CompactOrbit2D::with_apoapsis(
    // Initialize using apoapsis in place of eccentricity

    // Apoapsis
    2.0,

    // Periapsis
    1.0,

    // Argument of periapsis
    0.0,

    // Mean anomaly at epoch
    0.0,

    // Gravitational parameter of the parent body
    1.0,
);

See CompactOrbit2D::new and CompactOrbit2D::with_apoapsis for more information.

Fields§

§eccentricity: f64

The eccentricity of the orbit.
e < 1: ellipse
e = 1: parabola
e > 1: hyperbola\

See more: https://en.wikipedia.org/wiki/Orbital_eccentricity

§periapsis: f64

The periapsis of the orbit, in meters.

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.

§arg_pe: f64

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”.

§mean_anomaly: f64

The mean anomaly at orbit epoch, in radians.

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.

§mu: f64

The gravitational parameter of the parent body.

This is a constant value that represents the mass of the parent body multiplied by the gravitational constant.

In other words, mu = GM.

Implementations§

Source§

impl CompactOrbit2D

Source

pub fn new( eccentricity: f64, periapsis: f64, arg_pe: f64, mean_anomaly: f64, mu: f64, ) -> Self

Creates a new CompactOrbit2D instance with the given parameters.

Note: This function uses eccentricity instead of apoapsis. If you want to provide an apoapsis instead, consider using the CompactOrbit2D::with_apoapsis function instead.

§Parameters
  • eccentricity: The eccentricity of the orbit.
  • periapsis: The periapsis of the orbit, in meters.
  • arg_pe: The argument of periapsis of the orbit, in radians.
  • mean_anomaly: The mean anomaly of the orbit at epoch, in radians.
  • mu: The gravitational parameter of the parent body, in m^3 s^-2.
§Example
use keplerian_sim::{CompactOrbit2D, OrbitTrait2D};

let eccentricity = 0.2;
let periapsis = 2.8;
let argument_of_periapsis = 0.7;
let mean_anomaly_at_epoch = 2.9;
let gravitational_parameter = 9.2;

let orbit = CompactOrbit2D::new(
    eccentricity,
    periapsis,
    argument_of_periapsis,
    mean_anomaly_at_epoch,
    gravitational_parameter,
);

assert_eq!(orbit.eccentricity, eccentricity);
assert_eq!(orbit.periapsis, periapsis);
assert_eq!(orbit.arg_pe, argument_of_periapsis);
assert_eq!(orbit.mean_anomaly, mean_anomaly_at_epoch);
assert_eq!(orbit.mu, gravitational_parameter);
Source

pub fn with_apoapsis( apoapsis: f64, periapsis: f64, arg_pe: f64, mean_anomaly: f64, mu: f64, ) -> Self

Creates a new CompactOrbit2D instance with the given parameters.

Note: This function uses apoapsis instead of eccentricity.
Because of this, it’s not recommended to create parabolic or hyperbolic trajectories with this function.
If you’re looking to initialize a parabolic or hyperbolic trajectory, consider using the CompactOrbit2D::new function instead.

§Parameters
  • apoapsis: The apoapsis of the orbit, in meters. Must be more than the periapsis.
  • periapsis: The periapsis of the orbit, in meters.
  • arg_pe: The argument of periapsis of the orbit, in radians.
  • mean_anomaly: The mean anomaly of the orbit at epoch, in radians.
  • mu: The gravitational parameter of the parent body, in m^3 s^-2.
§Example
use keplerian_sim::{CompactOrbit2D, OrbitTrait2D};

let apoapsis = 4.1;
let periapsis = 2.8;
let argument_of_periapsis = 0.7;
let mean_anomaly_at_epoch = 2.9;
let gravitational_parameter = 9.2;

let orbit = CompactOrbit2D::with_apoapsis(
    apoapsis,
    periapsis,
    argument_of_periapsis,
    mean_anomaly_at_epoch,
    gravitational_parameter
);

let eccentricity = (apoapsis - periapsis) / (apoapsis + periapsis);

assert_eq!(orbit.eccentricity, eccentricity);
assert_eq!(orbit.periapsis, periapsis);
assert_eq!(orbit.arg_pe, argument_of_periapsis);
assert_eq!(orbit.mean_anomaly, mean_anomaly_at_epoch);
assert_eq!(orbit.mu, gravitational_parameter);
Source

pub fn new_circular(radius: f64, mean_anomaly: f64, mu: f64) -> Self

Creates a new circular CompactOrbit2D instance with the given parameters.

§Parameters
  • radius: The radius of the orbit, in meters.
  • mean_anomaly: The mean anomaly of hte orbit at epoch, in radians.
  • mu: The gravitational parameter of the parent body, in m^3 s^-2.
§Example
use keplerian_sim::{CompactOrbit2D, OrbitTrait2D};

let radius = 4.2;
let mean_anomaly_at_epoch = 1.5;
let gravitational_parameter = 5.0;

let orbit = CompactOrbit2D::new_circular(
    radius,
    mean_anomaly_at_epoch,
    gravitational_parameter,
);

assert_eq!(orbit.eccentricity, 0.0);
assert_eq!(orbit.periapsis, radius);
assert_eq!(orbit.arg_pe, 0.0);
assert_eq!(orbit.mean_anomaly, mean_anomaly_at_epoch);
assert_eq!(orbit.mu, gravitational_parameter);

Trait Implementations§

Source§

impl Clone for CompactOrbit2D

Source§

fn clone(&self) -> CompactOrbit2D

Returns a duplicate of the value. Read more
1.0.0 · Source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
Source§

impl Debug for CompactOrbit2D

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
Source§

impl Default for CompactOrbit2D

Source§

fn default() -> Self

Creates a unit orbit.

The unit orbit is a perfect circle of radius 1 and zero mean anomaly at epoch.

It also uses a gravitational parameter of 1.

Source§

impl From<CompactOrbit2D> for Orbit2D

Source§

fn from(compact: CompactOrbit2D) -> Self

Converts to this type from the input type.
Source§

impl From<Orbit2D> for CompactOrbit2D

Source§

fn from(cached: Orbit2D) -> Self

Converts to this type from the input type.
Source§

impl OrbitTrait2D for CompactOrbit2D

Source§

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

Sets the apoapsis of the orbit.\ Read more
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: Read more
Source§

fn get_transformation_matrix(&self) -> DMat2

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

fn get_pqw_basis_vector_p(&self) -> DVec2

Gets the p basis vector for the perifocal coordinate (PQW) system. Read more
Source§

fn get_pqw_basis_vector_q(&self) -> DVec2

Gets the q basis vector for the perifocal coordinate (PQW) system. Read more
Source§

fn get_eccentricity(&self) -> f64

Gets the eccentricity of the orbit. Read more
Source§

fn set_eccentricity(&mut self, eccentricity: f64)

Sets the eccentricity of the orbit. Read more
Source§

fn get_periapsis(&self) -> f64

Gets the periapsis of the orbit. Read more
Source§

fn set_periapsis(&mut self, periapsis: f64)

Sets the periapsis of the orbit. Read more
Source§

fn get_arg_pe(&self) -> f64

Gets the argument of periapsis of the orbit in radians. Read more
Source§

fn set_arg_pe(&mut self, arg_pe: f64)

Sets the argument of periapsis of the orbit in radians. Read more
Source§

fn get_mean_anomaly_at_epoch(&self) -> f64

Gets the mean anomaly of the orbit at a certain epoch. Read more
Source§

fn set_mean_anomaly_at_epoch(&mut self, mean_anomaly: f64)

Sets the mean anomaly of the orbit at a certain epoch. Read more
Source§

fn get_gravitational_parameter(&self) -> f64

Gets the gravitational parameter of the parent body. Read more
Source§

fn set_gravitational_parameter( &mut self, gravitational_parameter: f64, mode: MuSetterMode2D, )

Sets the gravitational parameter of the parent body. Read more
Source§

fn get_semi_major_axis(&self) -> f64

Gets the semi-major axis of the orbit. Read more
Source§

fn get_semi_minor_axis(&self) -> f64

Gets the semi-minor axis of the orbit. Read more
Source§

fn get_semi_latus_rectum(&self) -> f64

Gets the semi-latus rectum of the orbit, in meters. Read more
Source§

fn get_linear_eccentricity(&self) -> f64

Gets the linear eccentricity of the orbit, in meters. Read more
Source§

fn get_true_anomaly_at_asymptote(&self) -> f64

Gets the true anomaly asymptote (f_∞) of the hyperbolic trajectory. Read more
Source§

fn get_position_at_periapsis(&self) -> DVec2

Gets the position of the orbit at periapsis. Read more
Source§

fn get_position_at_periapsis_unchecked(&self, p_vector: DVec2) -> DVec2

Gets the position of the orbit at periapsis based on a known P basis vector. Read more
Source§

fn get_position_at_apoapsis(&self) -> DVec2

Gets the position of the orbit at apoapsis. Read more
Source§

fn get_position_at_apoapsis_unchecked(&self, p_vector: DVec2) -> DVec2

Gets the position of the orbit at apoapsis based on a known P basis vector. Read more
Source§

fn get_apoapsis(&self) -> f64

Gets the apoapsis of the orbit.
Returns infinity for parabolic orbits.
Returns negative values for hyperbolic orbits.\ Read more
Source§

fn get_mean_motion(&self) -> f64

Gets the mean motion of the orbit, in radians per second. Read more
Source§

fn get_focal_parameter(&self) -> f64

Gets the focal parameter of the orbit, in meters. Read more
Source§

fn get_specific_angular_momentum(&self) -> f64

Gets the specific angular momentum of the orbit, in square meters per second (m^2/s). Read more
Source§

fn get_specific_orbital_energy(&self) -> f64

Gets the specific orbital energy ε of the orbit, in joules per kilogram (J/kg, equiv. to m^2 ⋅ s^-2). Read more
Source§

fn get_area_sweep_rate(&self) -> f64

Gets the area swept out by the orbit in square meters per second (m^2/s). Read more
Source§

fn get_time_of_periapsis(&self) -> f64

Gets the time when the orbit is in periapsis, in seconds since epoch. Read more
Source§

fn get_time_of_apoapsis(&self) -> f64

Gets the time when the orbit is in apoapsis, in seconds since epoch. Read more
Source§

fn get_pqw_basis_vectors(&self) -> (DVec2, DVec2)

Gets the basis vectors for the perifocal coordinate (PQW) system. Read more
Source§

fn get_eccentricity_vector(&self) -> DVec2

Gets the eccentricity vector of this orbit. Read more
Source§

fn get_eccentricity_vector_unchecked(&self, p_vector: DVec2) -> DVec2

Gets the eccentricity vector of this orbit. Read more
Source§

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

Gets the eccentric anomaly at a given mean anomaly in the orbit. Read more
Source§

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

Get an initial guess for the hyperbolic eccentric anomaly of an orbit. Read more
Source§

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

Gets the hyperbolic eccentric anomaly of the orbit. Read more
Source§

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

Gets the elliptic eccentric anomaly of the orbit. Read more
Source§

fn get_eccentric_anomaly_at_true_anomaly(&self, true_anomaly: f64) -> f64

Gets the eccentric anomaly at a given true anomaly in the orbit. Read more
Source§

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

Gets the eccentric anomaly at a given time in the orbit. Read more
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. Read more
Source§

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

Gets the true anomaly at a given mean anomaly in the orbit. Read more
Source§

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

Gets the true anomaly at a given time in the orbit. Read more
Source§

fn get_true_anomaly_at_altitude(&self, altitude: f64) -> f64

Gets the true anomaly where a certain altitude is reached. Read more
Source§

fn get_true_anomaly_at_altitude_unchecked( &self, semi_latus_rectum: f64, altitude: f64, eccentricity_recip: f64, ) -> f64

Gets the true anomaly where a certain altitude is reached. Read more
Source§

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

Gets the mean anomaly at a given time in the orbit. Read more
Source§

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

Gets the mean anomaly at a given eccentric anomaly in the orbit. Read more
Source§

fn get_mean_anomaly_at_elliptic_eccentric_anomaly( &self, eccentric_anomaly: f64, sin_eccentric_anomaly: f64, ) -> f64

Gets the mean anomaly at a given eccentric anomaly in the orbit and its precomputed sine. Read more
Source§

fn get_mean_anomaly_at_hyperbolic_eccentric_anomaly( &self, eccentric_anomaly: f64, sinh_eccentric_anomaly: f64, ) -> f64

Gets the mean anomaly at a given eccentric anomaly in the orbit and its precomputed sine. Read more
Source§

fn get_mean_anomaly_at_true_anomaly(&self, true_anomaly: f64) -> f64

Gets the mean anomaly at a given true anomaly in the orbit. Read more
Source§

fn get_position_at_true_anomaly(&self, angle: f64) -> DVec2

Gets the position at a given angle (true anomaly) in the orbit. Read more
Source§

fn get_position_at_eccentric_anomaly(&self, eccentric_anomaly: f64) -> DVec2

Gets the position at a given eccentric anomaly in the orbit. Read more
Source§

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

Gets the speed at a given angle (true anomaly) in the orbit. Read more
Source§

fn get_speed_at_altitude(&self, altitude: f64) -> f64

Gets the speed at a given altitude in the orbit. Read more
Source§

fn get_speed_at_periapsis(&self) -> f64

Gets the speed at the periapsis of the orbit. Read more
Source§

fn get_speed_at_apoapsis(&self) -> f64

Gets the speed at the apoapsis of the orbit. Read more
Source§

fn get_speed_at_infinity(&self) -> f64

Gets the hyperbolic excess speed (v_∞) of the trajectory. Read more
Source§

fn get_speed_at_time(&self, time: f64) -> f64

Gets the speed at a given time in the orbit. Read more
Source§

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

Gets the speed at a given eccentric anomaly in the orbit. Read more
Source§

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

Gets the speed at a given mean anomaly in the orbit. Read more
Source§

fn get_pqw_velocity_at_true_anomaly(&self, angle: f64) -> DVec2

Gets the velocity at a given angle (true anomaly) in the orbit in the perifocal coordinate system. Read more
Source§

fn get_pqw_velocity_at_periapsis(&self) -> DVec2

Gets the velocity at the periapsis of the orbit in the perifocal coordinate system. Read more
Source§

fn get_pqw_velocity_at_apoapsis(&self) -> DVec2

Gets the velocity at the apoapsis of the orbit in the perifocal coordinate system. Read more
Source§

fn get_pqw_velocity_at_incoming_asymptote(&self) -> DVec2

Gets the velocity at the incoming asymptote of the trajectory in the perifocal coordinate system. Read more
Source§

fn get_pqw_velocity_at_outgoing_asymptote(&self) -> DVec2

Gets the velocity at the outgoing asymptote of the trajectory in the perifocal coordinate system. Read more
Source§

fn get_pqw_velocity_at_eccentric_anomaly(&self, eccentric_anomaly: f64) -> DVec2

Gets the velocity at a given eccentric anomaly in the orbit in the perifocal coordinate system. Read more
Source§

fn get_pqw_velocity_at_mean_anomaly(&self, mean_anomaly: f64) -> DVec2

Gets the velocity at a given mean anomaly in the orbit in the perifocal coordinate system. Read more
Source§

fn get_pqw_velocity_at_eccentric_anomaly_unchecked( &self, outer_mult: f64, q_mult: f64, trig_ecc_anom: (f64, f64), ) -> DVec2

Gets the velocity at a given eccentric anomaly in the orbit in the perifocal coordinate system. Read more
Source§

fn get_pqw_velocity_at_time(&self, time: f64) -> DVec2

Gets the velocity at a given time in the orbit in the perifocal coordinate system. Read more
Source§

fn get_pqw_position_at_true_anomaly(&self, angle: f64) -> DVec2

Gets the position at a given angle (true anomaly) in the orbit in the perifocal coordinate system. Read more
Source§

fn get_pqw_position_at_true_anomaly_unchecked( &self, altitude: f64, sincos_angle: (f64, f64), ) -> DVec2

Gets the position at a certain point in the orbit in the perifocal coordinate system. Read more
Source§

fn get_pqw_position_at_time(&self, time: f64) -> DVec2

Gets the position at a given time in the orbit in the perifocal coordinate system. Read more
Source§

fn get_pqw_position_at_eccentric_anomaly(&self, eccentric_anomaly: f64) -> DVec2

Gets the position at a given eccentric anomaly in the orbit in the perifocal coordinate system. Read more
Source§

fn get_velocity_at_true_anomaly(&self, angle: f64) -> DVec2

Gets the velocity at a given angle (true anomaly) in the orbit. Read more
Source§

fn get_velocity_at_periapsis(&self) -> DVec2

Gets the velocity at the periapsis of the orbit. Read more
Source§

fn get_velocity_at_apoapsis(&self) -> DVec2

Gets the velocity at the apoapsis of the orbit. Read more
Source§

fn get_velocity_at_incoming_asymptote(&self) -> DVec2

Gets the velocity at the incoming asymptote of the trajectory. Read more
Source§

fn get_velocity_at_outgoing_asymptote(&self) -> DVec2

Gets the velocity at the outgoing asymptote of the trajectory. Read more
Source§

fn get_velocity_at_eccentric_anomaly(&self, eccentric_anomaly: f64) -> DVec2

Gets the velocity at a given eccentric anomaly in the orbit. Read more
Source§

fn get_velocity_at_time(&self, time: f64) -> DVec2

Gets the velocity at a given time in the orbit. Read more
Source§

fn get_velocity_at_mean_anomaly(&self, mean_anomaly: f64) -> DVec2

Gets the velocity at a given mean anomaly in the orbit. Read more
Source§

fn get_altitude_at_true_anomaly(&self, true_anomaly: f64) -> f64

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

fn get_altitude_at_true_anomaly_unchecked( &self, semi_latus_rectum: f64, cos_true_anomaly: f64, ) -> f64

Gets the altitude of the body from its parent given the cosine of the true anomaly. Read more
Source§

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

Gets the altitude at a given eccentric anomaly in the orbit. Read more
Source§

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

Gets the altitude of the body from its parent at a given time in the orbit. Read more
Source§

fn get_position_at_time(&self, time: f64) -> DVec2

Gets the position at a given time in the orbit. Read more
Source§

fn get_position_at_mean_anomaly(&self, mean_anomaly: f64) -> DVec2

Gets the position at a given mean anomaly in the orbit. Read more
Source§

fn get_state_vectors_at_eccentric_anomaly( &self, eccentric_anomaly: f64, ) -> StateVectors2D

Gets the position and velocity at a given eccentric anomaly in the orbit. Read more
Source§

fn get_state_vectors_at_true_anomaly(&self, true_anomaly: f64) -> StateVectors2D

Gets the position and velocity at a given angle (true anomaly) in the orbit. Read more
Source§

fn get_state_vectors_at_mean_anomaly(&self, mean_anomaly: f64) -> StateVectors2D

Gets the position and velocity at a given mean anomaly in the orbit. Read more
Source§

fn get_state_vectors_at_time(&self, time: f64) -> StateVectors2D

Gets the position and velocity at a given time in the orbit. Read more
Source§

fn get_state_vectors_from_unchecked_parts( &self, sqrt_abs_gm_a: f64, altitude: f64, q_mult: f64, trig_ecc_anom: (f64, f64), sincos_angle: (f64, f64), ) -> StateVectors2D

Gets the position and velocity at a certain point in the orbit. Read more
Source§

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

Gets the time of the orbit at a certain mean anomaly. Read more
Source§

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

Gets the time of the orbit at a certain eccentric anomaly. Read more
Source§

fn get_time_at_true_anomaly(&self, true_anomaly: f64) -> f64

Gets the time of the orbit at a certain true anomaly. Read more
Source§

fn transform_pqw_vector(&self, position: DVec2) -> DVec2

Transforms a position from the perifocal coordinate (PQW) system into world space, using the orbital parameters. Read more
Source§

fn is_circular(&self) -> bool

Gets whether or not the orbit is circular. Read more
Source§

fn is_elliptic(&self) -> bool

Gets whether or not the orbit is strictly elliptic. Read more
Source§

fn is_closed(&self) -> bool

Gets whether or not the orbit is closed. Read more
Source§

fn is_parabolic(&self) -> bool

Gets whether or not the trajectory is parabolic. Read more
Source§

fn is_hyperbolic(&self) -> bool

Gets whether or not the trajectory is hyperbolic. Read more
Source§

fn is_open(&self) -> bool

Gets whether or not the trajectory is open. Read more
Source§

fn get_orbital_period(&self) -> f64

Gets the time it takes to complete one revolution of the orbit. Read more
Source§

impl PartialEq for CompactOrbit2D

Source§

fn eq(&self, other: &CompactOrbit2D) -> bool

Tests for self and other values to be equal, and is used by ==.
1.0.0 · Source§

fn ne(&self, other: &Rhs) -> bool

Tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
Source§

impl Copy for CompactOrbit2D

Source§

impl StructuralPartialEq for CompactOrbit2D

Auto Trait Implementations§

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> CloneToUninit for T
where T: Clone,

Source§

unsafe fn clone_to_uninit(&self, dest: *mut u8)

🔬This is a nightly-only experimental API. (clone_to_uninit)
Performs copy-assignment from self to dest. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> ToOwned for T
where T: Clone,

Source§

type Owned = T

The resulting type after obtaining ownership.
Source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
Source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.