![]() |
SuperNOVAS v1.5
The NOVAS C library, made better
|
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. | |
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.:
Once the orbital is defined, you can use it to define a generic object via make_orbital_object()
to enable the the same astrometric calculations as for ephemeris sources, e.g.:
Or, you can calculate geometric orbital positions and velocities directly via novas_orbit_posvel()
, while novas_orbit_native_posvel()
will do the same, but in the native coordinate system in which the orbital is defined.
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:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
orbit | Orbital 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). |
References TWOPI.
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:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
orbit | Orbital parameters | |
accuracy | NOVAS_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. |
References novas_orbit_native_posvel().
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:
type | Coordinate reference system in which ra and dec are defined (e.g. NOVAS_GCRS). | |
ra | [h] the R.A. of the pole of the oribtal reference plane. | |
dec | [deg] the declination of the pole of the oribtal reference plane. | |
[out] | sys | Orbital system |
sys
pointer is NULL.References NOVAS_EQUATORIAL_PLANE, and NOVAS_TIRS.