Struct nyx_space::md::ui::OrbitDual

source ·
pub struct OrbitDual {
    pub x: OHyperdual<f64, U7>,
    pub y: OHyperdual<f64, U7>,
    pub z: OHyperdual<f64, U7>,
    pub vx: OHyperdual<f64, U7>,
    pub vy: OHyperdual<f64, U7>,
    pub vz: OHyperdual<f64, U7>,
    pub dt: Epoch,
    pub frame: Frame,
}
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: OHyperdual<f64, U7>

in km

§y: OHyperdual<f64, U7>

in km

§z: OHyperdual<f64, U7>

in km

§vx: OHyperdual<f64, U7>

in km/s

§vy: OHyperdual<f64, U7>

in km/s

§vz: OHyperdual<f64, U7>

in km/s

§dt: Epoch§frame: Frame

Frame contains everything we need to compute state information

Implementations§

source§

impl OrbitDual

source

pub fn partial_for( &self, param: &StateParameter ) -> Result<OrbitPartial, NyxError>

source

pub fn rmag(&self) -> OrbitPartial

Returns the magnitude of the radius vector in km

source

pub fn vmag(&self) -> OrbitPartial

Returns the magnitude of the velocity vector in km/s

source

pub fn hx(&self) -> OrbitPartial

Returns the orbital momentum value on the X axis

source

pub fn hy(&self) -> OrbitPartial

Returns the orbital momentum value on the Y axis

source

pub fn hz(&self) -> OrbitPartial

Returns the orbital momentum value on the Z axis

source

pub fn hmag(&self) -> OrbitPartial

Returns the norm of the orbital momentum

source

pub fn energy(&self) -> OrbitPartial

Returns the specific mechanical energy

source

pub fn sma(&self) -> OrbitPartial

Returns the semi-major axis in km

source

pub fn ecc(&self) -> OrbitPartial

Returns the eccentricity (no unit)

source

pub fn inc(&self) -> OrbitPartial

Returns the inclination in degrees

source

pub fn aop(&self) -> OrbitPartial

Returns the argument of periapsis in degrees

source

pub fn raan(&self) -> OrbitPartial

Returns the right ascension of ther ascending node in degrees

source

pub fn ta(&self) -> OrbitPartial

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

source

pub fn tlong(&self) -> OrbitPartial

Returns the true longitude in degrees

source

pub fn aol(&self) -> OrbitPartial

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) -> OrbitPartial

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

source

pub fn apoapsis(&self) -> OrbitPartial

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

source

pub fn ea(&self) -> OrbitPartial

Returns the eccentric anomaly in degrees

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

source

pub fn fpa(&self) -> OrbitPartial

Returns the flight path angle in degrees

source

pub fn ma(&self) -> OrbitPartial

Returns the mean anomaly in degrees

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

source

pub fn semi_parameter(&self) -> OrbitPartial

Returns the semi parameter (or semilatus rectum)

source

pub fn geodetic_longitude(&self) -> OrbitPartial

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) -> OrbitPartial

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) -> OrbitPartial

Returns the geodetic height in km.

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

source

pub fn right_ascension(&self) -> OrbitPartial

Returns the right ascension of this orbit in degrees

source

pub fn declination(&self) -> OrbitPartial

Returns the declination of this orbit in degrees

source

pub fn semi_minor_axis(&self) -> OrbitPartial

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

source

pub fn velocity_declination(&self) -> OrbitPartial

Returns the velocity declination of this orbit in degrees

source

pub fn c3(&self) -> OrbitPartial

Returns the $C_3$ of this orbit

source

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

Returns the hyperbolic anomaly in degrees between 0 and 360.0

Trait Implementations§

source§

impl Clone for OrbitDual

source§

fn clone(&self) -> OrbitDual

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 OrbitDual

source§

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

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

impl From<Orbit> for OrbitDual

source§

fn from(orbit: Orbit) -> Self

Initialize a new OrbitDual from an orbit, no other initializers

source§

impl TimeTagged for OrbitDual

source§

fn epoch(&self) -> Epoch

Retrieve the Epoch
source§

fn set_epoch(&mut self, epoch: Epoch)

Set the Epoch
source§

fn shift_by(&mut self, duration: Duration)

Shift this epoch by a duration (can be negative)
source§

impl Copy for OrbitDual

Auto Trait Implementations§

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
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, 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