![]() |
SuperNOVAS v1.3
The NOVAS C library, made better
|
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) |
SuperNOVAS only functions, which are not integral to the functionality of novas.c, and thus can live in a separate, more manageably sized, module.
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.
jd_tt | [day] Terrestrial Time (TT) based Julian date |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) |
ra | [h] the apparent R.A. coordinate measured from the true equinox of date. |
References cio_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.
jd_tt | [day] Terrestrial Time (TT) based Julian date |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) |
ra | [h] The CIRS right ascension coordinate, measured from the CIO. |
References cio_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 | ||
) |
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:
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 | |
accuracy | NOVAS_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. | |
in | Position vector, geocentric equatorial rectangular coordinates, referred to CIRS axes (celestial system). | |
[out] | out | Position vector, geocentric equatorial rectangular coordinates, referred to ITRS axes (terrestrial system). |
References cel2ter(), EROT_ERA, and NOVAS_DYNAMICAL_CLASS.
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.
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 | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
in | CIRS Input (x, y, z) position or velocity vector | |
[out] | out | Output position or velocity 3-vector in the True of Date (TOD) frame. It can be the same vector as the input. |
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'.
jd_tt | [day] Terrestrial Time (TT) based Julian date. (Unused if 'coord_sys' is NOVAS_GCRS_EQUATOR[2]) | |
coord_sys | The astrometric reference system of the coordinates. If 'coord_sys' is NOVAS_GCRS_EQUATOR(2), the input GCRS coordinates are converted to J2000 ecliptic coordinates. | |
accuracy | NOVAS_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. |
References ecl2equ_vec().
int gal2equ | ( | double | glon, |
double | glat, | ||
double *restrict | ra, | ||
double *restrict | dec | ||
) |
Converts galactic longitude and latitude to ICRS right ascension and declination.
REFERENCES:
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. |
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:
leap_seconds | [s] Leap seconds at the time of observations |
dut1 | [s] UT1 - UTC time difference [-0.5:0.5] |
ut1_to_tt
argument.References get_utc_to_tt().
double get_utc_to_tt | ( | int | leap_seconds | ) |
Returns the difference between Terrestrial Time (TT) and Universal Coordinated Time (UTC)
leap_seconds | [s] The current leap seconds (see IERS Bulletins) |
References NOVAS_TAI_TO_TT.
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.
M_kg | [kg] Mass of gravitating body that is contained inside the emitting radius. |
r_m | [m] Radius at which light is emitted. |
References C.
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:
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date | |
accuracy | NOVAS_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. |
References grav_bodies_full_accuracy, grav_bodies_reduced_accuracy, grav_undo_planets(), NOVAS_FULL_ACCURACY, and obs_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:
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. | |
planets | Apparent 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. |
References grav_planets(), and novas_inv_max_iter.
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.
location | Observer 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] | itrs | Unit 3-vector direction in Earth-fixed ITRS frame |
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:
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 | |
accuracy | NOVAS_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. | |
in | Position vector, geocentric equatorial rectangular coordinates, referred to ITRS axes (terrestrial system) | |
[out] | out | Position vector, geocentric equatorial rectangular coordinates, referred to CIRS axes (celestial system). |
References EROT_ERA, NOVAS_DYNAMICAL_CLASS, and ter2cel().
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.
location | Observer location on Earth | |
itrs | 3-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. |
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:
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 | |
accuracy | NOVAS_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. | |
in | Position vector, geocentric equatorial rectangular coordinates, referred to ITRS axes (terrestrial system) | |
[out] | out | Position vector, geocentric equatorial rectangular coordinates, referred to True of Date (TOD) axes (celestial system) |
References EROT_GST, NOVAS_DYNAMICAL_CLASS, and ter2cel().
int j2000_to_gcrs | ( | const double * | in, |
double * | out | ||
) |
Change J2000 coordinates to GCRS coordinates. Same as frame_tie() called with J2000_TO_ICRS
in | J2000 input 3-vector | |
[out] | out | GCRS output 3-vector |
References frame_tie(), and J2000_TO_ICRS.
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.
location | Current longitude, latitude and altitude, and local weather (temperature and pressure) | |
vel | [km/s] Surface velocity. | |
[out] | obs | Pointer to data structure to populate. |
References IN_SPACE_INIT, make_observer(), NOVAS_AIRBORNE_OBSERVER, and in_space::sc_vel.
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.
star | Pointer to structure to populate with the catalog data for a celestial object located outside the solar system, specified with ICRS coordinates. | |
[out] | source | Pointer to the celestial object data structure to be populated. |
References make_object(), NOVAS_CATALOG_OBJECT, cat_entry::starname, and cat_entry::starnumber.
Populates and object data structure with the data for a catalog source for a given system of catalog coordinates.
star | Pointer to structure to populate with the catalog data for a celestial object located outside the solar system. | |
system | Input 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] | source | Pointer to the celestial object data structure to be populated with the corresponding ICRS catalog coordinates, after appying proper-motion and precession corrections as appropriate. |
References make_cat_object(), and object::star.
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.
name | Name 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. | |
num | Solar-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] | body | Pointer to structure to populate. |
References make_object(), and NOVAS_EPHEM_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.
name | Name of object. It may be NULL if not relevant. | |
num | Solar-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. | |
orbit | The orbital parameters to adopt. The data will be copied, not referenced. | |
[out] | body | Pointer to structure to populate. |
References make_object(), NOVAS_ORBITAL_OBJECT, and object::orbit.
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.
name | Object 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). | |
z | Redhift value (λobs / λrest - 1 = frest / fobs - 1). | |
[out] | source | Pointer to structure to populate. |
References make_cat_entry(), and novas_z2v().
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.
name | Object 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). | |
z | Redhift value (λobs / λrest - 1 = frest / fobs - 1). | |
[out] | source | Pointer to structure to populate. |
References make_cat_object(), and make_redshifted_cat_entry().
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.
name | Object 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). | |
system | Input 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. | |
z | Redhift value (λobs / λrest - 1 = frest / fobs - 1). | |
[out] | source | Pointer to the celestial object data structure to be populated with the corresponding ICRS catalog coordinates. |
References make_redshifted_object(), and object::star.
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.
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] | obs | Pointer to the data structure to populate |
References make_in_space(), make_observer(), and NOVAS_SOLAR_SYSTEM_OBSERVER.
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:
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. |
References novas_h2e_offset().
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.
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 |
double novas_epoch | ( | const char *restrict | system | ) |
Returns the Julian day corresponding to an astronomical coordinate epoch.
system | Coordinate 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. |
References NOVAS_JD_B1950, NOVAS_JD_HIP, NOVAS_JD_J2000, NOVAS_SYSTEM_FK4, NOVAS_SYSTEM_FK5, NOVAS_SYSTEM_HIP, and NOVAS_SYSTEM_ICRS.
double novas_equ_sep | ( | double | ra1, |
double | dec1, | ||
double | ra2, | ||
double | dec2 | ||
) |
Returns the angular separation of two equatorial locations on a sphere.
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 |
References novas_sep().
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:
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. |
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.
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. | |
source | Observed Solar-system source | |
[out] | rate | [AU/day] (optional) Returned rate of recession from Sun |
References ephemeris(), NOVAS_CATALOG_OBJECT, NOVAS_HELIOCENTER, and NOVAS_REDUCED_ACCURACY.
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.
az | [deg] Azimuth angle |
el | [deg] Elevation angle |
lat | [deg] Geodetic latitude of observer |
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:
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. |
References NOVAS_JD_J2000, precession(), and radec2vector().
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.
name | The 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. |
References NOVAS_PLANET_NAMES_INIT, NOVAS_PLANETS, and NOVAS_SSB.
double novas_sep | ( | double | lon1, |
double | lat1, | ||
double | lon2, | ||
double | lat2 | ||
) |
Returns the angular separation of two locations on a sphere.
lon1 | [deg] longitude of first location |
lat1 | [deg] latitude of first location |
lon2 | [deg] longitude of second location |
lat2 | [deg] latitude of second location |
int novas_set_orbsys_pole | ( | enum novas_reference_system | type, |
double | ra, | ||
double | dec, | ||
novas_orbital_system *restrict | sys | ||
) |
Sets the orientation of an orbital system using the RA and DEC coordinates of the pole of the Laplace (or else equatorial) plane relative to which the orbital elements are defined. Orbital parameters of planetary satellites normally include the R.A. and declination of the pole of the local Laplace plane in which the Keplerian orbital elements are referenced.
The system will become referenced to the equatorial plane, the relative obliquity is set to (90° - dec
), while the argument of the ascending node ('Omega') is set to (90° + ra
).
NOTES:
type | Coordinate reference system in which ra and dec are defined (e.g. NOVAS_GCRS). | |
ra | [h] the R.A. of the pole of the oribtal reference plane. | |
dec | [deg] the declination of the pole of the oribtal reference plane. | |
[out] | sys | Orbital system |
sys
pointer is NULL.References NOVAS_EQUATORIAL_PLANE.
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.
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. |
source | Observed Solar-system source |
References novas_helio_dist(), and NOVAS_SOLAR_CONSTANT.
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:
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. |
References NOVAS_JD_J2000, precession(), and radec2vector().
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.
vel | [km/s] velocity (i.e. rate) of recession. |
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).
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. |
References novas_xyz_to_los().
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).
z1 | One of the redshift values |
z2 | The other redshift value |
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.
z | A redhift value |
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.
jd_tt | [day] Terrestrial Time (TT) based Julian date of observation. | |
source | Catalog source or solar_system body. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos | Structure to populate with the calculated CIRS position data |
References NOVAS_CIRS, and place().
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.
jd_tt | [day] Terrestrial Time (TT) based Julian date of observation. | |
source | Catalog source or solar_system body. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos | Structure to populate with the calculated GCRS position data |
References NOVAS_GCRS, and place().
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.
jd_tt | [day] Terrestrial Time (TT) based Julian date of observation. | |
source | Catalog source or solar_system body. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos | Structure to populate with the calculated geocentric ICRS position data (Unlike place_gcrs(), the calculated coordinates do not account for aberration or gravitational deflection). |
References NOVAS_ICRS, and place().
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.
jd_tt | [day] Terrestrial Time (TT) based Julian date of observation. | |
source | Catalog source or solar_system body. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos | Structure to populate with the calculated CIRS position data |
References NOVAS_J2000, and place().
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.
jd_tt | [day] Terrestrial Time (TT) based Julian date of observation. | |
source | Catalog source or solar_system body. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos | Structure to populate with the calculated CIRS position data |
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.
jd_tt | [day] Terrestrial Time (TT) based Julian date of observation. | |
source | Catalog source or solar_system body. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos | Structure to populate with the calculated CIRS position data |
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.
vrad | [km/s] Radial velocity |
z | Redshift correction to apply |
References novas_v2z(), and novas_z2v().
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:
NOVAS_TOD
as the system, as well as all legacy NOVAS C functions that produce TOD coordinates. 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 | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
in | CIRS Input (x, y, z) position or velocity vector | |
[out] | out | Output position or velocity 3-vector in the True of Date (TOD) frame. It can be the same vector as the input. |
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:
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. | |
accuracy | NOVAS_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. | |
in | Position vector, geocentric equatorial rectangular coordinates, referred to True of Date (TOD) axes (celestial system). | |
[out] | out | Position vector, geocentric equatorial rectangular coordinates, referred to ITRS axes (terrestrial system). |
References cel2ter(), EROT_GST, and NOVAS_DYNAMICAL_CLASS.
double unredshift_vrad | ( | double | vrad, |
double | z | ||
) |
Undoes an incremental redshift correction that was applied to radial velocity.
vrad | [km/s] Radial velocity |
z | Redshift correction to apply |
References novas_v2z(), and novas_z2v().