PCL
|
Osculating conic orbital elements. More...
#include <OsculatingElements.h>
Public Member Functions | |
OsculatingElements ()=default | |
OsculatingElements (const OsculatingElements &)=default | |
OsculatingElements (const Vector &r, const Vector &v, const TimePoint &t, double GM=0, double m=0) | |
bool | IsNearParabolic () const |
double | MeanAnomalyFromTimeOfPerifocalPassage (const TimePoint &t, double GM=0, double m=0) const |
double | MeanAnomalyFromTimeOfPerihelionPassage (const TimePoint &t, double GM=0, double m=0) const |
double | MeanMotion (double GM=0, double m=0) const |
double | OrbitalPeriod (double GM=0, double m=0) const |
double | PerifocalDistanceFromSemimajorAxis () const |
double | PerihelionDistanceFromSemimajorAxis () const |
double | SemimajorAxisFromPerifocalDistance () const |
double | SemimajorAxisFromPerihelionDistance () const |
TimePoint | TimeOfPerifocalPassageFromMeanAnomaly (const TimePoint &t, double GM=0, double m=0) const |
TimePoint | TimeOfPerihelionPassageFromMeanAnomaly (const TimePoint &t, double GM=0, double m=0) const |
void | ToPerifocalPosition (double &xp, double &yp, const TimePoint &t, double GM=0, double m=0) const |
void | ToPerifocalPositionAndVelocity (double &xp, double &yp, double &xv, double &yv, const TimePoint &t, double GM=0, double m=0) const |
void | ToStateVectors (Vector &r, Vector &v, const TimePoint &t, double GM=0, double m=0) const |
Static Public Member Functions | |
static OsculatingElements | FromStateVectors (const Vector &r, const Vector &v, const TimePoint &t, double GM=0, double m=0) |
static void | GetOrbitOrientationFromStateVectors (double &i, double &O, double &w, const Vector &r, const Vector &v, double GM=0, double m=0) |
Public Attributes | |
double | a = 0.0 |
Semimajor axis. Unused for (near-)parabolic orbits. | |
double | e = 0.0 |
Eccentricity. | |
double | i = 0.0 |
Inclination over the reference plane (radians). | |
double | M = 0.0 |
Mean anomaly (radians). | |
double | O = 0.0 |
Longitude of the ascending node (radians). Undefined if i = 0. | |
double | q = 0.0 |
Perifocal distance. Required for (near-)parabolic orbits. | |
TimePoint | T = 0.0 |
Time of perifocal passage. Undefined for circular orbits. | |
double | w = 0.0 |
Argument of periapsis (radians). Undefined for circular orbits or when i = 0. | |
This class represents a set of osculating conic orbital elements that define the size and orientation of an orbit in space at a given epoch.
The length and time units are implicitly defined by the gravitational mass parameters of the central body used upon construction from state vectors (position and velocity vectors). If the default parameters are used, a heliocentric orbit is implicitly assumed where the units of length and time are the astronomical unit (au) and the day of 86400 SI seconds in the TDB timescale, respectively. All data members representing angular elements (mean and true anomalies, inclination, longitude of the ascending node and argument of periapsis) are stored in radians.
The orbital elements describe the orbit of a celestial body about a dynamical center located at the origin of an intertial rectangular coordinate system. The same elements are used to describe the three types of conic orbits: elliptic, hyperbolic, and (near-)parabolic.
Definition at line 90 of file OsculatingElements.h.
|
default |
Default constructor. Constructs an uninitialized object where all orbital elements are set to zero.
|
default |
Copy constructor.
pcl::OsculatingElements::OsculatingElements | ( | const Vector & | r, |
const Vector & | v, | ||
const TimePoint & | t, | ||
double | GM = 0 , |
||
double | m = 0 |
||
) |
Constructor from state vectors.
r | Position vector. |
v | Velocity vector. |
t | Epoch of state vectors (TDB). |
GM | Gravitational mass parameter of the central body in units of length^3/time^2. If not specified or ≤ 0, a default value for the Sun in au^3/day^2 units, consistent with current fundamental ephemerides, will be used. |
m | Mass of the body in central body mass units. Zero by default. |
|
inlinestatic |
Returns a new OsculatingElements object with orbital elements calculated from the specified state vectors.
r | Position vector. |
v | Velocity vector. |
t | Epoch of state vectors (TDB). |
GM | Gravitational mass parameter of the central body in units of length^3/time^2. If not specified or ≤ 0, a default value for the Sun in au^3/day^2 units, consistent with current fundamental ephemerides, will be used. |
m | Mass of the body in central body mass units. Zero by default. |
Definition at line 161 of file OsculatingElements.h.
|
static |
Computes the orbit orientation angles from given state vectors.
[out] | i | Inclination over the reference plane (radians). |
[out] | O | Longitude of the ascending node (radians). |
[out] | w | Argument of periapsis (radians). |
r | Position vector. | |
v | Velocity vector. | |
GM | Gravitational mass parameter of the central body in units of length^3/time^2. If not specified or ≤ 0, a default value for the Sun in au^3/day^2 units, consistent with current fundamental ephemerides, will be used. | |
m | Mass of the body in central body mass units. Zero by default. |
This function is useful when only the orbit orientation in space is required, saving the calculation of a complete set of orbital elements.
|
inline |
Returns true iff this object defines a parabolic or near-parabolic orbit.
With the current implementation, an orbit is considered near-parabolic if its eccentricity is > 0.999.
For parabolic and near-parabolic orbits, we expect to have the perifocal distance q instead of the semimajor axis a.
Definition at line 199 of file OsculatingElements.h.
References pcl::Abs().
|
inline |
The mean anomaly calculated from the time of perifocal passage.
t | Epoch of osculating elements, TDB. |
GM | Gravitational mass parameter of the central body in units of length^3/time^2. If not specified or ≤ 0, a default value for the Sun in au^3/day^2 units, consistent with current fundamental ephemerides, will be used. |
m | Mass of the body in central body mass units. Zero by default. |
Returns the mean anomaly in radians, in the range [0,2pi).
Definition at line 249 of file OsculatingElements.h.
References pcl::Norm2Pi().
|
inline |
The mean anomaly calculated from the time of perihelion passage.
This function is a convenience synonym for MeanAnomalyFromTimeOfPerifocalPassage().
Definition at line 260 of file OsculatingElements.h.
double pcl::OsculatingElements::MeanMotion | ( | double | GM = 0 , |
double | m = 0 |
||
) | const |
Mean orbital motion.
GM | Gravitational mass parameter of the central body in units of length^3/time^2. If not specified or ≤ 0, a default value for the Sun in au^3/day^2 units, consistent with current fundamental ephemerides, will be used. |
m | Mass of the body in central body mass units. Zero by default. |
Returns the mean motion in radians, computed from the semimajor axis for elliptical and hyperbolic orbits, or from the perifocal distance for near-parabolic orbits.
double pcl::OsculatingElements::OrbitalPeriod | ( | double | GM = 0 , |
double | m = 0 |
||
) | const |
Orbital period of the elliptic orbit.
GM | Gravitational mass parameter of the central body in units of length^3/time^2. If not specified or ≤ 0, a default value for the Sun in au^3/day^2 units, consistent with current fundamental ephemerides, will be used. |
m | Mass of the body in central body mass units. Zero by default. |
Returns the orbital period of an elliptic orbit (e < 1) in time units. For parabolic and hyperbolic orbits, this function returns zero.
|
inline |
The perifocal distance calculated from semimajor axis and eccentricity.
Definition at line 298 of file OsculatingElements.h.
|
inline |
The perihelion distance calculated from semimajor axis and eccentricity.
This function is a convenience synonym for PerifocalDistanceFromSemimajorAxis().
Definition at line 309 of file OsculatingElements.h.
|
inline |
The semimajor axis calculated from perifocal distance and eccentricity.
Definition at line 317 of file OsculatingElements.h.
References pcl::Abs().
|
inline |
The semimajor axis calculated from perihelion distance and eccentricity.
This function is a convenience synonym for SemimajorAxisFromPerifocalDistance().
Definition at line 331 of file OsculatingElements.h.
|
inline |
The time of perifocal passage calculated from the mean anomaly.
t | Epoch of osculating elements, TDB. |
GM | Gravitational mass parameter of the central body in units of length^3/time^2. If not specified or ≤ 0, a default value for the Sun in au^3/day^2 units, consistent with current fundamental ephemerides, will be used. |
m | Mass of the body in central body mass units. Zero by default. |
Returns the time of perifocal passage in the TDB timescale.
Definition at line 279 of file OsculatingElements.h.
|
inline |
The time of perihelion passage calculated from the mean anomaly.
This function is a convenience synonym for TimeOfPerifocalPassageFromMeanAnomaly().
Definition at line 290 of file OsculatingElements.h.
|
inline |
Position components in the perifocal reference frame.
[out] | xp,yp | Perifocal position vector components. |
t | Epoch of osculating elements (TDB). | |
GM | Gravitational mass parameter of the central body in units of length^3/time^2. If not specified or ≤ 0, a default value for the Sun in au^3/day^2 units, consistent with current fundamental ephemerides, will be used. | |
m | Mass of the body in central body mass units. Zero by default. |
The output position components are referred to a right-handed Cartesian coordinate system where the XY plane coincides with the orbit plane and the X axis is aligned with the semimajor axis, positive towards the perifocus. The Z axis points in the same direction as the angular momentum vector. This system is known as the perifocal reference frame.
Definition at line 386 of file OsculatingElements.h.
void pcl::OsculatingElements::ToPerifocalPositionAndVelocity | ( | double & | xp, |
double & | yp, | ||
double & | xv, | ||
double & | yv, | ||
const TimePoint & | t, | ||
double | GM = 0 , |
||
double | m = 0 |
||
) | const |
Position and velocity components in the perifocal reference frame.
[out] | xp,yp | Perifocal position vector components. |
[out] | xv,yv | Perifocal velocity vector components. |
t | Epoch of osculating elements (TDB). | |
GM | Gravitational mass parameter of the central body in units of length^3/time^2. If not specified or ≤ 0, a default value for the Sun in au^3/day^2 units, consistent with current fundamental ephemerides, will be used. | |
m | Mass of the body in central body mass units. Zero by default. |
The output position and velocity components are referred to a right-handed Cartesian coordinate system where the XY plane coincides with the orbit plane and the X axis is aligned with the semimajor axis, positive towards the perifocus. The Z axis points in the same direction as the angular momentum vector. This system is known as the perifocal reference frame.
void pcl::OsculatingElements::ToStateVectors | ( | Vector & | r, |
Vector & | v, | ||
const TimePoint & | t, | ||
double | GM = 0 , |
||
double | m = 0 |
||
) | const |
Computes position and velocity vectors from conic orbital elements.
[out] | r | Position vector. |
[out] | v | Velocity vector. |
t | Epoch of state vectors (TDB). | |
GM | Gravitational mass parameter of the central body in units of length^3/time^2. If not specified or ≤ 0, a default value for the Sun in au^3/day^2 units, consistent with current fundamental ephemerides, will be used. | |
m | Mass of the body in central body mass units. Zero by default. |