Struct nyx_space::celestia::OrbitDual [−][src]
pub struct OrbitDual { pub x: Hyperdual<f64, U7>, pub y: Hyperdual<f64, U7>, pub z: Hyperdual<f64, U7>, pub vx: Hyperdual<f64, U7>, pub vy: Hyperdual<f64, U7>, pub vz: Hyperdual<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: Hyperdual<f64, U7>
in km
y: Hyperdual<f64, U7>
in km
z: Hyperdual<f64, U7>
in km
vx: Hyperdual<f64, U7>
in km/s
vy: Hyperdual<f64, U7>
in km/s
vz: Hyperdual<f64, U7>
in km/s
dt: Epoch
frame: Frame
Frame contains everything we need to compute state information
Implementations
Returns the magnitude of the radius vector in km
Returns the magnitude of the velocity vector in km/s
Returns the orbital momentum value on the X axis
Returns the orbital momentum value on the Y axis
Returns the orbital momentum value on the Z axis
Returns the norm of the orbital momentum
Returns the specific mechanical energy
Returns the semi-major axis in km
Returns the eccentricity (no unit)
Returns the inclination in degrees
Returns the argument of periapsis in degrees
Returns the right ascension of ther ascending node in degrees
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
Returns the true longitude in degrees
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.
Returns the radius of periapsis (or perigee around Earth), in kilometers.
Returns the radius of apoapsis (or apogee around Earth), in kilometers.
Returns the eccentric anomaly in degrees
This is a conversion from GMAT’s StateConversionUtil::TrueToEccentricAnomaly
Returns the flight path angle in degrees
Returns the mean anomaly in degrees
This is a conversion from GMAT’s StateConversionUtil::TrueToMeanAnomaly
Returns the semi parameter (or semilatus rectum)
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
Returns the geodetic latitude (φ) in degrees. Value is between -180 and +180 degrees.
Reference: Vallado, 4th Ed., Algorithm 12 page 172.
Returns the geodetic height in km.
Reference: Vallado, 4th Ed., Algorithm 12 page 172.
Returns the right ascension of this orbit in degrees
Returns the declination of this orbit in degrees
Returns the semi minor axis in km, includes code for a hyperbolic orbit
Returns the velocity declination of this orbit in degrees
Returns the $C_3$ of this orbit
Returns the hyperbolic anomaly in degrees between 0 and 360.0
Trait Implementations
Auto Trait Implementations
impl RefUnwindSafe for OrbitDual
impl UnwindSafe for OrbitDual
Blanket Implementations
Mutably borrows from an owned value. Read more
type Output = T
type Output = T
Should always be Self
The inverse inclusion map: attempts to construct self
from the equivalent element of its
superset. Read more
pub fn is_in_subset(&self) -> bool
pub fn is_in_subset(&self) -> bool
Checks if self
is actually part of its subset T
(and can be converted to it).
pub fn to_subset_unchecked(&self) -> SS
pub fn to_subset_unchecked(&self) -> SS
Use with care! Same as self.to_subset
but without any property checks. Always succeeds.
pub fn from_subset(element: &SS) -> SP
pub fn from_subset(element: &SS) -> SP
The inclusion map: converts self
to the equivalent element of its superset.
pub fn vzip(self) -> V