Astro  0.0.0
A C++ library for space dynamics
Loading...
Searching...
No Matches
Astro::OrbitalElements::Equinoctial Struct Reference

Struct container for the (modified) equinoctial orbital elements. More...

#include <OrbitalElements.hh>

Public Member Functions

 Equinoctial ()
 Equinoctial (Real const t_p, Real const t_f, Real const t_g, Real const t_h, Real const t_k)
 Equinoctial (Vector5 const &t_equi)
 Equinoctial (Equinoctial const &)=default
 Equinoctial (Equinoctial &&)=default
Equinoctialoperator= (const Equinoctial &)=default
Equinoctialoperator= (Equinoctial &&)=default
Vector5 vector () const
std::string info () const
void info (std::ostream &os)
void reset ()
bool sanity_check () const
bool is_singular (Real tol_i=EPSILON_LOW) const
bool is_nonsingular (Real tol_i=EPSILON_LOW) const
Real u (Real L) const
Real a () const
Real e () const
Real i () const
Real omega () const
Real Omega () const

Public Attributes

Real p {QUIET_NAN}
Real f {QUIET_NAN}
Real g {QUIET_NAN}
Real h {QUIET_NAN}
Real k {QUIET_NAN}

Detailed Description

The (modified) equinoctial orbital elements are defined as follows:

  • the semi-latus rectum \( p \) (UA).
  • the \( x \)-axis component of the eccentricity vector in the orbital frame \( f \) (-).
  • the \( y \)-axis component of the eccentricity vector in the orbital frame \( g \) (-).
  • the \( x \)-axis component of the node vector in the orbital frame \( h \) (-).
  • the \( y \)-axis component of the node vector in the orbital frame \( k \) (-).
  • the posigrade (+1)/retrograde (-1) factor \( I = +1 \) (posigrade) or \( I = -1 \) (retrograde) (-). Singular at \( i = \pi \).
Note
For more information on orbital elements, refer to "Survey of Orbital Elements", by G. R. Hintz, Journal of Guidance, Control, and Dynamics, Vol. 31, No. 3, May-June 2008.

Constructor & Destructor Documentation

◆ Equinoctial() [1/5]

Astro::OrbitalElements::Equinoctial::Equinoctial ( )
inline

Struct constructor for the (modified) equinoctial orbit parameters.

◆ Equinoctial() [2/5]

Astro::OrbitalElements::Equinoctial::Equinoctial ( Real const t_p,
Real const t_f,
Real const t_g,
Real const t_h,
Real const t_k )
inline

Struct constructor for the (modified) equinoctial orbit parameters.

Parameters
[in]t_pThe semi-latus rectum \( p \).
[in]t_fThe \( x \)-axis component of the eccentricity vector in the orbital frame.
[in]t_gThe \( y \)-axis component of the eccentricity vector in the orbital frame.
[in]t_hThe \( x \)-axis component of the node vector in the orbital frame.
[in]t_kThe \( y \)-axis component of the node vector in the orbital frame.

◆ Equinoctial() [3/5]

Astro::OrbitalElements::Equinoctial::Equinoctial ( Vector5 const & t_equi)
inline

Struct constructor for the (modified) equinoctial orbit parameters.

Parameters
[in]t_equiThe vector of equinoctial orbit parameters \( [p, f, g, h, k]^\top \).
Note
The vector \( t_{\text{equinoctial}} \) must have 5 elements.

◆ Equinoctial() [4/5]

Astro::OrbitalElements::Equinoctial::Equinoctial ( Equinoctial const & )
default

Enable the default equinoctial orbit parameters copy constructor.

◆ Equinoctial() [5/5]

Astro::OrbitalElements::Equinoctial::Equinoctial ( Equinoctial && )
default

Enable the default equinoctial orbit parameters move constructor.

Member Function Documentation

◆ a()

Real Astro::OrbitalElements::Equinoctial::a ( ) const
inline

Compute the semi-major axis \( a = \displaystyle\frac{p}{1 - f^2 - g^2} \) (UA).

Returns
The semi-major axis \( a \).

◆ e()

Real Astro::OrbitalElements::Equinoctial::e ( ) const
inline

Compute the eccentricity \( \mathbf{e} = [f, g, 0]^\top \) (-).

Returns
The eccentricity \( \mathbf{e} \).

◆ i()

Real Astro::OrbitalElements::Equinoctial::i ( ) const
inline

Compute the inclination \( i = \arctan\left(2\sqrt{h^2 + k^2}, 1 - h^2 - k^2\right) \) (rad).

Returns
The inclination \( i \).

◆ info() [1/2]

std::string Astro::OrbitalElements::Equinoctial::info ( ) const
inline

Print the equinoctial orbit parameters on a string.

Returns
Equinoctial orbit parameters string.

◆ info() [2/2]

void Astro::OrbitalElements::Equinoctial::info ( std::ostream & os)
inline

Print the equinoctial orbit parameters on a stream.

Parameters
[in,out]osOutput stream.

◆ is_nonsingular()

bool Astro::OrbitalElements::Equinoctial::is_nonsingular ( Real tol_i = EPSILON_LOW) const
inline

Check if the equinoctial orbit is singular, i.e., * \( \pi-\varepsilon_i < i < \pi+\varepsilon_i\).

Parameters
[in]tol_iTolerance \( \varepsilon_i \) for the singularity check on the inclination.
Returns
True if the equinoctial orbit is nonsingular, false otherwise.

◆ is_singular()

bool Astro::OrbitalElements::Equinoctial::is_singular ( Real tol_i = EPSILON_LOW) const
inline

Check if the equinoctial orbit is singular, i.e., * \( i < \pi-\varepsilon_i \), or \( i > \pi+\varepsilon_i\).

Parameters
[in]tol_iTolerance \( \varepsilon_i \) for the singularity check on the inclination.
Returns
True if the equinoctial orbit is singular, false otherwise.

◆ Omega()

Real Astro::OrbitalElements::Equinoctial::Omega ( ) const
inline

Compute the right ascension of the ascending node \( \Omega = \arctan\left(k, h\right) \) (rad).

Returns
The right ascension of the ascending node \( \Omega \).

◆ omega()

Real Astro::OrbitalElements::Equinoctial::omega ( ) const
inline

Compute the argument of periapsis \( \omega = \arctan\left(gh - fk, fh + gk\right) \) (rad).

Returns
The argument of periapsis \( \omega \).

◆ operator=() [1/2]

Equinoctial & Astro::OrbitalElements::Equinoctial::operator= ( const Equinoctial & )
default

Enable the default equinoctial orbit parameters assignment operator.

◆ operator=() [2/2]

Equinoctial & Astro::OrbitalElements::Equinoctial::operator= ( Equinoctial && )
default

Enable the default equinoctial orbit parameters move assignment operator.

◆ reset()

void Astro::OrbitalElements::Equinoctial::reset ( )
inline

Reset the equinoctial orbit parameters to NaN.

◆ sanity_check()

bool Astro::OrbitalElements::Equinoctial::sanity_check ( ) const
inline

Check if the equinoctial orbit parameters are valid, i.e., finite, and with \( p > 0 \).

Returns
True if the equinoctial orbit parameters are valid, false otherwise.

◆ u()

Real Astro::OrbitalElements::Equinoctial::u ( Real L) const
inline

Compute the argument of latitude \( u = \omega + \Omega = \arctan\left(h\sin(L) - k\cos(L), h\cos(L) + k\sin(L)\right) \).

Parameters
[in]LThe true longitude \( L \).
Returns
The argument of latitude \( u \).

◆ vector()

Vector5 Astro::OrbitalElements::Equinoctial::vector ( ) const
inline

Get the equinoctial orbit parameters as a vector.

Returns
The vector of equinoctial orbit parameters \( [p, f, g, h, k]^\top \).

Member Data Documentation

◆ f

Real Astro::OrbitalElements::Equinoctial::f {QUIET_NAN}

\( X \)-axis component of the eccentricity vector in the orbital frame (-).

◆ g

Real Astro::OrbitalElements::Equinoctial::g {QUIET_NAN}

\( Y \)-axis component of the eccentricity vector in the orbital frame (-).

◆ h

Real Astro::OrbitalElements::Equinoctial::h {QUIET_NAN}

\( X \)-axis component of the node vector in the orbital frame (-).

◆ k

Real Astro::OrbitalElements::Equinoctial::k {QUIET_NAN}

\( Y \)-axis component of the node vector in the orbital frame (-).

◆ p

Real Astro::OrbitalElements::Equinoctial::p {QUIET_NAN}

Semi-latus rectum \( p \) (UA).


The documentation for this struct was generated from the following file: