Struct nyx_space::cosmic::eclipse::Orbit

source ·
pub struct Orbit {
    pub x: f64,
    pub y: f64,
    pub z: f64,
    pub vx: f64,
    pub vy: f64,
    pub vz: f64,
    pub dt: Epoch,
    pub frame: Frame,
    pub stm: Option<Matrix6<f64>>,
}
Expand description

Orbit defines an orbital state

Unless noted otherwise, algorithms are from GMAT 2016a StateConversionUtil.cpp. Regardless of the constructor used, this struct stores all the state information in Cartesian coordinates as these are always non singular. Note: although not yet supported, this struct may change once True of Date or other nutation frames are added to the toolkit.

Fields§

§x: f64

in km

§y: f64

in km

§z: f64

in km

§vx: f64

in km/s

§vy: f64

in km/s

§vz: f64

in km/s

§dt: Epoch§frame: Frame

Frame contains everything we need to compute state information

§stm: Option<Matrix6<f64>>

Optionally stores the state transition matrix from the start of the propagation until the current time (i.e. trajectory STM, not step-size STM)

Implementations§

source§

impl Orbit

source

pub fn cartesian( x: f64, y: f64, z: f64, vx: f64, vy: f64, vz: f64, dt: Epoch, frame: Frame ) -> Self

Creates a new Orbit in the provided frame at the provided Epoch.

Units: km, km, km, km/s, km/s, km/s

source

pub fn cartesian_stm( x: f64, y: f64, z: f64, vx: f64, vy: f64, vz: f64, dt: Epoch, frame: Frame ) -> Self

Creates a new Orbit and initializes its STM.

source

pub fn from_position(x: f64, y: f64, z: f64, dt: Epoch, frame: Frame) -> Self

Creates a new Orbit in the provided frame at the provided Epoch in time with 0.0 velocity.

Units: km, km, km

source

pub fn cartesian_vec(state: &Vector6<f64>, dt: Epoch, frame: Frame) -> Self

Creates a new Orbit around in the provided frame from the borrowed state vector

The state vector must be x, y, z, vx, vy, vz. This function is a shortcut to cartesian and as such it has the same unit requirements.

source

pub fn rmag(&self) -> f64

Returns the magnitude of the radius vector in km

source

pub fn vmag(&self) -> f64

Returns the magnitude of the velocity vector in km/s

source

pub fn radius(&self) -> Vector3<f64>

Returns the radius vector of this Orbit in [km, km, km]

source

pub fn velocity(&self) -> Vector3<f64>

Returns the velocity vector of this Orbit in [km/s, km/s, km/s]

source

pub fn to_cartesian_vec(self) -> Vector6<f64>

Returns this state as a Cartesian Vector6 in [km, km, km, km/s, km/s, km/s]

Note that the time is not returned in the vector.

source

pub fn distance_to(&self, other: &Orbit) -> f64

Returns the distance in kilometers between this state and another state. Will panic is the frames are different

source

pub fn distance_to_point(&self, other: &Vector3<f64>) -> f64

Returns the distance in kilometers between this state and a point assumed to be in the same frame.

source

pub fn r_hat(&self) -> Vector3<f64>

Returns the unit vector in the direction of the state radius

source

pub fn v_hat(&self) -> Vector3<f64>

Returns the unit vector in the direction of the state velocity

source

pub fn keplerian( sma: f64, ecc: f64, inc: f64, raan: f64, aop: f64, ta: f64, dt: Epoch, frame: Frame ) -> Self

Creates a new Orbit around the provided Celestial or Geoid frame from the Keplerian orbital elements.

Units: km, none, degrees, degrees, degrees, degrees

WARNING: This function will panic if the singularities in the conversion are expected. NOTE: The state is defined in Cartesian coordinates as they are non-singular. This causes rounding errors when creating a state from its Keplerian orbital elements (cf. the state tests). One should expect these errors to be on the order of 1e-12.

source

pub fn keplerian_altitude( sma_altitude: f64, ecc: f64, inc: f64, raan: f64, aop: f64, ta: f64, dt: Epoch, frame: Frame ) -> Self

Creates a new Orbit from the provided semi-major axis altitude in kilometers

source

pub fn keplerian_apsis_radii( r_a: f64, r_p: f64, inc: f64, raan: f64, aop: f64, ta: f64, dt: Epoch, frame: Frame ) -> Self

Creates a new Orbit from the provided radii of apoapsis and periapsis, in kilometers

source

pub fn keplerian_apsis_altitude( a_a: f64, a_p: f64, inc: f64, raan: f64, aop: f64, ta: f64, dt: Epoch, frame: Frame ) -> Self

Creates a new Orbit from the provided altitudes of apoapsis and periapsis, in kilometers

source

pub fn keplerian_vec(state: &Vector6<f64>, dt: Epoch, frame: Frame) -> Self

Creates a new Orbit around the provided frame from the borrowed state vector

The state vector must be sma, ecc, inc, raan, aop, ta. This function is a shortcut to cartesian and as such it has the same unit requirements.

source

pub fn from_geodesic( latitude: f64, longitude: f64, height: f64, dt: Epoch, frame: Frame ) -> Self

Creates a new Orbit from the geodetic latitude (φ), longitude (λ) and height with respect to the ellipsoid of the frame.

Units: degrees, degrees, km NOTE: This computation differs from the spherical coordinates because we consider the flattening of body. Reference: G. Xu and Y. Xu, “GPS”, DOI 10.1007/978-3-662-50367-6_2, 2016 WARNING: This uses the rotational rates known to Nyx. For other objects, use from_altlatlong for other celestial bodies.

source

pub fn from_altlatlong( latitude: f64, longitude: f64, height: f64, angular_velocity: f64, dt: Epoch, frame: Frame ) -> Self

Creates a new Orbit from the latitude (φ), longitude (λ) and height with respect to the frame’s ellipsoid.

Units: degrees, degrees, km, rad/s NOTE: This computation differs from the spherical coordinates because we consider the flattening of body. Reference: G. Xu and Y. Xu, “GPS”, DOI 10.1007/978-3-662-50367-6_2, 2016

source

pub fn to_keplerian_vec(self) -> Vector6<f64>

Returns this state as a Keplerian Vector6 in [km, none, degrees, degrees, degrees, degrees]

Note that the time is not returned in the vector.

source

pub fn hvec(&self) -> Vector3<f64>

Returns the orbital momentum vector

source

pub fn hx(&self) -> f64

Returns the orbital momentum value on the X axis

source

pub fn hy(&self) -> f64

Returns the orbital momentum value on the Y axis

source

pub fn hz(&self) -> f64

Returns the orbital momentum value on the Z axis

source

pub fn hmag(&self) -> f64

Returns the norm of the orbital momentum

source

pub fn energy(&self) -> f64

Returns the specific mechanical energy in km^2/s^2

source

pub fn with_radius(self, new_radius: &Vector3<f64>) -> Self

Returns a copy of the state with a new radius

source

pub fn with_velocity(self, new_velocity: &Vector3<f64>) -> Self

Returns a copy of the state with a new radius

source

pub fn sma(&self) -> f64

Returns the semi-major axis in km

source

pub fn set_sma(&mut self, new_sma_km: f64)

Mutates this orbit to change the SMA

source

pub fn with_sma(self, new_sma_km: f64) -> Self

Returns a copy of the state with a new SMA

source

pub fn add_sma(self, delta_sma: f64) -> Self

Returns a copy of the state with a provided SMA added to the current one

source

pub fn sma_altitude(&self) -> f64

Returns the SMA altitude in km

source

pub fn period(&self) -> Duration

Returns the period in seconds

source

pub fn evec(&self) -> Vector3<f64>

Returns the eccentricity vector (no unit)

source

pub fn ecc(&self) -> f64

Returns the eccentricity (no unit)

source

pub fn set_ecc(&mut self, new_ecc: f64)

Mutates this orbit to change the ECC

source

pub fn with_ecc(self, new_ecc: f64) -> Self

Returns a copy of the state with a new ECC

source

pub fn add_ecc(self, delta_ecc: f64) -> Self

Returns a copy of the state with a provided ECC added to the current one

source

pub fn inc(&self) -> f64

Returns the inclination in degrees

source

pub fn set_inc(&mut self, new_inc: f64)

Mutates this orbit to change the INC

source

pub fn with_inc(self, new_inc: f64) -> Self

Returns a copy of the state with a new INC

source

pub fn add_inc(self, delta_inc: f64) -> Self

Returns a copy of the state with a provided INC added to the current one

source

pub fn aop(&self) -> f64

Returns the argument of periapsis in degrees

source

pub fn set_aop(&mut self, new_aop: f64)

Mutates this orbit to change the AOP

source

pub fn with_aop(self, new_aop: f64) -> Self

Returns a copy of the state with a new AOP

source

pub fn add_aop(self, delta_aop: f64) -> Self

Returns a copy of the state with a provided AOP added to the current one

source

pub fn raan(&self) -> f64

Returns the right ascension of ther ascending node in degrees

source

pub fn set_raan(&mut self, new_raan: f64)

Mutates this orbit to change the RAAN

source

pub fn with_raan(self, new_raan: f64) -> Self

Returns a copy of the state with a new RAAN

source

pub fn add_raan(self, delta_raan: f64) -> Self

Returns a copy of the state with a provided RAAN added to the current one

source

pub fn ta(&self) -> f64

Returns the true anomaly in degrees between 0 and 360.0

NOTE: This function will emit a warning stating that the TA should be avoided if in a very near circular orbit Code from https://github.com/ChristopherRabotin/GMAT/blob/80bde040e12946a61dae90d9fc3538f16df34190/src/gmatutil/util/StateConversionUtil.cpp#L6835

LIMITATION: For an orbit whose true anomaly is (very nearly) 0.0 or 180.0, this function may return either 0.0 or 180.0 with a very small time increment. This is due to the precision of the cosine calculation: if the arccosine calculation is out of bounds, the sign of the cosine of the true anomaly is used to determine whether the true anomaly should be 0.0 or 180.0. In other words, there is an ambiguity in the computation in the true anomaly exactly at 180.0 and 0.0.

source

pub fn set_ta(&mut self, new_ta: f64)

Mutates this orbit to change the TA

source

pub fn with_ta(self, new_ta: f64) -> Self

Returns a copy of the state with a new TA

source

pub fn add_ta(self, delta_ta: f64) -> Self

Returns a copy of the state with a provided TA added to the current one

source

pub fn with_apoapsis_periapsis(self, new_ra: f64, new_rp: f64) -> Self

Returns a copy of this state with the provided apoasis and periapse

source

pub fn add_apoapsis_periapsis(self, delta_ra: f64, delta_rp: f64) -> Self

Returns a copy of this state with the provided apoasis and periapse added to the current values

source

pub fn tlong(&self) -> f64

Returns the true longitude in degrees

source

pub fn aol(&self) -> f64

Returns the argument of latitude in degrees

NOTE: If the orbit is near circular, the AoL will be computed from the true longitude instead of relying on the ill-defined true anomaly.

source

pub fn periapsis(&self) -> f64

Returns the radius of periapsis (or perigee around Earth), in kilometers.

source

pub fn apoapsis(&self) -> f64

Returns the radius of apoapsis (or apogee around Earth), in kilometers.

source

pub fn periapsis_altitude(&self) -> f64

Returns the altitude of periapsis (or perigee around Earth), in kilometers.

source

pub fn apoapsis_altitude(&self) -> f64

Returns the altitude of apoapsis (or apogee around Earth), in kilometers.

source

pub fn ea(&self) -> f64

Returns the eccentric anomaly in degrees

This is a conversion from GMAT’s StateConversionUtil::TrueToEccentricAnomaly

source

pub fn fpa(&self) -> f64

Returns the flight path angle in degrees

source

pub fn ma(&self) -> f64

Returns the mean anomaly in degrees

This is a conversion from GMAT’s StateConversionUtil::TrueToMeanAnomaly

source

pub fn semi_parameter(&self) -> f64

Returns the semi parameter (or semilatus rectum)

source

pub fn is_brouwer_short_valid(&self) -> bool

Returns whether this state satisfies the requirement to compute the Mean Brouwer Short orbital element set.

This is a conversion from GMAT’s StateConversionUtil::CartesianToBrouwerMeanShort. The details are at the log level info. NOTE: Mean Brouwer Short are only defined around Earth. However, nyx does not check the main celestial body around which the state is defined (GMAT does perform this verification).

source

pub fn geodetic_longitude(&self) -> f64

Returns the geodetic longitude (λ) in degrees. Value is between 0 and 360 degrees.

Although the reference is not Vallado, the math from Vallado proves to be equivalent. Reference: G. Xu and Y. Xu, “GPS”, DOI 10.1007/978-3-662-50367-6_2, 2016

source

pub fn geodetic_latitude(&self) -> f64

Returns the geodetic latitude (φ) in degrees. Value is between -180 and +180 degrees.

Reference: Vallado, 4th Ed., Algorithm 12 page 172.

source

pub fn geodetic_height(&self) -> f64

Returns the geodetic height in km.

Reference: Vallado, 4th Ed., Algorithm 12 page 172.

source

pub fn right_ascension(&self) -> f64

Returns the right ascension of this orbit in degrees

source

pub fn declination(&self) -> f64

Returns the declination of this orbit in degrees

source

pub fn semi_minor_axis(&self) -> f64

Returns the semi minor axis in km, includes code for a hyperbolic orbit

source

pub fn velocity_declination(&self) -> f64

Returns the velocity declination of this orbit in degrees

source

pub fn b_plane(&self) -> Result<BPlane, NyxError>

source

pub fn c3(&self) -> f64

Returns the $C_3$ of this orbit

source

pub fn vinf_periapsis(&self, turn_angle_degrees: f64) -> Result<f64, NyxError>

Returns the radius of periapse in kilometers for the provided turn angle of this hyperbolic orbit.

source

pub fn vinf_turn_angle(&self, periapsis_km: f64) -> Result<f64, NyxError>

Returns the turn angle in degrees for the provided radius of periapse passage of this hyperbolic orbit

source

pub fn hyperbolic_anomaly(&self) -> Result<f64, NyxError>

Returns the hyperbolic anomaly in degrees between 0 and 360.0

source

pub fn dcm_from_traj_frame(&self, from: Frame) -> Result<Matrix3<f64>, NyxError>

Returns the direct cosine rotation matrix to convert to this state’s frame (inertial or otherwise).

Example

let dcm_vnc2inertial = orbit.dcm_from_traj_frame(Frame::VNC)?; let vector_inertial = dcm_vnc2inertial * vector_vnc;

source

pub fn dcm6x6_from_traj_frame( &self, from: Frame ) -> Result<Matrix6<f64>, NyxError>

Returns a 6x6 DCM to convert to this inertial state. WARNING: This DCM does NOT contain the corrections needed for the transport theorem, and therefore the velocity rotation is wrong.

source

pub fn apply_dv(&mut self, dv: Vector3<f64>)

Apply the provided delta-v (in km/s)

source

pub fn with_dv(self, dv: Vector3<f64>) -> Self

Copies this orbit after applying the provided delta-v (in km/s)

source

pub fn rotate_by(&mut self, dcm: Matrix6<f64>)

Rotate this state provided a direct cosine matrix of position and velocity

source

pub fn with_rotation_by(&self, dcm: Matrix6<f64>) -> Self

Rotate this state provided a direct cosine matrix of position and velocity

source

pub fn position_rotated_by(&mut self, dcm: Matrix3<f64>)

Rotate the position and the velocity of this state provided a direct cosine matrix of position and velocity WARNING: You only want to use this if you’ll only be using the position components of the rotated state. This does not account for the transport theorem and therefore is physically WRONG.

source

pub fn with_position_rotated_by(&self, dcm: Matrix3<f64>) -> Self

Rotate the position of this state provided a direct cosine matrix of position and velocity WARNING: You only want to use this if you’ll only be using the position components of the rotated state. This does not account for the transport theorem and therefore is physically WRONG.

source

pub fn enable_stm(&mut self)

Sets the STM of this state of identity, which also enables computation of the STM for spacecraft navigation

source

pub fn disable_stm(&mut self)

Disable the STM of this state

source

pub fn with_stm(self) -> Self

Copies the current state but sets the STM to identity

source

pub fn without_stm(self) -> Self

Copies the current state but disables the STM

source

pub fn rss(&self, other: &Self) -> (f64, f64)

Returns the root sum square error between this state and the other, in kilometers for the position and kilometers per second in velocity

source

pub fn eq_within(&self, other: &Self, radial_tol: f64, velocity_tol: f64) -> bool

Returns whether this orbit and another are equal within the specified radial and velocity absolute tolerances

source

pub fn to_objectives( &self, params: &[StateParameter] ) -> Result<Vec<Objective>, NyxError>

Use the current orbit as a template to generate mission design objectives. Note: this sets the objective tolerances to be quite tight, so consider modifying them.

source

pub fn disperse( &self, mean: Vector6<f64>, cov: Matrix6<f64> ) -> Result<MultivariateNormal<Orbit>, NyxError>

Create a multivariate normal dispersion structure from this orbit with the provided mean and covariance, specified as {X, Y, Z, VX, VY, VZ} in km and km/s

source

pub fn disperse_zero_mean( &self, cov: Matrix6<f64> ) -> Result<MultivariateNormal<Orbit>, NyxError>

Create a multivariate normal dispersion structure from this orbit with the provided covariance, specified as {X, Y, Z, VX, VY, VZ} in km and km/s

Trait Implementations§

source§

impl Add<&Orbit> for &Orbit

source§

fn add(self, other: &Orbit) -> Orbit

Add one state from another. Frame must be manually changed if needed. STM will be copied from &self.

§

type Output = Orbit

The resulting type after applying the + operator.
source§

impl Add<Matrix<f64, Const<6>, Const<1>, <DefaultAllocator as Allocator<f64, Const<6>, Const<1>>>::Buffer>> for Orbit

source§

fn add(self, other: OVector<f64, Const<6>>) -> Self

Adds the provided state deviation to this orbit

§

type Output = Orbit

The resulting type after applying the + operator.
source§

impl Add<Orbit> for Orbit

source§

fn add(self, other: Orbit) -> Orbit

Add one state from another. Frame must be manually changed if needed. STM will be copied from &self.

§

type Output = Orbit

The resulting type after applying the + operator.
source§

impl AddAssign<Orbit> for Orbit

source§

fn add_assign(&mut self, other: Self)

Performs the += operation. Read more
source§

impl Clone for Orbit

source§

fn clone(&self) -> Orbit

Returns a copy 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 Orbit

source§

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

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

impl Default for Orbit

source§

fn default() -> Self

Returns the “default value” for a type. Read more
source§

impl Display for Orbit

source§

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

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

impl EstimateFrom<BaseSpacecraft<GuidanceMode>> for Orbit

source§

fn extract(from: Spacecraft) -> Self

From the state extract the state to be estimated
source§

impl EstimateFrom<Orbit> for Orbit

source§

fn extract(from: Orbit) -> Self

From the state extract the state to be estimated
source§

impl EventEvaluator<Orbit> for Event

source§

fn epoch_precision(&self) -> Duration

source§

fn value_precision(&self) -> f64

source§

fn eval(&self, state: &Orbit) -> f64

source§

fn eval_crossing(&self, prev_state: &S, next_state: &S) -> bool

source§

impl EventEvaluator<Orbit> for PenumbraEvent

source§

fn epoch_precision(&self) -> Duration

Stop searching when the time has converged to less than 0.1 seconds

source§

fn value_precision(&self) -> f64

Finds the slightest penumbra within 2%(i.e. 98% in visibility)

source§

fn eval(&self, observer: &Orbit) -> f64

source§

fn eval_crossing(&self, prev_state: &S, next_state: &S) -> bool

source§

impl EventEvaluator<Orbit> for UmbraEvent

source§

fn epoch_precision(&self) -> Duration

Stop searching when the time has converged to less than 0.1 seconds

source§

fn value_precision(&self) -> f64

Finds the darkest part of an eclipse within 2% of penumbra (i.e. 98% in shadow)

source§

fn eval(&self, observer: &Orbit) -> f64

source§

fn eval_crossing(&self, prev_state: &S, next_state: &S) -> bool

source§

impl From<Orbit> for OrbitDual

source§

fn from(orbit: Orbit) -> Self

Initialize a new OrbitDual from an orbit, no other initializers

source§

impl InterpState for Orbit

source§

fn params() -> Vec<StateParameter>

Return the parameters in order
source§

fn value_and_deriv(&self, param: &StateParameter) -> Result<(f64, f64), NyxError>

Return the requested parameter and its time derivative
source§

fn set_value_and_deriv( &mut self, param: &StateParameter, value: f64, _: f64 ) -> Result<(), NyxError>

Sets the requested parameter
source§

fn deriv(&self, param: &StateParameter) -> Result<f64, NyxError>

Return the time derivative requested parameter
source§

impl LowerExp for Orbit

source§

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

Formats the value using the given formatter.
source§

impl LowerHex for Orbit

source§

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

Formats the value using the given formatter.
source§

impl MdHdlr<Orbit> for OrbitStateOutput

source§

fn handle(&mut self, state: &Orbit)

source§

impl MeasurementDevice<Orbit, StdMeasurement> for GroundStation

source§

fn measure(&self, rx: &Orbit) -> Option<StdMeasurement>

Perform a measurement from the ground station to the receiver (rx).

source§

impl NavSolution<Orbit> for KfEstimate<Orbit>

source§

fn orbital_state(&self) -> Orbit

source§

fn expected_state(&self) -> Orbit

Returns the nominal state as computed by the dynamics
source§

impl Neg for &Orbit

source§

fn neg(self) -> Self::Output

Subtract one state from another. STM will be copied form &self.

§

type Output = Orbit

The resulting type after applying the - operator.
source§

impl Neg for Orbit

source§

fn neg(self) -> Self::Output

Subtract one state from another. STM will be copied from &self.

§

type Output = Orbit

The resulting type after applying the - operator.
source§

impl PartialEq<Orbit> for Orbit

source§

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

Two states are equal if their position are equal within one centimeter and their velocities within one centimeter per second.

1.0.0 · source§

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

This method tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
source§

impl Serialize for Orbit

source§

fn serialize<S>(&self, serializer: S) -> Result<S::Ok, S::Error>where S: Serializer,

NOTE: This is not part of unit testing because there is no deseralization of Orbit (yet)

source§

impl State for Orbit

Implementation of Orbit as a State for orbital dynamics with STM

source§

fn zeros() -> Self

Returns a state whose position, velocity and frame are zero, and STM is I_{6x6}.

§

type Size = Const<6>

Size of the state and its STM
§

type VecLength = Const<42>

source§

fn reset_stm(&mut self)

Return this state as a vector for the propagation/estimation By default, this is not implemented. This function must be implemented when filtering on this state.
source§

fn as_vector(&self) -> Result<OVector<f64, Const<42>>, NyxError>

Return this state as a vector for the propagation/estimation
source§

fn set( &mut self, epoch: Epoch, vector: &OVector<f64, Const<42>> ) -> Result<(), NyxError>

Set this state
source§

fn stm(&self) -> Result<Matrix6<f64>, NyxError>

Return this state as a vector for the propagation/estimation By default, this is not implemented. This function must be implemented when filtering on this state.
source§

fn epoch(&self) -> Epoch

Retrieve the Epoch
source§

fn set_epoch(&mut self, epoch: Epoch)

Set the Epoch
source§

fn add(self, other: OVector<f64, Self::Size>) -> Self

By default, this is not implemented. This function must be implemented when filtering on this state.
source§

fn value(&self, param: &StateParameter) -> Result<f64, NyxError>

Return the value of the parameter, returns an error by default
source§

fn set_value(&mut self, param: &StateParameter, val: f64) -> Result<(), NyxError>

Allows setting the value of the given parameter. NOTE: Most parameters where the value is available CANNOT be also set for that parameter (it’s a much harder problem!)
source§

fn set_with_delta_seconds( self, delta_t_s: f64, vector: &OVector<f64, Self::VecLength> ) -> Selfwhere DefaultAllocator: Allocator<f64, Self::VecLength>,

Reconstruct a new State from the provided delta time in seconds compared to the current state and with the provided vector.
source§

impl Sub<&Orbit> for &Orbit

source§

fn sub(self, other: &Orbit) -> Orbit

Subtract one state from another. STM will be copied from &self.

§

type Output = Orbit

The resulting type after applying the - operator.
source§

impl Sub<Orbit> for Orbit

source§

fn sub(self, other: Orbit) -> Orbit

Subtract one state from another. STM will be copied from &self.

§

type Output = Orbit

The resulting type after applying the - operator.
source§

impl SubAssign<Orbit> for Orbit

source§

fn sub_assign(&mut self, other: Self)

Performs the -= operation. Read more
source§

impl UpperExp for Orbit

source§

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

Formats the value using the given formatter.
source§

impl UpperHex for Orbit

source§

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

Formats the value using the given formatter.
source§

impl Copy for Orbit

Auto Trait Implementations§

§

impl RefUnwindSafe for Orbit

§

impl Send for Orbit

§

impl Sync for Orbit

§

impl Unpin for Orbit

§

impl UnwindSafe for Orbit

Blanket Implementations§

source§

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

source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
source§

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

const: unstable · source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
source§

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

const: unstable · source§

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

Mutably borrows from an owned value. Read more
source§

impl<T> From<T> for T

const: unstable · source§

fn from(t: T) -> T

Returns the argument unchanged.

source§

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

const: unstable · 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.

§

impl<T> Pointable for T

§

const ALIGN: usize = mem::align_of::<T>()

The alignment of pointer.
§

type Init = T

The type for initializers.
§

unsafe fn init(init: <T as Pointable>::Init) -> usize

Initializes a with the given initializer. Read more
§

unsafe fn deref<'a>(ptr: usize) -> &'a T

Dereferences the given pointer. Read more
§

unsafe fn deref_mut<'a>(ptr: usize) -> &'a mut T

Mutably dereferences the given pointer. Read more
§

unsafe fn drop(ptr: usize)

Drops the object pointed to by the given pointer. Read more
§

impl<T> Printing<T> for Twhere T: Display,

§

fn to_str(self) -> String

Method to serialize. Decorates Vecs with square brackets and tuples with round ones. Implementation code is in printing.rs.
§

fn to_plainstr(self) -> String

Method to serialize in minimal form (space separated, no brackets) Implementation code is in printing.rs.
§

fn rd(self) -> String

Printable in red
§

fn gr(self) -> String

Printable in green
§

fn bl(self) -> String

Printable in blue
§

fn yl(self) -> String

Printable in yellow
§

fn mg(self) -> String

Printable in magenta
§

fn cy(self) -> String

Printable in cyan
§

fn wvec(self, f: &mut File) -> Result<(), Error>

Method to write vector(s) to file f (space separated, without brackets). Passes up io errors
§

fn pvec(self)

Method to print vector(s) to stdout (space separated,without brackets).
source§

impl<T> Same<T> for T

§

type Output = T

Should always be Self
§

impl<SS, SP> SupersetOf<SS> for SPwhere SS: SubsetOf<SP>,

§

fn to_subset(&self) -> Option<SS>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
§

fn is_in_subset(&self) -> bool

Checks if self is actually part of its subset T (and can be converted to it).
§

fn to_subset_unchecked(&self) -> SS

Use with care! Same as self.to_subset but without any property checks. Always succeeds.
§

fn from_subset(element: &SS) -> SP

The inclusion map: converts self to the equivalent element of its superset.
source§

impl<T> ToOwned for Twhere T: Clone,

§

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> ToString for Twhere T: Display + ?Sized,

source§

default fn to_string(&self) -> String

Converts the given value to a String. Read more
source§

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

§

type Error = Infallible

The type returned in the event of a conversion error.
const: unstable · source§

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

Performs the conversion.
source§

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

§

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

The type returned in the event of a conversion error.
const: unstable · source§

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

Performs the conversion.
§

impl<V, T> VZip<V> for Twhere V: MultiLane<T>,

§

fn vzip(self) -> V

§

impl<T, Right> ClosedAdd<Right> for Twhere T: Add<Right, Output = T> + AddAssign<Right>,

§

impl<T> ClosedNeg for Twhere T: Neg<Output = T>,

§

impl<T, Right> ClosedSub<Right> for Twhere T: Sub<Right, Output = T> + SubAssign<Right>,

source§

impl<T> Scalar for Twhere T: 'static + Clone + PartialEq<T> + Debug,