SuperNOVAS v1.5
The NOVAS C library, made better
Loading...
Searching...
No Matches
Geometric equatorial positions and velocities

Functions

short astro_planet (double jd_tt, const object *restrict ss_body, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec, double *restrict dis)
 reference_system: ICRS / GCRS
observer location: geocenter
position_type: geometric

 
short astro_star (double jd_tt, const cat_entry *restrict star, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec)
 reference_system: ICRS / GCRS
observer location: geocenter
position_type: geometric

 
short ephemeris (const double *restrict jd_tdb, const object *restrict body, enum novas_origin origin, enum novas_accuracy accuracy, double *restrict pos, double *restrict vel)
 Retrieves the position and velocity of a solar system body using the currently configured plugins that provide them.
 
short geo_posvel (double jd_tt, double ut1_to_tt, enum novas_accuracy accuracy, const observer *restrict obs, double *restrict pos, double *restrict vel)
 Computes the geocentric GCRS position and velocity of an observer.
 
short light_time (double jd_tdb, const object *restrict body, const double *pos_obs, double tlight0, enum novas_accuracy accuracy, double *pos_src_obs, double *restrict tlight)
 Computes the geocentric position of a solar system body, as antedated for light-time.
 
int light_time2 (double jd_tdb, const object *restrict body, const double *restrict pos_obs, double tlight0, enum novas_accuracy accuracy, double *p_src_obs, double *restrict v_ssb, double *restrict tlight)
 Computes the geocentric position and velocity of a solar system body, as antedated for light-time.
 
int novas_approx_heliocentric (enum novas_planet id, double jd_tdb, double *restrict pos, double *restrict vel)
 Returns the approximate geometric heliocentric orbital positions and velocities for the major planets, Sun, or Earth-Moon Barycenter (EMB).
 
int novas_geom_posvel (const object *restrict source, const novas_frame *restrict frame, enum novas_reference_system sys, double *restrict pos, double *restrict vel)
 Calculates the geometric position and velocity vectors, relative to the observer, for a source in the given observing frame, in the specified coordinate system of choice.
 
int novas_geom_to_app (const novas_frame *restrict frame, const double *restrict pos, enum novas_reference_system sys, sky_pos *restrict out)
 Converts an geometric position in ICRS to an apparent position on sky, by applying appropriate corrections for aberration and gravitational deflection for the observer's frame.
 
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_transform_vector (const double *in, const novas_transform *restrict transform, double *out)
 Transforms a position or velocity 3-vector from one coordinate reference system to another.
 
int obs_posvel (double jd_tdb, double ut1_to_tt, enum novas_accuracy accuracy, const observer *restrict obs, const double *restrict geo_pos, const double *restrict geo_vel, double *restrict pos, double *restrict vel)
 Calculates the ICRS position and velocity of the observer relative to the Solar System Barycenter (SSB).
 

Detailed Description

Function Documentation

◆ astro_planet()

short astro_planet ( double jd_tt,
const object *restrict ss_body,
enum novas_accuracy accuracy,
double *restrict ra,
double *restrict dec,
double *restrict dis )

reference_system: ICRS / GCRS
observer location: geocenter
position_type: geometric

Computes the astrometric place of a solar system body, as would be seen by a non-moving observer located at the current position of the geocenter, referenced to the ICRS without light deflection or aberration. This is the same as calling place_icrs() for the body, except the different set of return values used.

REFERENCES:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Explanatory Supplement to the Astronomical Almanac (1992),Chapter 3.
Parameters
jd_tt[day] Terretrial Time (TT) based Julian date.
ss_bodyPointer to structure containing the body designation for the solar system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]ra[h] Astrometric (geometric) right ascension for a fictitous geocentric observer, referred to the ICRS, without light deflection or aberration. (It may be NULL if not required)
[out]dec[deg] Astrometric (geometric) declination, for a fictitous geocentric observer, referred to the ICRS, without light deflection or aberration. (It may be NULL if not required)
[out]dis[AU] Apparent distance from Earth to the body at 'jd_tt' (it may be NULL if not needed).
Returns
0 if successful, or -1 if the object is NULL, or else 1 if the value of 'type' in structure 'ss_body' is invalid, or 10 + the error from place().
See also
place_icrs(), app_planet(), local_planet(), topo_planet(), virtual_planet(), astro_star()
novas_geom_posvel(), make_observer_at_geocenter(), NOVAS_GCRS

References NOVAS_ICRS, and radec_planet().

◆ astro_star()

short astro_star ( double jd_tt,
const cat_entry *restrict star,
enum novas_accuracy accuracy,
double *restrict ra,
double *restrict dec )

reference_system: ICRS / GCRS
observer location: geocenter
position_type: geometric

Computes the astrometric place of a star, as would be seen by a non-moving observer at the current location of the geocenter, referred to the ICRS without light deflection or aberration, at date 'jd_tt', given its catalog mean place, proper motion, parallax, and radial velocity.

Notwithstanding the different set of return values, this is the same as calling place_star() with a NULL observer location and NOVAS_ICRS as the system, or place_icrs() for an object that specifies the star.

REFERENCES:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Explanatory Supplement to the Astronomical Almanac (1992), Chapter 3.
Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date.
starPointer to catalog entry structure containing catalog data for the object in the ICRS.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]ra[h] Astrometric (geometric) right ascension for a fictitous geocentric observer, referred to the ICRS, without light deflection or aberration. (It may be NULL if not required)
[out]dec[deg] Astrometric (geometric) declination for a fictitous geocentric observer, referred to the ICRS, without light deflection or aberration. (It may be NULL if not required)
Returns
0 if successful, or -1 if a required pointer argument is NULL, or 20 + the error from place().
See also
place_star(), place_icrs(), app_star(), local_star(), topo_star(), virtual_star(), astro_planet()
novas_geom_posvel(), make_observer_at_geocenter(), NOVAS_ICRS

References NOVAS_ICRS, and radec_star().

◆ ephemeris()

short ephemeris ( const double *restrict jd_tdb,
const object *restrict body,
enum novas_origin origin,
enum novas_accuracy accuracy,
double *restrict pos,
double *restrict vel )

Retrieves the position and velocity of a solar system body using the currently configured plugins that provide them.

It is recommended that the input structure 'cel_obj' be created using make_object()

Parameters
jd_tdb[day] Barycentric Dynamic Time (TDB) based Julian date
bodyPointer to structure containing the designation of the body of interest
originNOVAS_BARYCENTER (0) or NOVAS_HELIOCENTER (1)
accuracyNOCAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]pos[AU] Pointer to structure containing the designation of the body of interest
[out]vel[AU/day] Velocity vector of the body at 'jd_tdb'; equatorial rectangular coordinates in AU/day referred to the ICRS.
Returns
0 if successful, -1 if the 'jd_tdb' or input object argument is NULL, or else 1 if 'origin' is invalid, 2 if cel_obj->type is invalid, 10 + the error code from the currently configured novas_planet_provider_hp call, or 20 + the error code from readeph().
See also
set_planet_provider(), set_planet_provider_hp(), set_ephem_provider()
make_planet(), make_ephem_object()

References ephemeris(), get_ephem_provider(), make_planet(), NOVAS_BARYCENTER, NOVAS_EPHEM_OBJECT, NOVAS_FULL_ACCURACY, NOVAS_HELIOCENTER, novas_orbit_posvel(), NOVAS_ORBITAL_OBJECT, NOVAS_ORIGIN_TYPES, NOVAS_PLANET, NOVAS_SSB, NOVAS_SUN, and readeph().

◆ geo_posvel()

short geo_posvel ( double jd_tt,
double ut1_to_tt,
enum novas_accuracy accuracy,
const observer *restrict obs,
double *restrict pos,
double *restrict vel )

Computes the geocentric GCRS position and velocity of an observer.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date.
ut1_to_tt[s] TT - UT1 time difference in seconds
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
obsITRF / GRS80 geodetic observer location.
[out]pos[AU] Position 3-vector of observer, with respect to origin at geocenter, referred to GCRS axes, components in AU. (It may be NULL if not required.)
[out]vel[AU/day] Velocity 3-vector of observer, with respect to origin at geocenter, referred to GCRS axes, components in AU/day. (It must be distinct from the pos output vector, and may be NULL if not required)
Returns
0 if successful, -1 if the 'obs' is NULL or the two output vectors are the same, or else 1 if 'accuracy' is invalid, or 2 if 'obserrver->where' is invalid.
See also
novas_make_frame(), make_observer(), get_ut1_to_tt()

References ephemeris(), geo_posvel(), NOVAS_AIRBORNE_OBSERVER, NOVAS_BARYCENTER, NOVAS_EARTH_INIT, NOVAS_FULL_ACCURACY, novas_gast(), NOVAS_OBSERVER_AT_GEOCENTER, NOVAS_OBSERVER_IN_EARTH_ORBIT, NOVAS_OBSERVER_ON_EARTH, NOVAS_REDUCED_ACCURACY, NOVAS_SOLAR_SYSTEM_OBSERVER, terra(), tod_to_gcrs(), tt2tdb(), and observer::where.

◆ light_time()

short light_time ( double jd_tdb,
const object *restrict body,
const double * pos_obs,
double tlight0,
enum novas_accuracy accuracy,
double * pos_src_obs,
double *restrict tlight )

Computes the geocentric position of a solar system body, as antedated for light-time.

Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian date
bodyPointer to structure containing the designation for the solar system body
pos_obs[AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes.
tlight0[day] First approximation to light-time (can be set to 0.0 if not readily available – it will then be computed as needed).
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]pos_src_obs[AU] Position 3-vector of body, with respect to origin at observer (or the geocenter), referred to ICRS axes. It can be the same vector as either of the inputs.
[out]tlight[day] Calculated light time
Returns
0 if successful, -1 if any of the poiinter arguments is NULL, 1 if the algorithm failed to converge after 10 iterations, or 10 + the error from ephemeris().
See also
light_time2(), novas_make_frame()

References light_time2().

◆ light_time2()

int light_time2 ( double jd_tdb,
const object *restrict body,
const double *restrict pos_obs,
double tlight0,
enum novas_accuracy accuracy,
double * p_src_obs,
double *restrict v_ssb,
double *restrict tlight )

Computes the geocentric position and velocity of a solar system body, as antedated for light-time.

It is effectively the same as the original NOVAS C light_time(), except that this returns the antedated source velocity vector also.

NOTES:

  1. This function is called by novas_geom_posvel() or place() to calculate observed positions, radial velocity, and distance for the time when the observed light originated from the source.
Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian date
bodyPointer to structure containing the designation for the solar system body
pos_obs[AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes.
tlight0[day] First approximation to light-time (can be set to 0.0 if not readily avasilable – so it will be calculated as needed).
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]p_src_obs[AU] Position 3-vector of body, relative to observer, referred to ICRS axes, components in AU.
[out]v_ssb[AU/day] Velocity 3-vector of body, with respect to the Solar-system barycenter, referred to ICRS axes.
[out]tlight[day] Calculated light time, or NAN when returning with an error code.
Returns
0 if successful, -1 if any of the pointer arguments is NULL or if the output vectors are the same or if they are the same as pos_obs, 1 if the algorithm failed to converge after 10 iterations, or 10 + the error from ephemeris().
See also
light_time()
novas_sky_pos()
Since
1.0
Author
Attila Kovacs

References bary2obs(), ephemeris(), NOVAS_BARYCENTER, NOVAS_FULL_ACCURACY, and novas_inv_max_iter.

◆ novas_approx_heliocentric()

int novas_approx_heliocentric ( enum novas_planet id,
double jd_tdb,
double *restrict pos,
double *restrict vel )

Returns the approximate geometric heliocentric orbital positions and velocities for the major planets, Sun, or Earth-Moon Barycenter (EMB).

The returned positions and velocities are not suitable for precise calculations. However, they may be used to provide rough guidance about the approximate locations of the planets, e.g. for determining approximate rise or set times or the approximate azimuth / elevation angles from an observing site.

The orbitals can provide planet positions to arcmin-level precision for the rocky inner planets, and to a fraction of a degree precision for the gas and ice giants and Pluto. The accuracies for Uranus, Neptune, and Pluto are significantly improved (to the arcmin level) if used in the time range of 1800 AD to 2050 AD. For a more detailed summary of the typical accuracies, see either of the top two references below.

For accurate positions, you should use planetary ephemerides (such as the JPL ephemerides via the CALCEPH or CSPICE plugins) and novas_geom_posvel() instead.

While this function is generally similar to creating an orbital object with an orbit initialized with novas_make_planet_orbit() or novas_make_moon_orbit(), and then calling novas_geom_posvel(), there are a few important differences:

  1. This function returns geometric positions referenced to the Sun (i.e., heliocentric), whereas novas_geom_posvel() references the calculated positions to the Solar-system Barycenter (SSB).
  2. This function calculates Earth and Moon positions about the Keplerian orbital position of the Earth-Moon Barycenter (EMB). In constrast, novas_make_planet_orbit() does not provide orbitals for the Earth directly, and make_moot_orbit() references the Moon's orbital to the Earth position returned by the currently configured planet calculator function (see set_planet_provider()).

NOTES:

  1. The Sun's position w.r.t. the Solar-system Barycenter is calculated using earth_sun_calc(). All other orbitals are also referenced to the Sun's position calculated that way.

REFERENCES:

  1. E.M. Standish and J.G. Williams 1992.
  2. https://ssd.jpl.nasa.gov/planets/approx_pos.html
  3. Chapront, J. et al., 2002, A&A 387, 700–709
  4. Chapront-Touze, M, and Chapront, J. 1983, Astronomy and Astrophysics (ISSN 0004-6361), vol. 124, no. 1, July 1983, p. 50-62.
Parameters
idNOVAS major planet ID. All major planets, plus the Sun, Moon, Earth-Moon Barycenter (EMB), and Pluto system Barycenter are supported. (For Pluto, the Pluto System Barycenter value are returned.)
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian Date. Dates between 3000 BC and 3000 AD are supported. For dates between 1800 AD and 2050 AD the returned positions are somewhat more accurate.
[out]pos[AU] Output Heliocentric ICRS position vector, or NULL if not required.
[out]vel[AU/day] Output Heliocentric ICRS velocity vector, or NULL if not required.
Returns
0 if successful, or if the JD date is outside of the range with valid parameters, or else -1 if the planet ID is not supported or if both output vectors are NULL. In case of errors errno will be set to indicate the type of error.
Since
1.4
Author
Attila Kovacs
See also
novas_approx_sky_pos(), earth_sun_calc(), novas_geom_posvel()
novas_use_calceph(), novas_use_cspice()

References NOVAS_EARTH, NOVAS_EMB, novas_make_moon_orbit(), novas_make_planet_orbit(), NOVAS_MOON, NOVAS_ORBIT_INIT, novas_orbit_posvel(), NOVAS_REDUCED_ACCURACY, and NOVAS_SUN.

◆ novas_geom_posvel()

int novas_geom_posvel ( const object *restrict source,
const novas_frame *restrict frame,
enum novas_reference_system sys,
double *restrict pos,
double *restrict vel )

Calculates the geometric position and velocity vectors, relative to the observer, for a source in the given observing frame, in the specified coordinate system of choice.

The geometric position includes proper motion, and for solar-system bodies it is antedated for light travel time, so it effectively represents the geometric position as seen by the observer. However, the geometric does not include aberration correction, nor gravitational deflection.

If you want apparent positions, which account for aberration and gravitational deflection, use novas_skypos() instead.

You can also use novas_transform_vector() to convert the output position and velocity vectors to a dfferent coordinate system of choice afterwards if you want the results expressed in more than one coordinate system.

It implements the same geometric transformations as place() but at a reduced computational cost. See place() for references.

NOTES:

  1. As of SuperNOVAS v1.3, the returned velocity vector is a proper observer-based velocity measure. In prior releases, and in NOVAS C 3.1, this was inconsistent, with pseudo LSR-based measures being returned for catalog sources.
Parameters
sourcePointer to a celestial source data structure that is observed. Catalog sources should have coordinates and properties in ICRS. You can use transform_cat() to convert catalog entries to ICRS as necessary.
frameObserver frame, defining the location and time of observation
sysThe coordinate system in which to return positions and velocities.
[out]pos[AU] Calculated geometric position vector of the source relative to the observer location, in the designated coordinate system. It may be NULL if not required.
[out]vel[AU/day] The calculated velocity vector of the source relative to the observer in the designated coordinate system. It must be distinct from the pos output vector, and may be NULL if not required.
Returns
0 if successful, or else -1 if any of the arguments is invalid, 50–70 error is 50 + error from light_time2().
Since
1.1
Author
Attila Kovacs
See also
novas_geom_to_app(), novas_sky_pos(),
novas_transform_vector(), novas_make_frame()

References bary2obs(), d_light(), light_time2(), NOVAS_CATALOG_OBJECT, NOVAS_FULL_ACCURACY, novas_get_time(), NOVAS_JD_J2000, NOVAS_PLANET, NOVAS_REDUCED_ACCURACY, NOVAS_TDB, proper_motion(), and starvectors().

◆ novas_geom_to_app()

int novas_geom_to_app ( const novas_frame *restrict frame,
const double *restrict pos,
enum novas_reference_system sys,
sky_pos *restrict out )

Converts an geometric position in ICRS to an apparent position on sky, by applying appropriate corrections for aberration and gravitational deflection for the observer's frame.

Unlike place() the output reports the distance calculated from the parallax for sidereal sources. The radial velocity of the output is set to NAN (if needed use novas_sky_pos() instead).

Parameters
frameThe observer frame, defining the location and time of observation
pos[AU] Geometric position of source in ICRS coordinates
sysThe coordinate system in which to return the apparent sky location
[out]outPointer to the data structure which is populated with the calculated apparent location in the designated coordinate system. It may be the same pounter as the input position.
Returns
0 if successful, or an error from grav_def2(), or else -1 (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_app_to_geom()
novas_app_to_hor(), novas_sky_pos(), novas_geom_posvel()

References grav_planets(), NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, novas_vlen(), and vector2radec().

◆ 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.

◆ novas_transform_vector()

int novas_transform_vector ( const double * in,
const novas_transform *restrict transform,
double * out )

Transforms a position or velocity 3-vector from one coordinate reference system to another.

Parameters
inInput 3-vector in the original coordinate reference system
transformPointer to a coordinate transformation matrix
[out]outOutput 3-vector in the new coordinate reference system. It may be the same as the input.
Returns
0 if successful, or else -1 if there was an error (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_make_transform(), novas_transform_skypos()

◆ obs_posvel()

int obs_posvel ( double jd_tdb,
double ut1_to_tt,
enum novas_accuracy accuracy,
const observer *restrict obs,
const double *restrict geo_pos,
const double *restrict geo_vel,
double *restrict pos,
double *restrict vel )

Calculates the ICRS position and velocity of the observer relative to the Solar System Barycenter (SSB).

Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian date.
ut1_to_tt[s] TT - UT1 time difference. Used only when 'location->where' is NOVAS_OBSERVER_ON_EARTH (1) or NOVAS_OBSERVER_IN_EARTH_ORBIT (2).
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
obsThe observer location, relative to which the output positions and velocities are to be calculated
geo_pos[AU] ICRS position vector of the geocenter w.r.t. the Solar System Barycenter (SSB). If either geo_pos or geo_vel is NULL, it will be calculated when needed.
geo_vel[AU/day] ICRS velocity vector of the geocenter w.r.t. the Solar System Barycenter (SSB). If either geo_pos or geo_vel is NULL, it will be calculated when needed.
[out]pos[AU] Position 3-vector of the observer w.r.t. the Solar System Barycenter (SSB). It may be NULL if not required.
[out]vel[AU/day] Velocity 3-vector of the observer w.r.t. the Solar System Barycenter (SSB). It must be distinct from the pos output vector, and may be NULL if not required.
Returns
0 if successful, or the error from geo_posvel(), or else -1 (with errno indicating the type of error).
Author
Attila Kovacs
Since
1.3
See also
novas_make_frame()

References ephemeris(), geo_posvel(), NOVAS_AIRBORNE_OBSERVER, NOVAS_BARYCENTER, NOVAS_EARTH_INIT, NOVAS_OBSERVER_IN_EARTH_ORBIT, NOVAS_OBSERVER_ON_EARTH, NOVAS_OBSERVER_PLACES, NOVAS_SOLAR_SYSTEM_OBSERVER, and tt2tdb().