PCL
pcl::OsculatingElements Class Reference

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.
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ OsculatingElements() [1/3]

pcl::OsculatingElements::OsculatingElements ( )
default

Default constructor. Constructs an uninitialized object where all orbital elements are set to zero.

◆ OsculatingElements() [2/3]

pcl::OsculatingElements::OsculatingElements ( const OsculatingElements )
default

Copy constructor.

◆ OsculatingElements() [3/3]

pcl::OsculatingElements::OsculatingElements ( const Vector r,
const Vector v,
const TimePoint t,
double  GM = 0,
double  m = 0 
)

Constructor from state vectors.

Parameters
rPosition vector.
vVelocity vector.
tEpoch of state vectors (TDB).
GMGravitational 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.
mMass of the body in central body mass units. Zero by default.
Note
The calculated orbital elements will describe an orbit in the same inertial rectangular coordinate system where the specified r and v vectors are defined. For example, if the state vectors provide heliocentric ecliptic rectangular coordiantes, the orbit will be referred to the same fundamental plane (ecliptic) and origin (the Sun).

Member Function Documentation

◆ FromStateVectors()

static OsculatingElements pcl::OsculatingElements::FromStateVectors ( const Vector r,
const Vector v,
const TimePoint t,
double  GM = 0,
double  m = 0 
)
inlinestatic

Returns a new OsculatingElements object with orbital elements calculated from the specified state vectors.

Parameters
rPosition vector.
vVelocity vector.
tEpoch of state vectors (TDB).
GMGravitational 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.
mMass of the body in central body mass units. Zero by default.
Note
The calculated orbital elements will describe an orbit in the same inertial rectangular coordinate system where the specified r and v vectors are defined. For example, if the state vectors provide heliocentric ecliptic rectangular coordiantes, the orbit will be referred to the same fundamental plane (ecliptic) and origin (the Sun).

Definition at line 161 of file OsculatingElements.h.

◆ GetOrbitOrientationFromStateVectors()

static void pcl::OsculatingElements::GetOrbitOrientationFromStateVectors ( double &  i,
double &  O,
double &  w,
const Vector r,
const Vector v,
double  GM = 0,
double  m = 0 
)
static

Computes the orbit orientation angles from given state vectors.

Parameters
[out]iInclination over the reference plane (radians).
[out]OLongitude of the ascending node (radians).
[out]wArgument of periapsis (radians).
rPosition vector.
vVelocity vector.
GMGravitational 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.
mMass 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.

Note
The output orbital elements will describe the orientation of an orbit in the same inertial rectangular coordinate system where the specified r and v vectors are defined. For example, if the state vectors provide heliocentric ecliptic rectangular coordiantes, the output orientation angles will be referred to the same fundamental plane (ecliptic) and origin (the Sun).

◆ IsNearParabolic()

bool pcl::OsculatingElements::IsNearParabolic ( ) const
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().

◆ MeanAnomalyFromTimeOfPerifocalPassage()

double pcl::OsculatingElements::MeanAnomalyFromTimeOfPerifocalPassage ( const TimePoint t,
double  GM = 0,
double  m = 0 
) const
inline

The mean anomaly calculated from the time of perifocal passage.

Parameters
tEpoch of osculating elements, TDB.
GMGravitational 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.
mMass 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().

◆ MeanAnomalyFromTimeOfPerihelionPassage()

double pcl::OsculatingElements::MeanAnomalyFromTimeOfPerihelionPassage ( const TimePoint t,
double  GM = 0,
double  m = 0 
) const
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.

◆ MeanMotion()

double pcl::OsculatingElements::MeanMotion ( double  GM = 0,
double  m = 0 
) const

Mean orbital motion.

Parameters
GMGravitational 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.
mMass 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.

◆ OrbitalPeriod()

double pcl::OsculatingElements::OrbitalPeriod ( double  GM = 0,
double  m = 0 
) const

Orbital period of the elliptic orbit.

Parameters
GMGravitational 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.
mMass 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.

◆ PerifocalDistanceFromSemimajorAxis()

double pcl::OsculatingElements::PerifocalDistanceFromSemimajorAxis ( ) const
inline

The perifocal distance calculated from semimajor axis and eccentricity.

Definition at line 298 of file OsculatingElements.h.

◆ PerihelionDistanceFromSemimajorAxis()

double pcl::OsculatingElements::PerihelionDistanceFromSemimajorAxis ( ) const
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.

◆ SemimajorAxisFromPerifocalDistance()

double pcl::OsculatingElements::SemimajorAxisFromPerifocalDistance ( ) const
inline

The semimajor axis calculated from perifocal distance and eccentricity.

Definition at line 317 of file OsculatingElements.h.

References pcl::Abs().

◆ SemimajorAxisFromPerihelionDistance()

double pcl::OsculatingElements::SemimajorAxisFromPerihelionDistance ( ) const
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.

◆ TimeOfPerifocalPassageFromMeanAnomaly()

TimePoint pcl::OsculatingElements::TimeOfPerifocalPassageFromMeanAnomaly ( const TimePoint t,
double  GM = 0,
double  m = 0 
) const
inline

The time of perifocal passage calculated from the mean anomaly.

Parameters
tEpoch of osculating elements, TDB.
GMGravitational 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.
mMass 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.

◆ TimeOfPerihelionPassageFromMeanAnomaly()

TimePoint pcl::OsculatingElements::TimeOfPerihelionPassageFromMeanAnomaly ( const TimePoint t,
double  GM = 0,
double  m = 0 
) const
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.

◆ ToPerifocalPosition()

void pcl::OsculatingElements::ToPerifocalPosition ( double &  xp,
double &  yp,
const TimePoint t,
double  GM = 0,
double  m = 0 
) const
inline

Position components in the perifocal reference frame.

Parameters
[out]xp,ypPerifocal position vector components.
tEpoch of osculating elements (TDB).
GMGravitational 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.
mMass 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.

◆ ToPerifocalPositionAndVelocity()

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.

Parameters
[out]xp,ypPerifocal position vector components.
[out]xv,yvPerifocal velocity vector components.
tEpoch of osculating elements (TDB).
GMGravitational 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.
mMass 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.

◆ ToStateVectors()

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.

Parameters
[out]rPosition vector.
[out]vVelocity vector.
tEpoch of state vectors (TDB).
GMGravitational 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.
mMass of the body in central body mass units. Zero by default.
Note
The output r and v vectors will be referred to the same inertial rectangular coordinate system where orbital elements are defined. For example, in the case of a heliocentric orbit, the state vectors will give heliocentric ecliptic rectangular coordiantes.

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