|
Astro
0.0.0
A C++ library for space dynamics
|
Orbit class container. More...
#include <Orbit.hh>
Public Member Functions | |
| Orbit () | |
| Orbit (Orbit const &)=default | |
| Orbit (Orbit &&)=default | |
| Orbit & | operator= (const Orbit &)=default |
| Orbit & | operator= (Orbit &&)=default |
| Cartesian const & | cartesian () const |
| void | set_cartesian (Real r_x, Real r_y, Real r_z, Real v_x, Real v_y, Real v_z) |
| void | set_cartesian (Vector3 const &t_r, Vector3 const &t_v) |
| void | set_cartesian (Vector6 const &t_cart) |
| void | set_cartesian (Cartesian const &t_cart) |
| Keplerian const & | keplerian () const |
| void | set_keplerian (Real const t_a, Real const t_e, Real const t_i, Real const t_Omega, Real const t_omega, Real const nu) |
| void | set_keplerian (Vector5 const &t_kepl, Real const nu) |
| void | set_keplerian (Keplerian const &t_kepl, Real const nu) |
| Equinoctial const & | equinoctial () const |
| void | set_equinoctial (Real const t_p, Real const t_f, Real const t_g, Real const t_h, Real const t_k, Real const L) |
| void | set_equinoctial (Vector5 const &t_equi, Real const L) |
| void | set_equinoctial (Equinoctial const &t_equi, Real const L) |
| Quaternionic const & | quaternionic () const |
| void | set_quaternionic (Real const t_q_1, Real const t_q_2, Real const t_q_3, Real const t_q_4) |
| void | set_quaternionic (Quaternion const &t_quat) |
| Type | type () const |
| void | set_type (Real const e) |
| Factor | factor () const |
| void | set_factor (Factor t_factor) |
| Real | mu () const |
| void | set_mu (Real const t_mu) |
| std::string | info () const |
| void | info (std::ostream &os) |
| void | reset () |
| bool | sanity_check () const |
| Rotation | keplerian_to_reference () const |
| Rotation | equinoctial_to_reference () const |
| Vector3 | cartesian_rtn_to_xyz (Vector3 const &r, Vector3 const &v, Vector3 const &vec) const |
| Vector3 | equinoctial_rtn_to_xyz (Vector3 const &vec) const |
Protected Attributes | |
| Cartesian | m_cart |
| Keplerian | m_kepl |
| Equinoctial | m_equi |
| Quaternionic | m_quat |
| Type | m_type {Type::UNDEFINED} |
| Factor | m_factor {Factor::UNDEFINED} |
| Real | m_mu {QUIET_NAN} |
Private Types | |
| using | Cartesian = OrbitalElements::Cartesian |
| using | Keplerian = OrbitalElements::Keplerian |
| using | Equinoctial = OrbitalElements::Equinoctial |
| using | Quaternionic = OrbitalElements::Quaternionic |
| using | Anomaly = OrbitalElements::Anomaly |
The Orbit class container is used to store the orbital elements of an object in space. The class provides methods for setting and getting the orbital elements in various representations, such as cartesian, keplerian, equinoctial, and quaternionic. The class also provides methods for converting between these representations.
|
private |
|
private |
|
private |
|
private |
|
private |
|
inline |
Class constructor for the astro object.
|
default |
Enable the default astro copy constructor.
|
default |
Enable the default astro move constructor.
|
inline |
Get the cartesian orbital elements.
|
inline |
Compute the Frenet-Serret frame of the orbit (radial, tangential, normal) given the cartesian orbital elements.
| [in] | r | The position vector \( \mathbf{r} \). |
| [in] | v | The velocity vector \( \mathbf{v} \f $. \return The Frenet-Serret frame of the orbit. */ Rotation cartesian_to_frenet_rtn(Vector3 const & r, Vector3 const & v) const { // Check if the input vectors are valid ASTRO_ASSERT(r.allFinite() && v.allFinite(), "Astro::Orbit::cartesian_to_frenet_rtn(...): invalid input vectors detected."); // Initialize the Frenet-Serret frame Rotation rtn; // Compute the radial vector rtn.col(0) = r; rtn.col(0).normalize(); // Compute the binormal vector rtn.col(2) = r.cross(v); rtn.col(2).normalize(); // Compute the tangential vector rtn.col(1) = rtn.col(2).cross(r); rtn.col(1).normalize(); // Return the Frenet-Serret frame return rtn; } /** Compute the Frenet-Serret frame of the orbit (radial, tangential, normal) through the equinoctial orbital elements. \param[in] L The true longitude \) L \(. \return The Frenet-Serret frame of the orbit. */ Rotation equinoctial_to_frenet_rtn(Real const L) const { Real h{this->m_equi.h}; Real k{this->m_equi.k}; Real c_L{std::cos(L)}; Real s_L{std::sin(L)}; Real h2{h*h}; Real k2{k*k}; Real hk{h*k}; Real bf{1.0+h2+k2}; Real I{static_cast<Real>(this->m_factor)}; Rotation R; R << ((h2-k2+1)*c_L+2.0*I*hk*s_L)/bf, ((k2-h2-1.0)*s_L+2.0*hk*I*c_L)/bf, 2.0*k/bf, (2.0*hk*c_L-I*(h2-k2-1.0)*s_L)/bf, (I*(1.0+k2-h2)*c_L-2.0*hk*s_L)/bf, -2.0*h/bf, (2.0*(h*s_L-c_L*k*I))/bf, (2.0*k*I*s_L+2.0*h*c_L)/bf, I*(1-h2-k2)/bf; return R; } /** Transform a vector in the Frent-Serret frame of the orbit (radial, tangential, normal) to a vector in the cartesian frame given the cartesian orbital elements. \param[in] r The position vector \) \mathbf{r} \(. \param[in] v The velocity vector \) \mathbf{v} \f $. |
| [in] | vec | The vector in Frenet-Serret frame. |
|
inline |
Get the equinoctial orbital elements.
Transform a vector in the Frent-Serret frame of the orbit (radial, tangential, normal) to a vector in the cartesian frame through the equinoctial orbital elements.
| [in] | vec | The vector in Frenet-Serret frame. |
|
inline |
Compute the rotation matrix of the orbital plane through the equinoctial orbital elements.
|
inline |
Get the posigrade (+1)/retrograde (-1) factor \( I \).
|
inline |
Print the orbit information on a string.
|
inline |
Print the orbit information on a stream.
| [in,out] | os | Output stream. |
|
inline |
Get the keplerian orbital elements.
|
inline |
Compute the rotation matrix of the orbital plane through the keplerian orbital elements.
|
inline |
Get the gravitational constant of the central body.
Enable the default astro assignment operator.
Enable the default astro move assignment operator.
|
inline |
Get the quaternionic orbital elements.
|
inline |
Reset the orbit parameters to undefined values.
|
inline |
Check if the orbit can accept the cartesian orbital elements, i.e., if factor, type, an the gravitational constant are set.
|
inline |
Set the cartesian orbital elements.
| [in] | t_cart | The cartesian orbital elements. |
|
inline |
Set the cartesian orbital elements.
| [in] | r_x | Position vector \( x \)-axis component. |
| [in] | r_y | Position vector \( y \)-axis component. |
| [in] | r_z | Position vector \( z \)-axis component. |
| [in] | v_x | Velocity vector \( x \)-axis component. |
| [in] | v_y | Velocity vector \( y \)-axis component. |
| [in] | v_z | Velocity vector \( z \)-axis component. |
Set the cartesian orbital elements.
| [in] | t_r | The position vector \( \mathbf{r} \). |
| [in] | t_v | The velocity vector \( \mathbf{v} \). |
|
inline |
Set the cartesian orbital elements.
| [in] | t_cart | The vector of cartesian orbital elements \( [r_x, r_y, r_z, v_x, v_y, v_z]^\top \). |
|
inline |
Set the equinoctial orbital elements.
| [in] | t_equi | The equinoctial orbital elements. |
| [in] | L | The true longitude \( L \). |
|
inline |
Set the (modified) equinoctial orbit 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. |
| [in] | L | The true longitude \( L \). |
Set the (modified) equinoctial orbit parameters.
| [in] | t_equi | The vector of equinoctial orbit parameters \( [p, f, g, h, k]^\top \). |
| [in] | L | The true longitude \( L \). |
|
inline |
Set the posigrade (+1)/retrograde (-1) factor \( I \).
| [in] | t_factor | The posigrade (+1)/retrograde (-1) factor \( I \). |
Set the keplerian orbital elements.
| [in] | t_kepl | The keplerian orbital elements. |
| [in] | nu | The true anomaly \( \nu \). |
|
inline |
Set the keplerian orbital elements.
| [in] | t_a | The semi-major axis \( a \). |
| [in] | t_e | The eccentricity \( e \). |
| [in] | t_i | The inclination \( i \). |
| [in] | t_Omega | The longitude of the ascending node \( \Omega \). |
| [in] | t_omega | The argument of periapsis \( \omega \). |
| [in] | nu | The true anomaly \( \nu \). |
Set the keplerian orbital elements.
| [in] | t_kepl | The vector of keplerian orbital elements \( [a, e, i, \Omega, \omega]^\top \). |
| [in] | nu | The true anomaly \( \nu \). |
|
inline |
Set the gravitational constant of the central body.
| [in] | t_mu | The gravitational constant of the central body. |
|
inline |
Set the quaternionic orbital elements.
| [in] | t_quat | The quaternion vector. |
|
inline |
Set the quaternionic orbital elements.
| [in] | t_q_1 | The first quaternionic orbit parameter. |
| [in] | t_q_2 | The second quaternionic orbit parameter. |
| [in] | t_q_3 | The third quaternionic orbit parameter. |
| [in] | t_q_4 | The fourth quaternionic orbit parameter. |
|
inline |
Set the type of the orbit.
| [in] | e | The eccentricity of the orbit. |
|
inline |
Get the type of the orbit.
|
protected |
Equinoctial orbit parameters.
|
protected |
Orbit posigrade (+1)/retrograde (-1) factor.
|
protected |
Quaternionic orbit parameters.