SuperNOVAS C++ API v1.6
High-precision C/C++ astrometry library
Loading...
Searching...
No Matches
orbital.c File Reference

Function relating to the use of Keplerian orbital elements. More...

Functions

int novas_orbit_native_posvel (double jd_tdb, const novas_orbital *restrict orbit, double *restrict pos, double *restrict vel)
 Calculates a rectangular equatorial position and velocity vector for the given orbital elements for the specified time of observation, in the native coordinate system in which the orbital is defined (e.g.
int novas_orbit_posvel (double jd_tdb, const novas_orbital *restrict orbit, enum novas_accuracy accuracy, double *restrict pos, double *restrict vel)
 Calculates a rectangular equatorial position and velocity vector for the given orbital elements for the specified time of observation.
int novas_set_orbsys_pole (enum novas_reference_system type, double ra, double dec, novas_orbital_system *restrict sys)
 Sets the orientation of an orbital system using the RA and DEC coordinates of the pole of the Laplace (or else equatorial) plane relative to which the orbital elements are defined.

Detailed Description

Function relating to the use of Keplerian orbital elements.

For example, the IAU Minor Planet Center publishes up-to-date Keplerial orbital elements for all asteroids, comets, and Near-Earth Objects (NEOs) regularly. On short timescales these can provide accurate positions for such objects, that are as good as, or comparable to, the ephemeris data provided by NASA JPL.

For newly discovered objects, the MPC orbital elements may be the only source, or the most accurate source, of position information.

To use Keplerian orbital elements with SuperNOVAS, simply define the orbital parameters in a novas_orbital structure. By default heliocentric orbits are assumed, but you may also define orbitals around other bodies, such as a major planet, by defining the body (or barycenter) that is at the center of the orbital, and the orientation of the orbital system (w.r.t. the ecliptic or equator of date) with novas_set_orbsys_pole(). E.g.:

orb.a = ... // [AU] semi-major axis
... // populate the rest of the orbital parameters.
// Set Jupiter as the orbital center
// Set the orientation by defining the RA/Dec of the orbital pole, say in ICRS
novas_set_orbsys_pole(NOVAS_ICRS, ra_pole, dec_pole, &orb.system);
int novas_set_orbsys_pole(enum novas_reference_system type, double ra, double dec, novas_orbital_system *restrict sys)
Sets the orientation of an orbital system using the RA and DEC coordinates of the pole of the Laplace...
Definition orbital.c:414
@ NOVAS_ICRS
International Celestial Reference system.
Definition novas.h:876
#define NOVAS_ORBIT_INIT
Initializer for novas_orbital for heliocentric orbits using GCRS ecliptic parametrization.
Definition novas.h:1431
#define NOVAS_JUPITER_INIT
object initializer for the planet Jupiter
Definition novas.h:1538
enum novas_planet center
major planet or barycenter at the center of the orbit.
Definition novas.h:1346
Keplerian orbital elements for NOVAS_ORBITAL_OBJECT type.
Definition novas.h:1400
struct novas_orbital_system system
orbital reference system assumed for the parametrization
Definition novas.h:1401
double a
[AU] semi-major axis
Definition novas.h:1403

Function Documentation

◆ novas_orbit_native_posvel()

int novas_orbit_native_posvel ( double jd_tdb,
const novas_orbital *restrict orbit,
double *restrict pos,
double *restrict vel )

Calculates a rectangular equatorial position and velocity vector for the given orbital elements for the specified time of observation, in the native coordinate system in which the orbital is defined (e.g.

ecliptic for heliocentric orbitals).

REFERENCES:

  1. E.M. Standish and J.G. Williams 1992.
  2. https://ssd.jpl.nasa.gov/planets/approx_pos.html
  3. https://en.wikipedia.org/wiki/Orbital_elements
  4. https://orbitalofficial.com/
  5. https://downloads.rene-schwarz.com/download/M001-Keplerian_Orbit_Elements_to_Cartesian_State_Vectors.pdf
Parameters
jd_tdb[day] Barycentric Dynamic Time (TDB) based Julian date
orbitOrbital parameters
[out]pos[AU] Output ICRS equatorial position vector around orbital center, or NULL if not required.
[out]vel[AU/day] Output ICRS velocity vector rel. to orbital center, in the native system of the orbital, or NULL if not required. 0 if successful, or else -1 if the orbital parameters are NULL, or if the position and velocity output vectors are the same or the orbital system is ill defined (errno set to EINVAL), or if the calculation did not converge (errno set to ECANCELED).
Author
Attila Kovacs
Since
1.4
See also
novas_geom_posvel(), ephemeris(), make_orbital_object()

References TWOPI.

Referenced by novas_moon_phase(), and novas_orbit_posvel().

◆ novas_orbit_posvel()

int novas_orbit_posvel ( double jd_tdb,
const novas_orbital *restrict orbit,
enum novas_accuracy accuracy,
double *restrict pos,
double *restrict vel )

Calculates a rectangular equatorial position and velocity vector for the given orbital elements for the specified time of observation.

REFERENCES:

  1. E.M. Standish and J.G. Williams 1992.
  2. https://ssd.jpl.nasa.gov/planets/approx_pos.html
  3. https://en.wikipedia.org/wiki/Orbital_elements
  4. https://orbitalofficial.com/
  5. https://downloads.rene-schwarz.com/download/M001-Keplerian_Orbit_Elements_to_Cartesian_State_Vectors.pdf
Parameters
jd_tdb[day] Barycentric Dynamic Time (TDB) based Julian date
orbitOrbital parameters
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1).
[out]pos[AU] Output ICRS equatorial position vector around orbital center, or NULL if not required.
[out]vel[AU/day] Output ICRS equatorial velocity vector rel. to orbital center, or NULL if not required.
Returns
0 if successful, or else -1 if the orbital parameters are NULL, or if the position and velocity output vectors are the same or the orbital system is ill defined (errno set to EINVAL), or if the calculation did not converge (errno set to ECANCELED).
Author
Attila Kovacs
Since
1.2
See also
novas_geom_posvel(), ephemeris(), make_orbital_object()

References novas_orbit_native_posvel().

Referenced by ephemeris(), novas_approx_heliocentric(), supernovas::Orbital::position(), and supernovas::Orbital::velocity().

◆ novas_set_orbsys_pole()

int novas_set_orbsys_pole ( enum novas_reference_system type,
double ra,
double dec,
novas_orbital_system *restrict sys )

Sets the orientation of an orbital system using the RA and DEC coordinates of the pole of the Laplace (or else equatorial) plane relative to which the orbital elements are defined.

Orbital parameters of planetary satellites normally include the R.A. and declination of the pole of the local Laplace plane in which the Keplerian orbital elements are referenced.

The system will become referenced to the equatorial plane, the relative obliquity is set to (90° - dec), while the argument of the ascending node ('Omega') is set to (90° + ra).

NOTES:

  1. You should not expect much precision from the long-range orbital approximations for planetary satellites. For applications that require precision at any level, you should rely on appropriate ephemerides, or else on up-to-date short-term orbital elements.
Parameters
typeCoordinate reference system in which ra and dec are defined (e.g. NOVAS_GCRS).
ra[h] the R.A. of the pole of the orbital reference plane.
dec[deg] the declination of the pole of the oribtal reference plane.
[out]sysOrbital system
Returns
0 if successful, or else -1 (errno will be set to EINVAL) if the output sys pointer is NULL.
Author
Attila Kovacs
Since
1.2
See also
make_orbital_object()

References NOVAS_EQUATORIAL_PLANE, and NOVAS_TIRS.