![]() |
SuperNOVAS C++ API v1.6
High-precision C/C++ astrometry library
|
Keplerian orbital elements, for example, for a comet using parameters published by the IAU Minor Planet Center. More...
#include <supernovas.h>
Public Member Functions | |
| Orbital (const OrbitalSystem &system, const Time &ref_time, const Coordinate &semi_major, const Angle &mean_anomaly, const Interval &period) | |
| Instantiates a new Keplerian orbital in the specified orbital system and the basic circular orbital parameters. | |
| Orbital (const OrbitalSystem &system, double jd_tdb, double semi_major_m, double mean_anomaly_rad, double period_s) | |
| Instantiates a new Keplerian orbital in the specified orbital system and the basic circular orbital parameters. | |
| const novas_orbital * | _novas_orbital () const |
| (for internal use) Returns the underlying NOVAS C data structure containing the orbital parameters. | |
| Interval | apsis_period () const |
| Returns the rotation period of the apsis location in the orbital system, in which the orbital was defined (positive for counter-clockwise rotation, or negative for clockwise rotation, when viewed from the orbital system's pole.). | |
| Orbital & | apsis_period (const Interval &period) |
| Sets the apsis rotation period (positive for counter-clockwise rotation when viewed from the orbital system's pole). | |
| Orbital & | apsis_period (double seconds) |
| Sets the apsis rotation period (positive for counter-clockwise rotation when viewed from the orbital system's pole). | |
| double | apsis_rate () const |
| Returns the angular rate at which the apsis rotates in the orbital system (positive for counter-clockwise rotation, or negative for clockwise rotation, when viewed from the orbital system's pole.). | |
| Orbital & | apsis_rate (double rad_per_sec) |
| Sets the apsis rotation rate (positive for counter-clockwise rotation when viewed from the orbital system's pole). | |
| Angle | ascending_node () const |
| Returns the longitude of the ascending node of this orbit in the orbital system, in which the orbital was defined. | |
| double | eccentricity () const |
| Returns the eccentricity of this orbital. | |
| Orbital & | eccentricity (double e, const Angle &periapsis_angle) |
| Sets parameters for an elliptical orbit. | |
| Orbital & | eccentricity (double e, double periapsis_rad) |
| Sets parameters for an elliptical orbit. | |
| Angle | inclination () const |
| Returns the inclination angle of this orbit, relative to the orbital system's plane. | |
| Orbital & | inclination (const Angle &angle, const Angle &ascending_node_angle) |
| Sets parameters for an orbit that is inclined relative to the orbital system's native plane. | |
| Orbital & | inclination (double angle_rad, double ascending_node_rad) |
| Sets parameters for an orbit that is inclined relative to the orbital system's native plane. | |
| double | mean_motion () const |
| Returns the mean motion (circular angular velocity) of the object in this orbit. | |
| Interval | node_period () const |
| Returns the rotation period (due to precession) of the orbit's ascending node in the orbital system, in which the orbital was defined (positive for counter-clockwise rotation, or negative for clockwise rotation, when viewed from the orbital system's pole.). | |
| Orbital & | node_period (const Interval &period) |
| Sets the node precession period (positive for counter-clockwise rotation when viewed from the orbital system's pole). | |
| Orbital & | node_period (double seconds) |
| Sets the node precession period (positive for counter-clockwise rotation when viewed from the orbital system's pole). | |
| double | node_rate () const |
| Returns the angular rate at which the ascending node of the orbit rotates in the orbital system (positive for counter-clockwise rotation, or negative for clockwise rotation, when viewed from the orbital system's pole.). | |
| Orbital & | node_rate (double rad_per_sec) |
| Sets the node precession rate (positive for counter-clockwise rotation when viewed from the orbital system's pole). | |
| Angle | periapsis () const |
| Returns the periapsis angle of this orbit, in the orbital system, in which the orbit was defined. | |
| Interval | period () const |
| Returns the period of the object on this orbit. | |
| Spherical | pole () const |
| Returns the spherical coordinates of the orbit's pole in the orbital system, in which the orbital was defined. | |
| Orbital & | pole (const Angle &longitude, const Angle &latitude) |
| Sets the orbit's pole, in the orbital system in which the orbit is defined. | |
| Orbital & | pole (const Spherical &coords) |
| Sets the orbit's pole, in the orbital system in which the orbit is defined. | |
| Orbital & | pole (double longitude_rad, double latitude_rad) |
| Sets the orbit's pole, in the orbital system in which the orbit is defined. | |
| Position | position (const Time &time, enum novas_accuracy accuracy=NOVAS_FULL_ACCURACY) const |
| Calculates a rectangular equatorial position vector for this Keplerian orbital for the specified time of observation. | |
| double | reference_jd_tdb () const |
| Returns the reference time, as a Barycentric Dynamical Time (TDB) based Julian date. | |
| Angle | reference_mean_anomaly () const |
| Returns the mean anomaly (or longitude for circular orbits) of the object at the reference time, in the orbital system in which the orbit was defined. | |
| Coordinate | semi_major_axis () const |
| Returns the semi-major axis (that is the radius for circular orbits) of this orbital. | |
| OrbitalSystem | system () const |
| Returns a new instance of the orbital system in which this orbit is defined. | |
| OrbitalSource | to_source (const std::string &name) const |
| Returns a new orbital source from this orbital and with the specified name and number designations. | |
| std::string | to_string () const |
| Returns a basic human-readable description of this orbital, with just the major parameters. | |
| Velocity | velocity (const Time &time, enum novas_accuracy accuracy=NOVAS_FULL_ACCURACY) const |
| Calculates a rectangular equatorial velocity vector for this Keplerian orbital for the specified time of observation. | |
| Public Member Functions inherited from supernovas::Validating | |
| bool | is_valid () const |
| Returns the previously set 'valid' stae of the implementing instance. | |
| operator bool () const | |
| Objects that implement Validating can be used in conditionals directly, without explicitly calling is_valid(). | |
Static Public Member Functions | |
| static Orbital | from_mean_motion (const OrbitalSystem &system, const Time &ref_time, const Coordinate &semi_major, const Angle &mean_anomaly, double rad_per_s) |
| Return a new instance of a new Keplerian orbital in the specified orbital system and the basic circular orbital parameters, with mean motion used instead of a orbital period. | |
| static Orbital | from_mean_motion (const OrbitalSystem &system, double jd_tdb, double semi_majpr_m, double mean_anomaly_rad, double rad_per_s) |
| Return a new instance of a new Keplerian orbital in the specified orbital system and the basic circular orbital parameters, with mean motion used instead of a orbital period. | |
| static Orbital | from_novas_orbit (const novas_orbital *orbit) |
| (for internal use) Returns a new instance of a Keplerian orbital, using a copy of a NOVAS C orbital data structure, or std::nullopt if the argument is NULL. | |
| static Orbital | moon_mean_orbit_at (const Time &time) |
| Returns the mean Keplerian orbital for the Moon relative to the geocenter for the specified epoch of observation. | |
| static Orbital | moon_orbit_at (const Time &time) |
| Returns an approximation of the current Keplerian orbital of the Moon relative to the geocenter for the specified epoch of observation. | |
Additional Inherited Members | |
| Protected Member Functions inherited from supernovas::Validating | |
| Validating () | |
| dummy constructor; | |
| Protected Attributes inherited from supernovas::Validating | |
| bool | _valid = false |
| the state variable. | |
Keplerian orbital elements, for example, for a comet using parameters published by the IAU Minor Planet Center.
While Keplerian orbitals cannot provide accurate positions or velocities for Solar-system bodies over the long term (for that you need ephemeris data), they can be sufficiently accurate on the short term. And, in case of recently discovered objects, such as Near-Earth Objects (NEOs), orbital elements may be the only source of up-to-date positional data.
NOTES:
| supernovas::Orbital::Orbital | ( | const OrbitalSystem & | system, |
| double | jd_tdb, | ||
| double | semi_major_m, | ||
| double | mean_anomaly_rad, | ||
| double | period_s ) |
Instantiates a new Keplerian orbital in the specified orbital system and the basic circular orbital parameters.
You can further specify the parameters for elliptical orbits using a builder pattern after instantiation.
| system | the orbital system in which the orbit is defined. |
| jd_tdb | [day] reference date of the orbital parameters as a Barycentric Dynamical Time (TDB) based Julian date |
| semi_major_m | [m] semi-major axis (circular radius) of the orbit |
| mean_anomaly_rad | [rad] Mean anomaly (circular longitude) of the object at the reference time, in the orbital system. |
| period_s | [s] orbital period. |
References supernovas::Unit::au, supernovas::Unit::day, supernovas::Unit::deg, and system().
| supernovas::Orbital::Orbital | ( | const OrbitalSystem & | system, |
| const Time & | ref_time, | ||
| const Coordinate & | semi_major, | ||
| const Angle & | mean_anomaly, | ||
| const Interval & | period ) |
Instantiates a new Keplerian orbital in the specified orbital system and the basic circular orbital parameters.
You can further specify the parameters for elliptical orbits using a builder pattern after instantiation.
| system | the orbital system in which the orbit is defined. |
| ref_time | reference time of the orbital parameters. |
| semi_major | semi-major axis (circular radius) of the orbit |
| mean_anomaly | Mean anomaly (circular longitude) of the object at the reference time, in the orbital system. |
| period | orbital period. |
| const novas_orbital * supernovas::Orbital::_novas_orbital | ( | ) | const |
(for internal use) Returns the underlying NOVAS C data structure containing the orbital parameters.
Referenced by supernovas::OrbitalSource::OrbitalSource().
| Interval supernovas::Orbital::apsis_period | ( | ) | const |
Returns the rotation period of the apsis location in the orbital system, in which the orbital was defined (positive for counter-clockwise rotation, or negative for clockwise rotation, when viewed from the orbital system's pole.).
References supernovas::Unit::day, and supernovas::Validating::is_valid().
Referenced by apsis_period(), apsis_rate(), and apsis_rate().
| Orbital & supernovas::Orbital::apsis_period | ( | const Interval & | period | ) |
Sets the apsis rotation period (positive for counter-clockwise rotation when viewed from the orbital system's pole).
| period | time it takes for a full rotation of the apsis in the orbital system. It may be negative for clockwise (retrograde) rotation seen from the orbital system's pole. |
References apsis_period(), and period().
| Orbital & supernovas::Orbital::apsis_period | ( | double | seconds | ) |
Sets the apsis rotation period (positive for counter-clockwise rotation when viewed from the orbital system's pole).
| seconds | [s] counter-clockwise rotation period of the apsis. It may be negative for clockwise (retrograde) rotation seen from the orbital system's pole. |
References supernovas::Validating::_valid, and supernovas::Unit::day.
| double supernovas::Orbital::apsis_rate | ( | ) | const |
Returns the angular rate at which the apsis rotates in the orbital system (positive for counter-clockwise rotation, or negative for clockwise rotation, when viewed from the orbital system's pole.).
References apsis_period(), and supernovas::Constant::two_pi.
| Orbital & supernovas::Orbital::apsis_rate | ( | double | rad_per_sec | ) |
Sets the apsis rotation rate (positive for counter-clockwise rotation when viewed from the orbital system's pole).
| rad_per_sec | [rad/s] counter-clockwise rotation rate of the apsis. It may be negative for clockwise (retrograde) rotation seen from the orbital system's pole. |
References apsis_period(), and supernovas::Constant::two_pi.
| Angle supernovas::Orbital::ascending_node | ( | ) | const |
Returns the longitude of the ascending node of this orbit in the orbital system, in which the orbital was defined.
References supernovas::Unit::deg, and supernovas::Validating::is_valid().
| double supernovas::Orbital::eccentricity | ( | ) | const |
Returns the eccentricity of this orbital.
Referenced by eccentricity(), and to_string().
| Orbital & supernovas::Orbital::eccentricity | ( | double | e, |
| const Angle & | periapsis_angle ) |
Sets parameters for an elliptical orbit.
| e | eccenticity value (dimensionless). |
| periapsis_angle | longitude of the apsis (the point at which the elliptial orbit is closest to the center), in the orbital system, in which the orbit is defined. |
References eccentricity(), and supernovas::Angle::rad().
| Orbital & supernovas::Orbital::eccentricity | ( | double | e, |
| double | periapsis_rad ) |
Sets parameters for an elliptical orbit.
| e | eccenticity value (dimensionless). |
| periapsis_rad | [rad] longitude of the apsis (the point at which the elliptial orbit is closest to the center), in the orbital system, in which the orbit is defined. |
References supernovas::Validating::_valid, and supernovas::Unit::deg.
|
static |
Return a new instance of a new Keplerian orbital in the specified orbital system and the basic circular orbital parameters, with mean motion used instead of a orbital period.
You can further specify the parameters for elliptical orbits using a builder pattern after instantiation.
| system | the orbital system in which the orbit is defined. |
| ref_time | reference time of the orbital parameters. |
| semi_major | semi-major axis (circular radius; a) of the orbit |
| mean_anomaly | Mean anomaly (circular longitude; M0) of the object at the reference time, in the orbital system. |
| rad_per_sec | [rad/s] mean motion (circular angular velocity) on orbit. |
References from_mean_motion(), supernovas::Time::jd(), supernovas::Coordinate::m(), NOVAS_TDB, supernovas::Angle::rad(), and system().
|
static |
Return a new instance of a new Keplerian orbital in the specified orbital system and the basic circular orbital parameters, with mean motion used instead of a orbital period.
You can further specify the parameters for elliptical orbits using a builder pattern after instantiation.
| system | the orbital system in which the orbit is defined. |
| jd_tdb | [day] reference date of the orbital parameters as a Barycentric Dynamical Time (TDB) based Julian date |
| semi_major_m | [m] semi-major axis (circular radius; a) of the orbit |
| mean_anomaly_rad | [rad] Mean anomaly (circular longitude; M0) of the object at the reference time, in the orbital system. |
| rad_per_sec | [rad/s] mean motion (circular angular velocity) on orbit. |
References supernovas::Validating::is_valid(), system(), and supernovas::Constant::two_pi.
Referenced by from_mean_motion().
|
static |
(for internal use) Returns a new instance of a Keplerian orbital, using a copy of a NOVAS C orbital data structure, or std::nullopt if the argument is NULL.
It's best practice to call is_valid() after to check that the supplied parameters do in fact define a valid orbital system.
| orbit | The NOVAS C orbital data structure (copied) |
References supernovas::Validating::is_valid().
Referenced by moon_mean_orbit_at(), moon_orbit_at(), supernovas::Planet::orbit(), and supernovas::OrbitalSource::orbital().
| Angle supernovas::Orbital::inclination | ( | ) | const |
Returns the inclination angle of this orbit, relative to the orbital system's plane.
References supernovas::Unit::deg, and supernovas::Validating::is_valid().
Referenced by inclination(), and pole().
| Orbital & supernovas::Orbital::inclination | ( | const Angle & | angle, |
| const Angle & | ascending_node_angle ) |
Sets parameters for an orbit that is inclined relative to the orbital system's native plane.
| angle | inclination angle |
| ascending_node_angle | longitude of the ascending node in the orbital system, in which the orbit is defined. |
References inclination(), and supernovas::Angle::rad().
| Orbital & supernovas::Orbital::inclination | ( | double | angle_rad, |
| double | ascending_node_rad ) |
Sets parameters for an orbit that is inclined relative to the orbital system's native plane.
| angle_rad | [rad] inclination angle |
| ascending_node_rad | [rad] longitude of the ascending node in the orbital system, in which the orbit is defined. |
References supernovas::Validating::_valid, and supernovas::Unit::deg.
| double supernovas::Orbital::mean_motion | ( | ) | const |
Returns the mean motion (circular angular velocity) of the object in this orbit.
References supernovas::Unit::day, and supernovas::Unit::deg.
Referenced by period().
|
static |
Returns the mean Keplerian orbital for the Moon relative to the geocenter for the specified epoch of observation.
It is based on the secular parameters of the ELP2000-85 model, but does not include the harmonic series and the perturbation terms. As such it has accuracy at the few degrees level only, however it is 'valid' for long-term projections (i.e. for years around the orbit's reference epoch) at that coarse level.
For the short-term , Orbital::moon_orbit_at() can provide somewhat more accurate predictions for up to a day or so around the reference epoch of the orbit.
REFERENCES:
| Time of observation and astronomical timescales | Astrometric time for which to calculate the secular orbital parameters of the moon. |
References from_novas_orbit(), supernovas::Validating::is_valid(), supernovas::Time::jd(), novas_make_moon_mean_orbit(), and NOVAS_TDB.
|
static |
Returns an approximation of the current Keplerian orbital of the Moon relative to the geocenter for the specified epoch of observation.
The orbit includes the most dominant Solar perturbation terms to produce results with an accuracy at the ~10 arcmin level near (+- 0.5 days) the reference time argument of the orbit. The perturbed orbit is based on the ELP/MPP02 model.
While, the ELP/MPP02 model itself can be highly precise, the Moon's orbit is strongly non-Keplerian, and so any attempt to describe it in purely Keplerian terms is inherently flawed, which is the reason for the generally poor accuracy of this model.
REFERENCES:
| Time of observation and astronomical timescales | Astrometric time for which to calculate the secular orbital parameters of the moon. |
References from_novas_orbit(), supernovas::Validating::is_valid(), supernovas::Time::jd(), novas_make_moon_orbit(), and NOVAS_TDB.
| Interval supernovas::Orbital::node_period | ( | ) | const |
Returns the rotation period (due to precession) of the orbit's ascending node in the orbital system, in which the orbital was defined (positive for counter-clockwise rotation, or negative for clockwise rotation, when viewed from the orbital system's pole.).
References supernovas::Unit::day, and supernovas::Validating::is_valid().
Referenced by node_period(), node_rate(), and node_rate().
| Orbital & supernovas::Orbital::node_period | ( | const Interval & | period | ) |
Sets the node precession period (positive for counter-clockwise rotation when viewed from the orbital system's pole).
| period | counter-clockwise precession period of the node. It may be negative for clockwise (retrograde) rotation seen from the orbital system's pole. |
References node_period(), and period().
| Orbital & supernovas::Orbital::node_period | ( | double | seconds | ) |
Sets the node precession period (positive for counter-clockwise rotation when viewed from the orbital system's pole).
| seconds | [s] counter-clockwise precession period of the node. It may be negative for clockwise (retrograde) rotation seen from the orbital system's pole. |
References supernovas::Validating::_valid, and supernovas::Unit::day.
| double supernovas::Orbital::node_rate | ( | ) | const |
Returns the angular rate at which the ascending node of the orbit rotates in the orbital system (positive for counter-clockwise rotation, or negative for clockwise rotation, when viewed from the orbital system's pole.).
References node_period(), and supernovas::Constant::two_pi.
| Orbital & supernovas::Orbital::node_rate | ( | double | rad_per_sec | ) |
Sets the node precession rate (positive for counter-clockwise rotation when viewed from the orbital system's pole).
| rad_per_sec | [rad/sec] counter-clockwise precession rate of the node. It may be negative for clockwise (retrograde) rotation seen from the orbital system's pole. |
References node_period(), and supernovas::Constant::two_pi.
| Angle supernovas::Orbital::periapsis | ( | ) | const |
Returns the periapsis angle of this orbit, in the orbital system, in which the orbit was defined.
References supernovas::Unit::deg, and supernovas::Validating::is_valid().
| Interval supernovas::Orbital::period | ( | ) | const |
Returns the period of the object on this orbit.
References supernovas::Validating::is_valid(), mean_motion(), and supernovas::Constant::two_pi.
Referenced by Orbital(), apsis_period(), node_period(), and to_string().
| Spherical supernovas::Orbital::pole | ( | ) | const |
Returns the spherical coordinates of the orbit's pole in the orbital system, in which the orbital was defined.
References supernovas::Unit::deg, supernovas::Constant::half_pi, and supernovas::Validating::is_valid().
Sets the orbit's pole, in the orbital system in which the orbit is defined.
| longitude | longitude of orbit's pole in the orbital system. |
| latitude | latitude of the orbit's pole in the orbital system. |
References pole(), and supernovas::Angle::rad().
| Orbital & supernovas::Orbital::pole | ( | const Spherical & | coords | ) |
Sets the orbit's pole, in the orbital system in which the orbit is defined.
| coords | location of the pole in the orbital system. |
References supernovas::Spherical::latitude(), supernovas::Spherical::longitude(), and pole().
| Orbital & supernovas::Orbital::pole | ( | double | longitude_rad, |
| double | latitude_rad ) |
Sets the orbit's pole, in the orbital system in which the orbit is defined.
| longitude_rad | [rad] longitude of orbit's pole in the orbital system. |
| latitude_rad | [rad] latitude of the orbit's pole in the orbital system. |
References supernovas::Constant::half_pi, and inclination().
| Position supernovas::Orbital::position | ( | const Time & | time, |
| enum novas_accuracy | accuracy = NOVAS_FULL_ACCURACY ) const |
Calculates a rectangular equatorial position vector for this Keplerian orbital for the specified time of observation.
REFERENCES:
| Time of observation and astronomical timescales | Astrometric time of observation |
| accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1). |
References supernovas::Unit::au, supernovas::Time::jd(), novas_orbit_posvel(), and supernovas::Position::undefined().
| double supernovas::Orbital::reference_jd_tdb | ( | ) | const |
Returns the reference time, as a Barycentric Dynamical Time (TDB) based Julian date.
| Angle supernovas::Orbital::reference_mean_anomaly | ( | ) | const |
Returns the mean anomaly (or longitude for circular orbits) of the object at the reference time, in the orbital system in which the orbit was defined.
References supernovas::Unit::deg, and supernovas::Validating::is_valid().
| Coordinate supernovas::Orbital::semi_major_axis | ( | ) | const |
Returns the semi-major axis (that is the radius for circular orbits) of this orbital.
References supernovas::Unit::au, and supernovas::Validating::is_valid().
Referenced by to_string().
| OrbitalSystem supernovas::Orbital::system | ( | ) | const |
Returns a new instance of the orbital system in which this orbit is defined.
References supernovas::OrbitalSystem::from_novas_orbital_system(), and supernovas::Validating::is_valid().
Referenced by Orbital(), Orbital(), from_mean_motion(), from_mean_motion(), and to_string().
| std::string supernovas::Orbital::to_string | ( | ) | const |
Returns a basic human-readable description of this orbital, with just the major parameters.
References eccentricity(), period(), semi_major_axis(), system(), supernovas::Coordinate::to_string(), supernovas::Interval::to_string(), and supernovas::OrbitalSystem::to_string().
| Velocity supernovas::Orbital::velocity | ( | const Time & | time, |
| enum novas_accuracy | accuracy = NOVAS_FULL_ACCURACY ) const |
Calculates a rectangular equatorial velocity vector for this Keplerian orbital for the specified time of observation.
REFERENCES:
| Time of observation and astronomical timescales | Astrometric time of observation |
| accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1). |
References supernovas::Unit::au, supernovas::Unit::day, supernovas::Time::jd(), novas_orbit_posvel(), and supernovas::Velocity::undefined().