SuperNOVAS v1.1
The NOVAS C library, made better
Loading...
Searching...
No Matches
super.c File Reference

Functions

double app_to_cirs_ra (double jd_tt, enum novas_accuracy accuracy, double ra)
 
double cirs_to_app_ra (double jd_tt, enum novas_accuracy accuracy, double ra)
 
int cirs_to_itrs (double jd_tt_high, double jd_tt_low, double ut1_to_tt, enum novas_accuracy accuracy, double xp, double yp, const double *in, double *out)
 
int cirs_to_tod (double jd_tt, enum novas_accuracy accuracy, const double *in, double *out)
 
int ecl2equ (double jd_tt, enum novas_equator_type coord_sys, enum novas_accuracy accuracy, double elon, double elat, double *ra, double *dec)
 
int gal2equ (double glon, double glat, double *ra, double *dec)
 
double get_ut1_to_tt (int leap_seconds, double dut1)
 
double get_utc_to_tt (int leap_seconds)
 
int grav_undef (double jd_tdb, enum novas_accuracy accuracy, const double *pos_app, const double *pos_obs, double *out)
 
int grav_undo_planets (const double *pos_app, const double *pos_obs, const novas_planet_bundle *planets, double *out)
 
int hor_to_itrs (const on_surface *location, double az, double za, double *itrs)
 
int itrs_to_cirs (double jd_tt_high, double jd_tt_low, double ut1_to_tt, enum novas_accuracy accuracy, double xp, double yp, const double *in, double *out)
 
int itrs_to_hor (const on_surface *location, const double *itrs, double *az, double *za)
 
int itrs_to_tod (double jd_tt_high, double jd_tt_low, double ut1_to_tt, enum novas_accuracy accuracy, double xp, double yp, const double *in, double *out)
 
int j2000_to_gcrs (const double *in, double *out)
 
int make_airborne_observer (const on_surface *location, const double *vel, observer *obs)
 
int make_cat_object (const cat_entry *star, object *source)
 
int make_ephem_object (const char *name, long num, object *body)
 
int make_solar_system_observer (const double *sc_pos, const double *sc_vel, observer *obs)
 
int place_cirs (double jd_tt, const object *source, enum novas_accuracy accuracy, sky_pos *pos)
 
int place_gcrs (double jd_tt, const object *source, enum novas_accuracy accuracy, sky_pos *pos)
 
int place_icrs (double jd_tt, const object *source, enum novas_accuracy accuracy, sky_pos *pos)
 
int place_j2000 (double jd_tt, const object *source, enum novas_accuracy accuracy, sky_pos *pos)
 
int place_mod (double jd_tt, const object *source, enum novas_accuracy accuracy, sky_pos *pos)
 
int place_tod (double jd_tt, const object *source, enum novas_accuracy accuracy, sky_pos *pos)
 
int tod_to_cirs (double jd_tt, enum novas_accuracy accuracy, const double *in, double *out)
 
int tod_to_itrs (double jd_tt_high, double jd_tt_low, double ut1_to_tt, enum novas_accuracy accuracy, double xp, double yp, const double *in, double *out)
 

Detailed Description

Date
Created on Aug 24, 2024
Author
Attila Kovacs

SuperNOVAS only functions, which are not integral to the functionality of novas.c, and thus can live in a separate, more manageably sized, module.

Function Documentation

◆ app_to_cirs_ra()

double app_to_cirs_ra ( double  jd_tt,
enum novas_accuracy  accuracy,
double  ra 
)

Converts an apparent right ascension coordinate (measured from the true equinox of date) to a CIRS R.A., measured from the CIO.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
ra[h] the apparent R.A. coordinate measured from the true equinox of date.
Returns
[h] The CIRS right ascension coordinate, measured from the CIO [0:24], or NAN if the accuracy is invalid, or if there wan an error from cio_ra().
See also
cirs_to_app_ra()
tod_to_cirs()
Since
1.0.1
Author
Attila Kovacs

References cio_ra().

◆ cirs_to_app_ra()

double cirs_to_app_ra ( double  jd_tt,
enum novas_accuracy  accuracy,
double  ra 
)

Converts a CIRS right ascension coordinate (measured from the CIO) to an apparent R.A. measured from the true equinox of date.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
ra[h] The CIRS right ascension coordinate, measured from the CIO.
Returns
[h] the apparent R.A. coordinate measured from the true equinox of date [0:24], or NAN if the accuracy is invalid, or if there wan an error from cio_ra().
See also
app_to_cirs_ra()
cirs_to_tod()
Since
1.0.1
Author
Attila Kovacs

References cio_ra().

◆ cirs_to_itrs()

int cirs_to_itrs ( double  jd_tt_high,
double  jd_tt_low,
double  ut1_to_tt,
enum novas_accuracy  accuracy,
double  xp,
double  yp,
const double *  in,
double *  out 
)

Rotates a position vector from the dynamical CIRS frame of date to the Earth-fixed ITRS frame (IAU 2000 standard method).

If both 'xp' and 'yp' are set to 0 no polar motion is included in the transformation.

If extreme (sub-microarcsecond) accuracy is not required, you can use UT1-based Julian date instead of the TT-based Julian date and set the 'ut1_to_tt' argument to 0.0. and you can use UTC-based Julian date the same way.for arcsec-level precision also.

REFERENCES:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Kaplan, G. H. (2003), 'Another Look at Non-Rotating Origins', Proceedings of IAU XXV Joint Discussion 16.
Parameters
jd_tt_high[day] High-order part of Terrestrial Time (TT) based Julian date.
jd_tt_low[day] Low-order part of 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)
xp[arcsec] Conventionally-defined X coordinate of celestial intermediate pole with respect to ITRS pole, in arcseconds.
yp[arcsec] Conventionally-defined Y coordinate of celestial intermediate pole with respect to ITRS pole, in arcseconds.
inPosition vector, geocentric equatorial rectangular coordinates, referred to CIRS axes (celestial system).
[out]outPosition vector, geocentric equatorial rectangular coordinates, referred to ITRS axes (terrestrial system).
Returns
0 if successful, -1 if either of the vector arguments is NULL, 1 if 'accuracy' is invalid, 2 if 'method' is invalid 10–20, 3 if the method and option are mutually incompatible, or else 10 + the error from cio_location(), or 20 + error from cio_basis().
See also
tod_to_itrs()
itrs_to_cirs()
gcrs_to_cirs()
cirs_to_gcrs()
cirs_to_tod()
Since
1.0
Author
Attila Kovacs

References cel2ter(), EROT_ERA, and NOVAS_DYNAMICAL_CLASS.

◆ cirs_to_tod()

int cirs_to_tod ( double  jd_tt,
enum novas_accuracy  accuracy,
const double *  in,
double *  out 
)

Transforms a rectangular equatorial (x, y, z) vector from the Celestial Intermediate Reference System (CIRS) at the given epoch to the True of Date (TOD) reference system.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date that defines the output epoch. Typically it does not require much precision, and Julian dates in other time measures will be unlikely to affect the result
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
inCIRS Input (x, y, z) position or velocity vector
[out]outOutput position or velocity 3-vector in the True of Date (TOD) frame. It can be the same vector as the input.
Returns
0 if successful, or -1 if either of the vector arguments is NULL or the accuracy is invalid, or 10 + the error from cio_location(), or else 20 + the error from cio_basis().
See also
tod_to_cirs()
cirs_to_app_ra()
cirs_to_gcrs()
cirs_to_itrs()
Since
1.1
Author
Attila Kovacs

References cio_ra(), and spin().

◆ ecl2equ()

int ecl2equ ( double  jd_tt,
enum novas_equator_type  coord_sys,
enum novas_accuracy  accuracy,
double  elon,
double  elat,
double *  ra,
double *  dec 
)

Convert ecliptic longitude and latitude to right ascension and declination. To convert GCRS ecliptic coordinates (mean ecliptic and equinox of J2000.0), set 'coord_sys' to NOVAS_GCRS_EQUATOR(2); in this case the value of 'jd_tt' can be set to anything, since J2000.0 is assumed. Otherwise, all input coordinates are dynamical at'jd_tt'.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date. (Unused if 'coord_sys' is NOVAS_GCRS_EQUATOR[2])
coord_sysThe astrometric reference system of the coordinates. If 'coord_sys' is NOVAS_GCRS_EQUATOR(2), the input GCRS coordinates are converted to J2000 ecliptic coordinates.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
elon[deg] Ecliptic longitude in degrees, referred to specified ecliptic and equinox of date.
elat[deg] Ecliptic latitude in degrees, referred to specified ecliptic and equinox of date.
[out]ra[h] Right ascension in hours, referred to specified equator and equinox of date.
[out]dec[deg] Declination in degrees, referred to specified equator and equinox of date.
Returns
0 if successful, or else 1 if the value of 'coord_sys' is invalid.
See also
ecl2equ_vec()
equ2ecl()
Since
1.0
Author
Attila Kovacs

References ecl2equ_vec().

◆ gal2equ()

int gal2equ ( double  glon,
double  glat,
double *  ra,
double *  dec 
)

Converts galactic longitude and latitude to ICRS right ascension and declination.

REFERENCES:

  1. Hipparcos and Tycho Catalogues, Vol. 1, Section 1.5.3.
Parameters
glon[deg] Galactic longitude in degrees.
glat[deg] Galactic latitude in degrees.
[out]ra[h] ICRS right ascension in hours.
[out]dec[deg] ICRS declination in degrees.
Returns
0 if successful, or -1 if either of the output pointer arguments are NULL.
See also
equ2gal()
Since
1.0
Author
Attila Kovacs

◆ get_ut1_to_tt()

double get_ut1_to_tt ( int  leap_seconds,
double  dut1 
)

Returns the TT - UT1 time difference given the leap seconds and the actual UT1 - UTC time difference as measured and published by IERS.

NOTES:

  1. The current UT1 - UTC time difference, and polar offsets, historical data and near-term projections are published in the <a href="https://www.iers.org/IERS/EN/Publications/Bulletins/bulletins.html>IERS Bulletins
Parameters
leap_seconds[s] Leap seconds at the time of observations
dut1[s] UT1 - UTC time difference [-0.5:0.5]
Returns
[s] The TT - UT1 time difference that is suitable for used with all calls in this library that require a ut1_to_tt argument.
See also
get_utc_to_tt()
place()
cel_pole()
Since
1.0
Author
Attila Kovacs

References get_utc_to_tt().

◆ get_utc_to_tt()

double get_utc_to_tt ( int  leap_seconds)

Returns the difference between Terrestrial Time (TT) and Universal Coordinated Time (UTC)

Parameters
leap_seconds[s] The current leap seconds (see IERS Bulletins)
Returns
[s] The TT - UTC time difference
See also
get_ut1_to_tt()
julian_date()
Since
1.0
Author
Attila Kovacs

References NOVAS_TAI_TO_TT.

◆ grav_undef()

int grav_undef ( double  jd_tdb,
enum novas_accuracy  accuracy,
const double *  pos_app,
const double *  pos_obs,
double *  out 
)

Computes the gravitationally undeflected position of an observed source position due to the major gravitating bodies in the solar system. This function valid for an observed body within the solar system as well as for a star.

If 'accuracy' is set to zero (full accuracy), three bodies (Sun, Jupiter, and Saturn) are used in the calculation. If the reduced-accuracy option is set, only the Sun is used in the calculation. In both cases, if the observer is not at the geocenter, the deflection due to the Earth is included.

The number of bodies used at full and reduced accuracy can be set by making a change to the code in this function as indicated in the comments.

REFERENCES:

  1. Klioner, S. (2003), Astronomical Journal 125, 1580-1597, Section 6.
Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian date
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
pos_app[AU] Apparent position 3-vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, components in AU.
pos_obs[AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes, components in AU.
[out]out[AU] Nominal position vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, without gravitational deflection, components in AU. It can be the same vector as the input, but not the same as pos_obs.
Returns
0 if successful, -1 if any of the pointer arguments is NULL (errno = EINVAL) or if the result did not converge (errno = ECANCELED), or else an error from obs_planets().
See also
grav_def()
novas_app_to_geom()
set_planet_provider()
set_planet_provider_hp()
grav_bodies_full_accuracy
grav_bodies_reduced_accuracy
Since
1.1
Author
Attila Kovacs

References grav_bodies_full_accuracy, grav_bodies_reduced_accuracy, grav_undo_planets(), NOVAS_FULL_ACCURACY, and obs_planets().

◆ grav_undo_planets()

int grav_undo_planets ( const double *  pos_app,
const double *  pos_obs,
const novas_planet_bundle planets,
double *  out 
)

Computes the gravitationally undeflected position of an observed source position due to the specified Solar-system bodies.

REFERENCES:

  1. Klioner, S. (2003), Astronomical Journal 125, 1580-1597, Section 6.
Parameters
pos_app[AU] Apparent position 3-vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, components in AU.
pos_obs[AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes, components in AU.
planetsApparent planet data containing positions and velocities for the major gravitating bodies in the solar-system.
[out]out[AU] Nominal position vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, without gravitational deflection, components in AU. It can be the same vector as the input, but not the same as pos_obs.
Returns
0 if successful, -1 if any of the pointer arguments is NULL.
See also
obs_planets()
grav_planets()
novas_app_to_geom()
Since
1.1
Author
Attila Kovacs

References grav_planets(), and novas_inv_max_iter.

◆ hor_to_itrs()

int hor_to_itrs ( const on_surface location,
double  az,
double  za,
double *  itrs 
)

Converts astrometric (unrefracted) azimuth and zenith angles at the specified observer location to a unit position vector in the Earth-fixed ITRS frame.

Parameters
locationObserver location on Earth
az[deg] astrometric azimuth angle at observer location [0:360]. It may be NULL if not required.
za[deg] astrometric zenith angle at observer location [0:180]. It may be NULL if not required.
[out]itrsUnit 3-vector direction in Earth-fixed ITRS frame
Returns
0 if successful, or else -1 if the location or the input vector is NULL.
See also
itrs_to_hor()
itrs_to_cirs()
itrs_to_tod()
refract()
Since
1.0
Author
Attila Kovacs

References on_surface::latitude, and on_surface::longitude.

◆ itrs_to_cirs()

int itrs_to_cirs ( double  jd_tt_high,
double  jd_tt_low,
double  ut1_to_tt,
enum novas_accuracy  accuracy,
double  xp,
double  yp,
const double *  in,
double *  out 
)

Rotates a position vector from the Earth-fixed ITRS frame to the dynamical CIRS frame of date (IAU 2000 standard method).

If both 'xp' and 'yp' are set to 0 no polar motion is included in the transformation.

If extreme (sub-microarcsecond) accuracy is not required, you can use UT1-based Julian date instead of the TT-based Julian date and set the 'ut1_to_tt' argument to 0.0. and you can use UTC-based Julian date the same way.for arcsec-level precision also.

REFERENCES:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Kaplan, G. H. (2003), 'Another Look at Non-Rotating Origins', Proceedings of IAU XXV Joint Discussion 16.
Parameters
jd_tt_high[day] High-order part of Terrestrial Time (TT) based Julian date.
jd_tt_low[day] Low-order part of 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)
xp[arcsec] Conventionally-defined X coordinate of celestial intermediate pole with respect to ITRS pole, in arcseconds.
yp[arcsec] Conventionally-defined Y coordinate of celestial intermediate pole with respect to ITRS pole, in arcseconds.
inPosition vector, geocentric equatorial rectangular coordinates, referred to ITRS axes (terrestrial system)
[out]outPosition vector, geocentric equatorial rectangular coordinates, referred to CIRS axes (celestial system).
Returns
0 if successful, -1 if either of the vector arguments is NULL, 1 if 'accuracy' is invalid, or else 10 + the error from cio_location(), or 20 + error from cio_basis().
See also
itrs_to_tod()
cirs_to_itrs()
cirs_to_gcrs()
Since
1.0
Author
Attila Kovacs

References EROT_ERA, NOVAS_DYNAMICAL_CLASS, and ter2cel().

◆ itrs_to_hor()

int itrs_to_hor ( const on_surface location,
const double *  itrs,
double *  az,
double *  za 
)

Converts a position vector in the Earth-fixed ITRS frame to astrometric (unrefracted) azimuth and zenith angles at the specified observer location.

Parameters
locationObserver location on Earth
itrs3-vector position in Earth-fixed ITRS frame
[out]az[deg] astrometric azimuth angle at observer location [0:360]. It may be NULL if not required.
[out]za[deg] astrometric zenith angle at observer location [0:180]. It may be NULL if not required.
Returns
0 if successful, or else -1 if the location or the input vector is NULL.
See also
hor_to_itrs()
cirs_to_itrs()
tod_to_itrs()
refract_astro()
Since
1.0
Author
Attila Kovacs

References on_surface::latitude, and on_surface::longitude.

◆ itrs_to_tod()

int itrs_to_tod ( double  jd_tt_high,
double  jd_tt_low,
double  ut1_to_tt,
enum novas_accuracy  accuracy,
double  xp,
double  yp,
const double *  in,
double *  out 
)

Rotates a position vector from the Earth-fixed ITRS frame to the dynamical True of Date (TOD) frame of date (pre IAU 2000 method).

If both 'xp' and 'yp' are set to 0 no polar motion is included in the transformation.

If extreme (sub-microarcsecond) accuracy is not required, you can use UT1-based Julian date instead of the TT-based Julian date and set the 'ut1_to_tt' argument to 0.0. and you can use UTC-based Julian date the same way.for arcsec-level precision also.

REFERENCES:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Kaplan, G. H. (2003), 'Another Look at Non-Rotating Origins', Proceedings of IAU XXV Joint Discussion 16.
Parameters
jd_tt_high[day] High-order part of Terrestrial Time (TT) based Julian date.
jd_tt_low[day] Low-order part of 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)
xp[arcsec] Conventionally-defined X coordinate of celestial intermediate pole with respect to ITRS pole, in arcseconds.
yp[arcsec] Conventionally-defined Y coordinate of celestial intermediate pole with respect to ITRS pole, in arcseconds.
inPosition vector, geocentric equatorial rectangular coordinates, referred to ITRS axes (terrestrial system)
[out]outPosition vector, geocentric equatorial rectangular coordinates, referred to True of Date (TOD) axes (celestial system)
Returns
0 if successful, -1 if either of the vector arguments is NULL, 1 if 'accuracy' is invalid, or else 10 + the error from cio_location(), or 20 + error from cio_basis().
See also
itrs_to_cirs()
tod_to_itrs()
tod_to_j2000()
Since
1.0
Author
Attila Kovacs

References EROT_GST, NOVAS_DYNAMICAL_CLASS, and ter2cel().

◆ j2000_to_gcrs()

int j2000_to_gcrs ( const double *  in,
double *  out 
)

Change J2000 coordinates to GCRS coordinates. Same as frame_tie() called with J2000_TO_ICRS

Parameters
inJ2000 input 3-vector
[out]outGCRS output 3-vector
Returns
0 if successful, or else an error from frame_tie()
See also
j2000_to_tod()
gcrs_to_j2000()
Since
1.0
Author
Attila Kovacs

References frame_tie(), and J2000_TO_ICRS.

◆ make_airborne_observer()

int make_airborne_observer ( const on_surface location,
const double *  vel,
observer obs 
)

Populates an 'observer' data structure for an observer moving relative to the surface of Earth, such as an airborne observer. Airborne observers have an earth fixed momentary location, defined by longitude, latitude, and altitude, the same was as for a stationary observer on Earth, but are moving relative to the surface, such as in an aircraft or balloon observatory.

Parameters
locationCurrent longitude, latitude and altitude, and local weather (temperature and pressure)
vel[km/s] Surface velocity.
[out]obsPointer to data structure to populate.
Returns
0 if successful, or -1 if the output argument is NULL.
See also
make_observer_at geocenter()
make_observer_in_space()
make_observer_on_surface()
make_solar_system_observer()
novas_calc_geometric_position()
place()
Since
1.1
Author
Attila Kovacs

References IN_SPACE_INIT, make_observer(), NOVAS_AIRBORNE_OBSERVER, and in_space::sc_vel.

◆ make_cat_object()

int make_cat_object ( const cat_entry star,
object source 
)

Populates and object data structure with the data for a catalog source.

Parameters
starPointer to structure to populate with the catalog data for a celestial object located outside the solar system.
[out]sourcePointer to the celestial object data structure to be populated.
Returns
0 if successful, or -1 if 'cel_obj' is NULL or when type is NOVAS_CATALOG_OBJECT and 'star' is NULL, or else 1 if 'type' is invalid, 2 if 'number' is out of legal range or 5 if 'name' is too long.
See also
make_cat_entry()
make_planet()
make_ephem_object()
place()
Since
1.1
Author
Attila Kovacs

References make_object(), NOVAS_CATALOG_OBJECT, cat_entry::starname, and cat_entry::starnumber.

◆ make_ephem_object()

int make_ephem_object ( const char *  name,
long  num,
object body 
)

Sets a celestial object to be a Solar-system ephemeris body. Typically this would be used to define minor planets, asteroids, comets and planetary satellites.

Parameters
nameName of object. By default converted to upper-case, unless novas_case_sensitive() was called with a non-zero argument. Max. SIZE_OF_OBJ_NAME long, including termination.
numSolar-system body ID number (e.g. NAIF)
[out]bodyPointer to structure to populate.
Returns
0 if successful, or else -1 if the 'planet' pointer is NULL or the name is too long.
See also
make_planet()
make_cat_entry()
place()
Since
1.0
Author
Attila Kovacs

References make_object(), and NOVAS_EPHEM_OBJECT.

◆ make_solar_system_observer()

int make_solar_system_observer ( const double *  sc_pos,
const double *  sc_vel,
observer obs 
)

Populates an 'observer' data structure, for an observer situated on a near-Earth spacecraft, with the specified geocentric position and velocity vectors. Solar-system observers are similar to observers in Earth-orbit but their momentary position and velocity is defined relative to the Solar System Barycenter, instead of the geocenter.

Parameters
sc_pos[AU] Solar-system barycentric (x, y, z) position vector in ICRS.
sc_vel[AU/day] Solar-system barycentric (x, y, z) velocity vector in ICRS.
[out]obsPointer to the data structure to populate
Returns
0 if successful, or -1 if the output argument is NULL.
See also
make_observer_in_space()
make_observer_on_surface()
make_observer_at_geocenter()
make_airborne_observer()
novas_calc_geometric_position()
place()
Since
1.1
Author
Attila Kovacs

References make_in_space(), make_observer(), and NOVAS_SOLAR_SYSTEM_OBSERVER.

◆ place_cirs()

int place_cirs ( double  jd_tt,
const object source,
enum novas_accuracy  accuracy,
sky_pos pos 
)

Computes the Celestial Intermediate Reference System (CIRS) dynamical position position of a source as 'seen' from the geocenter at the given time of observation. See place() for more information.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated CIRS position data
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
See also
place_tod()
place_gcrs()
Since
1.0
Author
Attila Kovacs

References NOVAS_CIRS, and place().

◆ place_gcrs()

int place_gcrs ( double  jd_tt,
const object source,
enum novas_accuracy  accuracy,
sky_pos pos 
)

Computes the Geocentric Celestial Reference System (GCRS) position of a source (as 'seen' from the geocenter) at the given time of observation. Unlike place_icrs(), this includes aberration for the moving frame of the geocenter as well as gravitational deflections calculated for a virtual observer located at the geocenter. See place() for more information.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated GCRS position data
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
See also
place_icrs()
place_cirs()
place_tod()
virtual_star()
virtual_planet()
Since
1.0
Author
Attila Kovacs

References NOVAS_GCRS, and place().

◆ place_icrs()

int place_icrs ( double  jd_tt,
const object source,
enum novas_accuracy  accuracy,
sky_pos pos 
)

Computes the International Celestial Reference System (ICRS) position of a source. (from the geocenter). Unlike place_gcrs(), this version does not include aberration or gravitational deflection corrections.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated geocentric ICRS position data (Unlike place_gcrs(), the calculated coordinates do not account for aberration or gravitational deflection).
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
See also
place_gcrs()
place_cirs()
place_tod()
mean_star()
Since
1.0
Author
Attila Kovacs

References NOVAS_ICRS, and place().

◆ place_j2000()

int place_j2000 ( double  jd_tt,
const object source,
enum novas_accuracy  accuracy,
sky_pos pos 
)

Computes the J2000 dynamical position position of a source as 'seen' from the geocenter at the given time of observation. See place() for more information.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated CIRS position data
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
See also
place_cirs()
place_gcrs()
app_star()
app_planet()
Since
1.1
Author
Attila Kovacs

References NOVAS_J2000, and place().

◆ place_mod()

int place_mod ( double  jd_tt,
const object source,
enum novas_accuracy  accuracy,
sky_pos pos 
)

Computes the Mean of Date (MOD) dynamical position position of a source as 'seen' from the geocenter at the given time of observation. See place() for more information.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated CIRS position data
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
See also
place_cirs()
place_gcrs()
app_star()
app_planet()
Since
1.1
Author
Attila Kovacs

References NOVAS_MOD, and place().

◆ place_tod()

int place_tod ( double  jd_tt,
const object source,
enum novas_accuracy  accuracy,
sky_pos pos 
)

Computes the True of Date (TOD) dynamical position position of a source as 'seen' from the geocenter at the given time of observation. See place() for more information.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated CIRS position data
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
See also
place_cirs()
place_gcrs()
app_star()
app_planet()
Since
1.0
Author
Attila Kovacs

References NOVAS_TOD, and place().

◆ tod_to_cirs()

int tod_to_cirs ( double  jd_tt,
enum novas_accuracy  accuracy,
const double *  in,
double *  out 
)

Transforms a rectangular equatorial (x, y, z) vector from the True of Date (TOD) reference system to the Celestial Intermediate Reference System (CIRS) at the given epoch to the .

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date that defines the output epoch. Typically it does not require much precision, and Julian dates in other time measures will be unlikely to affect the result
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
inCIRS Input (x, y, z) position or velocity vector
[out]outOutput position or velocity 3-vector in the True of Date (TOD) frame. It can be the same vector as the input.
Returns
0 if successful, or -1 if either of the vector arguments is NULL or the accuracy is invalid, or 10 + the error from cio_location(), or else 20 + the error from cio_basis().
See also
cirs_to_tod()
app_to_cirs_ra()
tod_to_gcrs()
tod_to_j2000()
tod_to_itrs()
Since
1.1
Author
Attila Kovacs

References cio_ra(), and spin().

◆ tod_to_itrs()

int tod_to_itrs ( double  jd_tt_high,
double  jd_tt_low,
double  ut1_to_tt,
enum novas_accuracy  accuracy,
double  xp,
double  yp,
const double *  in,
double *  out 
)

Rotates a position vector from the dynamical True of Date (TOD) frame of date the Earth-fixed ITRS frame (pre IAU 2000 method).

If both 'xp' and 'yp' are set to 0 no polar motion is included in the transformation.

If extreme (sub-microarcsecond) accuracy is not required, you can use UT1-based Julian date instead of the TT-based Julian date and set the 'ut1_to_tt' argument to 0.0. and you can use UTC-based Julian date the same way.for arcsec-level precision also.

REFERENCES:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Kaplan, G. H. (2003), 'Another Look at Non-Rotating Origins', Proceedings of IAU XXV Joint Discussion 16.
Parameters
jd_tt_high[day] High-order part of Terrestrial Time (TT) based Julian date.
jd_tt_low[day] Low-order part of 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)
xp[arcsec] Conventionally-defined X coordinate of celestial intermediate pole with respect to ITRS pole, in arcseconds.
yp[arcsec] Conventionally-defined Y coordinate of celestial intermediate pole with respect to ITRS pole, in arcseconds.
inPosition vector, geocentric equatorial rectangular coordinates, referred to True of Date (TOD) axes (celestial system).
[out]outPosition vector, geocentric equatorial rectangular coordinates, referred to ITRS axes (terrestrial system).
Returns
0 if successful, -1 if either of the vector arguments is NULL, 1 if 'accuracy' is invalid, 2 if 'method' is invalid 10–20, 3 if the method and option are mutually incompatible, or else 10 + the error from cio_location(), or 20 + error from cio_basis().
See also
cirs_to_itrs()
itrs_to_tod()
j2000_to_tod()
tod_to_gcrs()
tod_to_j2000()
tod_to_cirs()
Since
1.0
Author
Attila Kovacs

References cel2ter(), EROT_GST, and NOVAS_DYNAMICAL_CLASS.