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

Macros

#define _DEFAULT_SOURCE
 strcasecmp() feature macro starting glibc 2.20 (2014-09-08)
 

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 *restrict ra, double *restrict dec)
 
int gal2equ (double glon, double glat, double *restrict ra, double *restrict dec)
 
double get_ut1_to_tt (int leap_seconds, double dut1)
 
double get_utc_to_tt (int leap_seconds)
 
double grav_redshift (double M_kg, double r_m)
 
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 *restrict planets, double *out)
 
int hor_to_itrs (const on_surface *restrict location, double az, double za, double *restrict 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 *restrict location, const double *restrict itrs, double *restrict az, double *restrict 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_cat_object_sys (const cat_entry *star, const char *restrict system, object *source)
 
int make_ephem_object (const char *name, long num, object *body)
 
int make_orbital_object (const char *name, long num, const novas_orbital *orbit, object *body)
 
int make_redshifted_cat_entry (const char *name, double ra, double dec, double z, cat_entry *source)
 
int make_redshifted_object (const char *name, double ra, double dec, double z, object *source)
 
int make_redshifted_object_sys (const char *name, double ra, double dec, const char *restrict system, double z, object *source)
 
int make_solar_system_observer (const double *sc_pos, const double *sc_vel, observer *obs)
 
int novas_e2h_offset (double dra, double ddec, double pa, double *restrict daz, double *restrict del)
 
double novas_epa (double ha, double dec, double lat)
 
double novas_epoch (const char *restrict system)
 
double novas_equ_sep (double ra1, double dec1, double ra2, double dec2)
 
int novas_h2e_offset (double daz, double del, double pa, double *restrict dra, double *restrict ddec)
 
double novas_helio_dist (double jd_tdb, const object *restrict source, double *restrict rate)
 
double novas_hpa (double az, double el, double lat)
 
double novas_lsr_to_ssb_vel (double epoch, double ra, double dec, double vLSR)
 
enum novas_planet novas_planet_for_name (const char *restrict name)
 
double novas_sep (double lon1, double lat1, double lon2, double lat2)
 
int novas_set_orbsys_pole (enum novas_reference_system type, double ra, double dec, novas_orbital_system *restrict sys)
 
double novas_solar_power (double jd_tdb, const object *restrict source)
 
double novas_ssb_to_lsr_vel (double epoch, double ra, double dec, double vLSR)
 
double novas_v2z (double vel)
 
int novas_xyz_to_uvw (const double *xyz, double ha, double dec, double *uvw)
 
double novas_z_add (double z1, double z2)
 
double novas_z_inv (double z)
 
int place_cirs (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos)
 
int place_gcrs (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos)
 
int place_icrs (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos)
 
int place_j2000 (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos)
 
int place_mod (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos)
 
int place_tod (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos)
 
double redshift_vrad (double vrad, double z)
 
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)
 
double unredshift_vrad (double vrad, double z)
 

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 *restrict  ra,
double *restrict  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 *restrict  ra,
double *restrict  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_redshift()

double grav_redshift ( double  M_kg,
double  r_m 
)

Returns the gravitational redshift (z) for light emitted near a massive spherical body at some distance from its center, and observed at some very large (infinite) distance away.

Parameters
M_kg[kg] Mass of gravitating body that is contained inside the emitting radius.
r_m[m] Radius at which light is emitted.
Returns
The gravitational redshift (z) for an observer at very large (infinite) distance from the gravitating body.
See also
redshift_vrad()
unredshift_vrad()
novas_z_add()
Since
1.2
Author
Attila Kovacs

References C.

◆ 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 *restrict  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 *restrict  location,
double  az,
double  za,
double *restrict  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

◆ 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 *restrict  location,
const double *restrict  itrs,
double *restrict  az,
double *restrict  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

◆ 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. The input source must be defined with ICRS coordinates. To create objects with other types of coordiantes use make_cat_object_epoch() instead.

Parameters
starPointer to structure to populate with the catalog data for a celestial object located outside the solar system, specified with ICRS coordinates.
[out]sourcePointer to the celestial object data structure to be populated.
Returns
0 if successful, or -1 if either argument is NULL, or else 5 if 'name' is too long.
See also
make_cat_object_sys()
make_cat_entry()
make_planet()
make_ephem_object()
novas_geom_posvel()
place()
Since
1.1
Author
Attila Kovacs

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

◆ make_cat_object_sys()

int make_cat_object_sys ( const cat_entry star,
const char *restrict  system,
object source 
)

Populates and object data structure with the data for a catalog source for a given system of catalog coordinates.

Parameters
starPointer to structure to populate with the catalog data for a celestial object located outside the solar system.
systemInput catalog coordinate system epoch, e.g. "ICRS", "B1950.0", "J2000.0", "FK4", "FK5", or "HIP". In general, any Besselian or Julian year epoch can be used by year (e.g. "B1933.193" or "J2022.033"), or else the fixed value listed. If 'B' or 'J' is ommitted in front of the epoch year, then Besselian epochs are assumed prior to 1984.0.
[out]sourcePointer to the celestial object data structure to be populated with the corresponding ICRS catalog coordinates, after appying proper-motion and precession corrections as appropriate.
Returns
0 if successful, or -1 if any argument is NULL or if the input 'system' is invalid, or else 5 if 'name' is too long.
See also
make_cat_object()
make_redshifted_object_sys()
novas_jd_for_epoch()
make_cat_entry()
place()
NOVAS_SYSTEM_ICRS
NOVAS_SYSTEM_HIP
NOVAS_SYSTEM_J2000
NOVAS_SYSTEM_B1950
Since
1.3
Author
Attila Kovacs

References make_cat_object(), and object::star.

◆ 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. If the ephemeris provider uses names, then the name should match those of the ephemeris provider – otherwise it is not important.
numSolar-system body ID number (e.g. NAIF). The number should match the needs of the ephemeris provider used with NOVAS. (If the ephemeris provider is by name and not ID number, then the number here is not important).
[out]bodyPointer to structure to populate.
Returns
0 if successful, or else -1 if the 'body' pointer is NULL or the name is too long.
See also
set_ephem_provider()
make_planet()
make_cat_entry()
novas_geom_posvel()
place()
Since
1.0
Author
Attila Kovacs

References make_object(), and NOVAS_EPHEM_OBJECT.

◆ make_orbital_object()

int make_orbital_object ( const char *  name,
long  num,
const novas_orbital orbit,
object body 
)

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

Parameters
nameName of object. It may be NULL if not relevant.
numSolar-system body ID number (e.g. NAIF). It is not required and can be set e.g. to -1 if not relevant to the caller.
orbitThe orbital parameters to adopt. The data will be copied, not referenced.
[out]bodyPointer to structure to populate.
Returns
0 if successful, or else -1 if the 'orbit' or 'body' pointer is NULL or the name is too long.
See also
novas_orbit_posvel()
make_planet()
make_ephem_object()
novas_geom_posvel()
place()
Since
1.2
Author
Attila Kovacs

References make_object(), NOVAS_ORBITAL_OBJECT, and object::orbit.

◆ make_redshifted_cat_entry()

int make_redshifted_cat_entry ( const char *  name,
double  ra,
double  dec,
double  z,
cat_entry source 
)

Populates a celestial object data structure with the parameters for a redhifted catalog source, such as a distant quasar or galaxy. It is similar to make_cat_object() except that it takes a Doppler-shift (z) instead of radial velocity and it assumes no parallax and no proper motion (appropriately for a distant redshifted source). The catalog name is set to EXT to indicate an extragalactic source, and the catalog number defaults to 0. The user may change these default field values as appropriate afterwards, if necessary.

Parameters
nameObject name (less than SIZE_OF_OBJ_NAME in length). It may be NULL.
ra[h] Right ascension of the object (hours).
dec[deg] Declination of the object (degrees).
zRedhift value (λobs / λrest - 1 = frest / fobs - 1).
[out]sourcePointer to structure to populate.
Returns
0 if successful, or 5 if 'name' is too long, else -1 if the 'source' pointer is NULL.
See also
make_redshifted_object_sys()
novas_v2z()
Since
1.2
Author
Attila Kovacs

References make_cat_entry(), and novas_z2v().

◆ make_redshifted_object()

int make_redshifted_object ( const char *  name,
double  ra,
double  dec,
double  z,
object source 
)

Populates a celestial object data structure with the parameters for a redhifted catalog source, such as a distant quasar or galaxy. It is similar to make_cat_object() except that it takes a Doppler-shift (z) instead of radial velocity and it assumes no parallax and no proper motion (appropriately for a distant redshifted source). The catalog name is set to EXT to indicate an extragalactic source, and the catalog number defaults to 0. The user may change these default field values as appropriate afterwards, if necessary.

Parameters
nameObject name (less than SIZE_OF_OBJ_NAME in length). It may be NULL.
ra[h] ICRS Right ascension of the object (hours).
dec[deg] ICRS Declination of the object (degrees).
zRedhift value (λobs / λrest - 1 = frest / fobs - 1).
[out]sourcePointer to structure to populate.
Returns
0 if successful, or 5 if 'name' is too long, else -1 if the 'source' pointer is NULL.
See also
make_redshifted_object_sys()
make_cat_object()
novas_v2z()
Since
1.2
Author
Attila Kovacs

References make_cat_object(), and make_redshifted_cat_entry().

◆ make_redshifted_object_sys()

int make_redshifted_object_sys ( const char *  name,
double  ra,
double  dec,
const char *restrict  system,
double  z,
object source 
)

Populates a celestial object data structure with the parameters for a redhifted catalog source, such as a distant quasar or galaxy, for a given system of catalog coordinates.

Parameters
nameObject name (less than SIZE_OF_OBJ_NAME in length). It may be NULL.
ra[h] ICRS Right ascension of the object (hours).
dec[deg] ICRS Declination of the object (degrees).
systemInput catalog coordinate system epoch, e.g. "ICRS", "B1950.0", "J2000.0", "FK4", "FK5", or "HIP". In general, any Besselian or Julian year epoch can be used by year (e.g. "B1933.193" or "J2022.033"), or else the fixed value listed. If 'B' or 'J' is ommitted in front of the epoch year, then Besselian epochs are assumed prior to 1984.0.
zRedhift value (λobs / λrest - 1 = frest / fobs - 1).
[out]sourcePointer to the celestial object data structure to be populated with the corresponding ICRS catalog coordinates.
Returns
0 if successful, or -1 if any of the pointer arguments is NULL or the 'system' is invalid, 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_redshifted_object()
make_cat_object_sys()
novas_jd_for_epoch()
place()
NOVAS_SYSTEM_ICRS
NOVAS_SYSTEM_HIP
NOVAS_SYSTEM_J2000
NOVAS_SYSTEM_B1950
Since
1.3
Author
Attila Kovacs

References make_redshifted_object(), and object::star.

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

◆ novas_e2h_offset()

int novas_e2h_offset ( double  dra,
double  ddec,
double  pa,
double *restrict  daz,
double *restrict  del 
)

Converts coordinate offsets, from the local equatorial system to local horizontal offsets. Converting between local flat projections and spherical coordinates usually requires a WCS projection.

REFERENCES:

  1. Calabretta, M.R., & Greisen, E.W., (2002), Astronomy & Astrophysics, 395, 1077-1122.
Parameters
dra[arcsec] Projected ffset position in the apparent true-of-date R.A. direction. E.g. The projected offset between two RA coordinates at a same reference declination, is δRA = (RA2 - RA1) * cos(Dec0)
ddec[arcsec] Projected offset position in the apparent true-of-date declination direction.
pa[deg] Parallactic Angle
[out]daz[arcsec] Output offset position in the local azimuth direction. It can be a pointer to one of the input coordinates, or NULL if not required.
[out]del[arcsec] Output offset position in the local elevation direction. It can be a pointer to one of the input coordinates, or NULL if not required.
Returns
0
Since
1.3
Author
Attila Kovacs
See also
novas_h2e_offset()
novas_epa()

References novas_h2e_offset().

◆ novas_epa()

double novas_epa ( double  ha,
double  dec,
double  lat 
)

Returns the equatorial Parallactic Angle (PA) calculated for an R.A./Dec location of the sky at a given sidereal time. The PA is the angle between the local horizontal coordinate directions and the local true-of-date equatorial coordinate directions, at the given location and time. The polar wobble is not included in the calculation.

The Parallactic Angle is sometimes referrred to as the Vertical Position Angle (VPA). Both define the same quantity.

Parameters
ha[h] Hour angle (LST - RA) i.e., the difference between the Local (apparent) Sidereal Time and the apparent (true-of-date) Right Ascension of observed source.
dec[deg] Apparent (true-of-date) declination of observed source
lat[deg] Geodetic latitude of observer
Returns
[deg] Parallactic Angle (PA). I.e., the clockwise position angle of the elevation direction w.r.t. the declination axis in the equatorial system. Same as the clockwise position angle of the declination direction w.r.t. the elevation axis, in the horizontal system.
Since
1.3
Author
Attila Kovacs
See also
novas_hpa()
novas_lst()
novas_e2h_offset()

◆ novas_epoch()

double novas_epoch ( const char *restrict  system)

Returns the Julian day corresponding to an astronomical coordinate epoch.

Parameters
systemCoordinate system, e.g. "ICRS", "B1950.0", "J2000.0", "FK4", "FK5", "1950", "2000", or "HIP". In general, any Besselian or Julian year epoch can be used by year (e.g. "B1933.193" or "J2022.033"), or else the fixed values listed. If 'B' or 'J' is ommitted in front of the epoch year, then Besselian epochs are assumed prior to 1984.0.
Returns
[day] The Julian day corresponding to the given coordinate epoch, or else NAN if the input string is NULL or the input is not recognised as a coordinate epoch specification (errno will be set to EINVAL).
Since
1.3
Author
Attila Kovacs
See also
make_cat_object_sys()
make_redshifted_object_sys()
transform_cat()
precession()
NOVAS_SYSTEM_ICRS
NOVAS_SYSTEM_B1950
NOVAS_SYSTEM_J2000
NOVAS_SYSTEM_HIP

References NOVAS_JD_B1950, NOVAS_JD_HIP, NOVAS_JD_J2000, NOVAS_SYSTEM_FK4, NOVAS_SYSTEM_FK5, NOVAS_SYSTEM_HIP, and NOVAS_SYSTEM_ICRS.

◆ novas_equ_sep()

double novas_equ_sep ( double  ra1,
double  dec1,
double  ra2,
double  dec2 
)

Returns the angular separation of two equatorial locations on a sphere.

Parameters
ra1[h] right ascension of first location
dec1[deg] declination of first location
ra2[h] right ascension of second location
dec2[deg] declination of second location
Returns
[deg] the angular separation of the two locations.
Since
1.3
Author
Attila Kovacs
See also
novas_sep()
novas_sun_angle()
novas_moon_angle()

References novas_sep().

◆ novas_h2e_offset()

int novas_h2e_offset ( double  daz,
double  del,
double  pa,
double *restrict  dra,
double *restrict  ddec 
)

Converts coordinate offsets, from the local horizontal system to local equatorial offsets. Converting between local flat projections and spherical coordinates usually requires a WCS projection.

REFERENCES:

  1. Calabretta, M.R., & Greisen, E.W., (2002), Astronomy & Astrophysics, 395, 1077-1122.
Parameters
daz[arcsec] Projected offset position in the azimuth direction. The projected offset between two azimuth positions at the same reference elevation is δAz = (Az2 - Az1) * cos(El0).
del[arcsec] projected offset position in the elevation direction
pa[deg] Parallactic Angle
[out]dra[arcsec] Output offset position in the local true-of-date R.A. direction. It can be a pointer to one of the input coordinates, or NULL if not required.
[out]ddec[arcsec] Output offset position in the local true-of-date declination direction. It can be a pointer to one of the input coordinates, or NULL if not required.
Returns
0
Since
1.3
Author
Attila Kovacs
See also
novas_e2h_offset()
novas_hpa()

◆ novas_helio_dist()

double novas_helio_dist ( double  jd_tdb,
const object *restrict  source,
double *restrict  rate 
)

Returns a Solar-system body's distance from the Sun, and optionally also the rate of recession. It may be useful, e.g. to calculate the body's heating from the Sun.

Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian date. You may want to use a time that is antedated to when the observed light originated from the source.
sourceObserved Solar-system source
[out]rate[AU/day] (optional) Returned rate of recession from Sun
Returns
[AU] Distance from the Sun, or NAN if not a Solar-system source.
Since
1.3
Author
Attila Kovacs
See also
novas_solar_power()
novas_solar_illum()

References ephemeris(), NOVAS_CATALOG_OBJECT, NOVAS_HELIOCENTER, and NOVAS_REDUCED_ACCURACY.

◆ novas_hpa()

double novas_hpa ( double  az,
double  el,
double  lat 
)

Returns the horizontal Parallactic Angle (PA) calculated for a gorizontal Az/El location of the sky. The PA is the angle between the local horizontal coordinate directions and the local true-of-date equatorial coordinate directions at the given location. The polar wobble is not included in the calculation.

The Parallactic Angle is sometimes referrred to as the Vertical Position Angle (VPA). Both define the same quantity.

Parameters
az[deg] Azimuth angle
el[deg] Elevation angle
lat[deg] Geodetic latitude of observer
Returns
[deg] Parallactic Angle (PA). I.e., the clockwise position angle of the declination direction w.r.t. the elevation axis in the horizontal system. Same as the the clockwise position angle of the elevation direction w.r.t. the declination axis in the equatorial system.
Since
1.3
Author
Attila Kovacs
See also
novas_epa()
novas_h2e_offset()

◆ novas_lsr_to_ssb_vel()

double novas_lsr_to_ssb_vel ( double  epoch,
double  ra,
double  dec,
double  vLSR 
)

Returns a Solar System Baricentric (SSB) radial velocity for a radial velocity that is referenced to the Local Standard of Rest (LSR). Internally, NOVAS always uses barycentric radial velocities, but it is just as common to have catalogs define radial velocities referenced to the LSR.

The SSB motion w.r.t. the barycenter is assumed to be (11.1, 12.24, 7.25) km/s in ICRS (Shoenrich et al. 2010).

REFERENCES:

  1. Ralph Schoenrich, James Binney, Walter Dehnen, Monthly Notices of the Royal Astronomical Society, Volume 403, Issue 4, April 2010, Pages 1829–1833, https://doi.org/10.1111/j.1365-2966.2010.16253.x
Parameters
epoch[yr] Coordinate epoch in which the coordinates and velocities are defined. E.g. 2000.0.
ra[h] Right-ascenscion of source at given epoch.
dec[deg] Declination of source at given epoch.
vLSR[km/s] radial velocity defined against the Local Standard of Rest (LSR), at given epoch.
Returns
[km/s] Equivalent Solar-System Barycentric radial velocity.
Since
1.3
Author
Attila Kovacs
See also
make_cat_entry()
novas_ssb_to_lsr_vel()

References NOVAS_JD_J2000, precession(), and radec2vector().

◆ novas_planet_for_name()

enum novas_planet novas_planet_for_name ( const char *restrict  name)

Returns the NOVAS planet ID for a given name (case insensitive), or -1 if no match is found.

Parameters
nameThe planet name, or that for the "Sun", "Moon" or "SSB" (case insensitive). The spelled out "Solar System Barycenter" is also recognized with either spaces, hyphens ('-') or underscores ('_') separating the case insensitive words.
Returns
The NOVAS major planet ID, or -1 (errno set to EINVAL) if the input name is NULL or if there is no match for the name provided.
Author
Attila Kovacs
Since
1.2
See also
make_planet()

References NOVAS_PLANET_NAMES_INIT, NOVAS_PLANETS, and NOVAS_SSB.

◆ novas_sep()

double novas_sep ( double  lon1,
double  lat1,
double  lon2,
double  lat2 
)

Returns the angular separation of two locations on a sphere.

Parameters
lon1[deg] longitude of first location
lat1[deg] latitude of first location
lon2[deg] longitude of second location
lat2[deg] latitude of second location
Returns
[deg] the angular separation of the two locations.
Since
1.3
Author
Attila Kovacs
See also
novas_equ_sep()
novas_sun_angle()
novas_moon_angle()

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

◆ novas_solar_power()

double novas_solar_power ( double  jd_tdb,
const object *restrict  source 
)

Returns the typical incident Solar power on a Solar-system body at the time of observation.

Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian date. You may want to use a time that is antedated to when the observed light originated ( was reflected) from the source.
sourceObserved Solar-system source
Returns
[W/m2] Incident Solar power on the illuminated side of the object, or NAN if not a Solar-system source or if the source is the Sun itself.
Since
1.3
Author
Attila Kovacs
See also
novas_solar_illum()

References novas_helio_dist(), and NOVAS_SOLAR_CONSTANT.

◆ novas_ssb_to_lsr_vel()

double novas_ssb_to_lsr_vel ( double  epoch,
double  ra,
double  dec,
double  vLSR 
)

Returns a radial-velocity referenced to the Local Standard of Rest (LSR) for a given Solar-System Barycentric (SSB) radial velocity. Internally, NOVAS always uses barycentric radial velocities, but it is just as common to have catalogs define radial velocities referenced to the LSR.

The SSB motion w.r.t. the barycenter is assumed to be (11.1, 12.24, 7.25) km/s in ICRS (Shoenrich et al. 2010).

REFERENCES:

  1. Ralph Schoenrich, James Binney, Walter Dehnen, Monthly Notices of the Royal Astronomical Society, Volume 403, Issue 4, April 2010, Pages 1829–1833, https://doi.org/10.1111/j.1365-2966.2010.16253.x
Parameters
epoch[yr] Coordinate epoch in which the coordinates and velocities are defined. E.g. 2000.0.
ra[h] Right-ascenscion of source at given epoch.
dec[deg] Declination of source at given epoch.
vLSR[km/s] radial velocity defined against the Local Standard of Rest (LSR), at given epoch.
Returns
[km/s] Equivalent Solar-System Barycentric radial velocity.
Since
1.3
Author
Attila Kovacs
See also
make_cat_entry()
novas_lsr_to_ssb_vel()

References NOVAS_JD_J2000, precession(), and radec2vector().

◆ novas_v2z()

double novas_v2z ( double  vel)

Converts a radial recession velocity to a redshift value (z = δf / frest). It is based on the relativistic formula:

 1 + z = sqrt((1 + β) / (1 - β))

where β = v / c.

Parameters
vel[km/s] velocity (i.e. rate) of recession.
Returns
the corresponding redshift value (δλ / λrest), or NAN if the input velocity is invalid (i.e., it exceeds the speed of light).
See also
novas_z2v()
novas_z_add()
Author
Attila Kovacs
Since
1.2

References C, and NOVAS_KMS.

◆ novas_xyz_to_uvw()

int novas_xyz_to_uvw ( const double *  xyz,
double  ha,
double  dec,
double *  uvw 
)

Converts rectangular telescope x,y,z (absolute or relative) coordinates (in ITRS) to equatorial u,v,w projected coordinates for a specified line of sight.

x,y,z are Cartesian coordinates w.r.t the Greenwich meridian. The directions are x: long=0, lat=0; y: long=90, lat=0; z: lat=90.

u,v,w are Cartesian coordinates (u,v) along the local equatorial R.A. and declination directions as seen from a direction on the sky (w).

Parameters
xyz[arb.u.] Absolute or relative x,y,z coordinates (double[3]).
ha[h] Hourangle (LST - RA) i.e., the difference between the Local (apparent) Sidereal Time and the apparent (true-of-date) Right Ascension of observed source.
dec[deg] Apparent (true-of-date) declination of source
[out]uvw[arb.u.] Converted u,v,w coordinates (double[3]) in same units as xyz. It may be the same vector as the input.
Returns
0 if successful, or else -1 if either vector argument is NULL (errno will be set to EINVAL)
Since
1.3
Author
Attila Kovacs

References novas_xyz_to_los().

◆ novas_z_add()

double novas_z_add ( double  z1,
double  z2 
)

Compounds two redshift corrections, e.g. to apply (or undo) a series gravitational redshift corrections and/or corrections for a moving observer. It's effectively using (1 + z) = (1 + z1) * (1 + z2).

Parameters
z1One of the redshift values
z2The other redshift value
Returns
The compound redshift value, ot NAN if either input redshift is invalid (errno will be set to EINVAL).
See also
grav_redshift()
redshift_vrad()
unredshift_vrad()
Since
1.2
Author
Attila Kovacs

◆ novas_z_inv()

double novas_z_inv ( double  z)

Returns the inverse of a redshift value, that is the redshift for a body moving with the same velocity as the original but in the opposite direction.

Parameters
zA redhift value
Returns
The redshift value for a body moving in the opposite direction with the same speed, or NAN if the input redshift is invalid.
See also
novas_z_add()
Since
1.2
Author
Attila Kovacs

◆ place_cirs()

int place_cirs ( double  jd_tt,
const object *restrict  source,
enum novas_accuracy  accuracy,
sky_pos *restrict  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 *restrict  source,
enum novas_accuracy  accuracy,
sky_pos *restrict  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 *restrict  source,
enum novas_accuracy  accuracy,
sky_pos *restrict  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 *restrict  source,
enum novas_accuracy  accuracy,
sky_pos *restrict  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 *restrict  source,
enum novas_accuracy  accuracy,
sky_pos *restrict  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 *restrict  source,
enum novas_accuracy  accuracy,
sky_pos *restrict  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().

◆ redshift_vrad()

double redshift_vrad ( double  vrad,
double  z 
)

Applies an incremental redshift correction to a radial velocity. For example, you may use this function to correct a radial velocity calculated by rad_vel() or rad_vel2() for a Solar-system body to account for the gravitational redshift for light originating at a specific distance away from the body. For the Sun, you may want to undo the redshift correction applied for the photosphere using unredshift_vrad() first.

Parameters
vrad[km/s] Radial velocity
zRedshift correction to apply
Returns
[km/s] The redshift corrected radial velocity or NAN if the redshift value is invalid (errno will be set to EINVAL).
See also
unredshift_vrad()
grav_redshift()
novas_z_add()
Since
1.2
Author
Attila Kovacs

References novas_v2z(), and novas_z2v().

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

NOTES:

  1. The accuracy of the output CIRS coordinates depends on how the input TOD coordinates were obtained. If TOD was calculated via the old (pre IAU 2006) method, using the Lieske et al. 1997 nutation model, then the limited accuracy of that model will affect the resulting coordinates. This is the case for the SuperNOVAS functions novas_geom_posvel() and novas_sky_pos() also, when called with NOVAS_TOD as the system, as well as all legacy NOVAS C functions that produce TOD coordinates.
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.

◆ unredshift_vrad()

double unredshift_vrad ( double  vrad,
double  z 
)

Undoes an incremental redshift correction that was applied to radial velocity.

Parameters
vrad[km/s] Radial velocity
zRedshift correction to apply
Returns
[km/s] The radial velocity without the redshift correction or NAN if the redshift value is invalid. (errno will be set to EINVAL)
See also
redshift_vrad()
grav_redshift()
Since
1.2
Author
Attila Kovacs

References novas_v2z(), and novas_z2v().