Struct container for the (modified) equinoctial orbital elements.
More...
#include <OrbitalElements.hxx>
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.
◆ Equinoctial() [1/4]
Astro::OrbitalElements::Equinoctial::Equinoctial |
( |
| ) |
|
|
inline |
Struct constructor for the (modified) equinoctial orbit parameters.
◆ Equinoctial() [2/4]
Astro::OrbitalElements::Equinoctial::Equinoctial |
( |
Real | t_p, |
|
|
Real | t_f, |
|
|
Real | t_g, |
|
|
Real | t_h, |
|
|
Real | t_k ) |
|
inline |
Struct constructor for the (modified) equinoctial orbit parameters.
- Parameters
-
[in] | t_p | The semi-latus rectum \( p \). |
[in] | t_f | The \( x \)-axis component of the eccentricity vector in the orbital frame. |
[in] | t_g | The \( y \)-axis component of the eccentricity vector in the orbital frame. |
[in] | t_h | The \( x \)-axis component of the node vector in the orbital frame. |
[in] | t_k | The \( y \)-axis component of the node vector in the orbital frame. |
◆ Equinoctial() [3/4]
Astro::OrbitalElements::Equinoctial::Equinoctial |
( |
Equinoctial const & | | ) |
|
|
default |
Enable the default equinoctial orbit parameters copy constructor.
◆ Equinoctial() [4/4]
Astro::OrbitalElements::Equinoctial::Equinoctial |
( |
Equinoctial && | | ) |
|
|
default |
Enable the default equinoctial orbit parameters move constructor.
◆ 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
-
◆ 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_i | Tolerance \( \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_i | Tolerance \( \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]
Enable the default equinoctial orbit parameters assignment operator.
◆ operator=() [2/2]
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] | L | The true longitude \( L \). |
- Returns
- The argument of latitude \( u \).
\( X \)-axis component of the eccentricity vector in the orbital frame (-).
\( Y \)-axis component of the eccentricity vector in the orbital frame (-).
\( X \)-axis component of the node vector in the orbital frame (-).
\( Y \)-axis component of the node vector in the orbital frame (-).
Semi-latus rectum \( p \) (UA).
The documentation for this struct was generated from the following file: