SuperNOVAS v1.2
The NOVAS C library, made better
|
Functions | |
int | aberration (const double *pos, const double *vobs, double lighttime, double *out) |
double | accum_prec (double t) |
short | app_planet (double jd_tt, const object *ss_body, enum novas_accuracy accuracy, double *ra, double *dec, double *dis) |
short | app_star (double jd_tt, const cat_entry *star, enum novas_accuracy accuracy, double *ra, double *dec) |
short | astro_planet (double jd_tt, const object *ss_body, enum novas_accuracy accuracy, double *ra, double *dec, double *dis) |
short | astro_star (double jd_tt, const cat_entry *star, enum novas_accuracy accuracy, double *ra, double *dec) |
int | bary2obs (const double *pos, const double *pos_obs, double *out, double *lighttime) |
int | cal_date (double tjd, short *year, short *month, short *day, double *hour) |
short | cel2ter (double jd_ut1_high, double jd_ut1_low, double ut1_to_tt, enum novas_earth_rotation_measure erot, enum novas_accuracy accuracy, enum novas_equatorial_class class, double xp, double yp, const double *in, double *out) |
short | cel_pole (double jd_tt, enum novas_pole_offset_type type, double dpole1, double dpole2) |
short | cio_array (double jd_tdb, long n_pts, ra_of_cio *cio) |
short | cio_basis (double jd_tdb, double ra_cio, enum novas_cio_location_type loc_type, enum novas_accuracy accuracy, double *x, double *y, double *z) |
short | cio_location (double jd_tdb, enum novas_accuracy accuracy, double *ra_cio, short *loc_type) |
short | cio_ra (double jd_tt, enum novas_accuracy accuracy, double *ra_cio) |
int | cirs_to_gcrs (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
double | d_light (const double *pos_src, const double *pos_body) |
int | e_tilt (double jd_tdb, enum novas_accuracy accuracy, double *mobl, double *tobl, double *ee, double *dpsi, double *deps) |
short | ecl2equ_vec (double jd_tt, enum novas_equator_type coord_sys, enum novas_accuracy accuracy, const double *in, double *out) |
double | ee_ct (double jd_tt_high, double jd_tt_low, enum novas_accuracy accuracy) |
short | ephemeris (const double *jd_tdb, const object *body, enum novas_origin origin, enum novas_accuracy accuracy, double *pos, double *vel) |
short | equ2ecl (double jd_tt, enum novas_equator_type coord_sys, enum novas_accuracy accuracy, double ra, double dec, double *elon, double *elat) |
short | equ2ecl_vec (double jd_tt, enum novas_equator_type coord_sys, enum novas_accuracy accuracy, const double *in, double *out) |
int | equ2gal (double ra, double dec, double *glon, double *glat) |
int | equ2hor (double jd_ut1, double ut1_to_tt, enum novas_accuracy accuracy, double xp, double yp, const on_surface *location, double ra, double dec, enum novas_refraction_model ref_option, double *zd, double *az, double *rar, double *decr) |
double | era (double jd_ut1_high, double jd_ut1_low) |
int | frame_tie (const double *in, enum novas_frametie_direction direction, double *out) |
int | fund_args (double t, novas_delaunay_args *a) |
short | gcrs2equ (double jd_tt, enum novas_dynamical_type sys, enum novas_accuracy accuracy, double rag, double decg, double *ra, double *dec) |
int | gcrs_to_cirs (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
int | gcrs_to_j2000 (const double *in, double *out) |
int | gcrs_to_mod (double jd_tdb, const double *in, double *out) |
int | gcrs_to_tod (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
short | geo_posvel (double jd_tt, double ut1_to_tt, enum novas_accuracy accuracy, const observer *obs, double *pos, double *vel) |
novas_ephem_provider | get_ephem_provider () |
novas_planet_provider | get_planet_provider () |
novas_planet_provider_hp | get_planet_provider_hp () |
short | grav_def (double jd_tdb, enum novas_observer_place unused, enum novas_accuracy accuracy, const double *pos_src, const double *pos_obs, double *out) |
int | grav_planets (const double *pos_src, const double *pos_obs, const novas_planet_bundle *planets, double *out) |
int | grav_vec (const double *pos_src, const double *pos_obs, const double *pos_body, double rmass, double *out) |
double | ira_equinox (double jd_tdb, enum novas_equinox_type equinox, enum novas_accuracy accuracy) |
int | j2000_to_tod (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
double | julian_date (short year, short month, short day, double hour) |
short | light_time (double jd_tdb, const object *body, const double *pos_obs, double tlight0, enum novas_accuracy accuracy, double *pos_src_obs, double *tlight) |
int | light_time2 (double jd_tdb, const object *body, const double *pos_obs, double tlight0, enum novas_accuracy accuracy, double *p_src_obs, double *v_ssb, double *tlight) |
int | limb_angle (const double *pos_src, const double *pos_obs, double *limb_ang, double *nadir_ang) |
short | local_planet (double jd_tt, const object *ss_body, double ut1_to_tt, const on_surface *position, enum novas_accuracy accuracy, double *ra, double *dec, double *dis) |
short | local_star (double jd_tt, double ut1_to_tt, const cat_entry *star, const on_surface *position, enum novas_accuracy accuracy, double *ra, double *dec) |
short | make_cat_entry (const char *star_name, const char *catalog, long cat_num, double ra, double dec, double pm_ra, double pm_dec, double parallax, double rad_vel, cat_entry *star) |
int | make_in_space (const double *sc_pos, const double *sc_vel, in_space *loc) |
short | make_object (enum novas_object_type type, long number, const char *name, const cat_entry *star, object *source) |
short | make_observer (enum novas_observer_place where, const on_surface *loc_surface, const in_space *loc_space, observer *obs) |
int | make_observer_at_geocenter (observer *obs) |
int | make_observer_in_space (const double *sc_pos, const double *sc_vel, observer *obs) |
int | make_observer_on_surface (double latitude, double longitude, double height, double temperature, double pressure, observer *obs) |
int | make_on_surface (double latitude, double longitude, double height, double temperature, double pressure, on_surface *loc) |
int | make_planet (enum novas_planet num, object *planet) |
double | mean_obliq (double jd_tdb) |
short | mean_star (double jd_tt, double tra, double tdec, enum novas_accuracy accuracy, double *ira, double *idec) |
int | mod_to_gcrs (double jd_tdb, const double *in, double *out) |
double | norm_ang (double angle) |
void | novas_case_sensitive (int value) |
void | novas_debug (enum novas_debug_mode mode) |
enum novas_debug_mode | novas_get_debug_mode () |
int | novas_orbit_posvel (double jd_tdb, const novas_orbital *orbit, enum novas_accuracy accuracy, double *pos, double *vel) |
double | novas_z2v (double z) |
int | nutation (double jd_tdb, enum novas_nutation_direction direction, enum novas_accuracy accuracy, const double *in, double *out) |
int | nutation_angles (double t, enum novas_accuracy accuracy, double *dpsi, double *deps) |
int | obs_planets (double jd_tdb, enum novas_accuracy accuracy, const double *pos_obs, int pl_mask, novas_planet_bundle *planets) |
int | obs_posvel (double jd_tdb, double ut1_to_tt, enum novas_accuracy accuracy, const observer *obs, const double *geo_pos, const double *geo_vel, double *pos, double *vel) |
short | place (double jd_tt, const object *source, const observer *location, double ut1_to_tt, enum novas_reference_system coord_sys, enum novas_accuracy accuracy, sky_pos *output) |
int | place_star (double jd_tt, const cat_entry *star, const observer *obs, double ut1_to_tt, enum novas_reference_system system, enum novas_accuracy accuracy, sky_pos *pos) |
double | planet_lon (double t, enum novas_planet planet) |
short | precession (double jd_tdb_in, const double *in, double jd_tdb_out, double *out) |
int | proper_motion (double jd_tdb_in, const double *pos, const double *vel, double jd_tdb_out, double *out) |
int | rad_vel (const object *source, const double *pos_src, const double *vel_src, const double *vel_obs, double d_obs_geo, double d_obs_sun, double d_src_sun, double *rv) |
double | rad_vel2 (const object *source, const double *pos_emit, const double *vel_src, const double *pos_det, const double *vel_obs, double d_obs_geo, double d_obs_sun, double d_src_sun) |
int | radec2vector (double ra, double dec, double dist, double *pos) |
int | radec_planet (double jd_tt, const object *ss_body, const observer *obs, double ut1_to_tt, enum novas_reference_system sys, enum novas_accuracy accuracy, double *ra, double *dec, double *dis, double *rv) |
int | radec_star (double jd_tt, const cat_entry *star, const observer *obs, double ut1_to_tt, enum novas_reference_system sys, enum novas_accuracy accuracy, double *ra, double *dec, double *rv) |
double | refract (const on_surface *location, enum novas_refraction_model option, double zd_obs) |
double | refract_astro (const on_surface *location, enum novas_refraction_model option, double zd_astro) |
int | set_cio_locator_file (const char *filename) |
int | set_ephem_provider (novas_ephem_provider func) |
int | set_nutation_lp_provider (novas_nutation_provider func) |
int | set_planet_provider (novas_planet_provider func) |
int | set_planet_provider_hp (novas_planet_provider_hp func) |
short | sidereal_time (double jd_ut1_high, double jd_ut1_low, double ut1_to_tt, enum novas_equinox_type gst_type, enum novas_earth_rotation_measure erot, enum novas_accuracy accuracy, double *gst) |
int | spin (double angle, const double *in, double *out) |
int | starvectors (const cat_entry *star, double *pos, double *vel) |
int | tdb2tt (double jd_tdb, double *jd_tt, double *secdiff) |
short | ter2cel (double jd_ut1_high, double jd_ut1_low, double ut1_to_tt, enum novas_earth_rotation_measure erot, enum novas_accuracy accuracy, enum novas_equatorial_class class, double xp, double yp, const double *in, double *out) |
int | terra (const on_surface *location, double lst, double *pos, double *vel) |
int | tod_to_gcrs (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
int | tod_to_j2000 (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
short | topo_planet (double jd_tt, const object *ss_body, double ut1_to_tt, const on_surface *position, enum novas_accuracy accuracy, double *ra, double *dec, double *dis) |
short | topo_star (double jd_tt, double ut1_to_tt, const cat_entry *star, const on_surface *position, enum novas_accuracy accuracy, double *ra, double *dec) |
short | transform_cat (enum novas_transform_type option, double jd_tt_in, const cat_entry *in, double jd_tt_out, const char *out_id, cat_entry *out) |
int | transform_hip (const cat_entry *hipparcos, cat_entry *hip_2000) |
double | tt2tdb (double jd_tt) |
short | vector2radec (const double *pos, double *ra, double *dec) |
short | virtual_planet (double jd_tt, const object *ss_body, enum novas_accuracy accuracy, double *ra, double *dec, double *dis) |
short | virtual_star (double jd_tt, const cat_entry *star, enum novas_accuracy accuracy, double *ra, double *dec) |
int | wobble (double jd_tt, enum novas_wobble_direction direction, double xp, double yp, const double *in, double *out) |
Variables | |
double | EPS_COR = 0.0 |
int | grav_bodies_full_accuracy = DEFAULT_GRAV_BODIES_FULL_ACCURACY |
int | grav_bodies_reduced_accuracy = DEFAULT_GRAV_BODIES_REDUCED_ACCURACY |
int | novas_inv_max_iter = 100 |
double | PSI_COR = 0.0 |
SuperNOVAS astrometry software based on the Naval Observatory Vector Astrometry Software (NOVAS). It has been modified to fix outstanding issues and to make it easier to use.
Based on the NOVAS C Edition, Version 3.1:
U. S. Naval Observatory
Astronomical Applications Dept.
Washington, DC
http://www.usno.navy.mil/USNO/astronomical-applications
int aberration | ( | const double * | pos, |
const double * | vobs, | ||
double | lighttime, | ||
double * | out | ||
) |
Corrects position vector for aberration of light. Algorithm includes relativistic terms.
NOTES:
REFERENCES:
pos | [AU] Position vector of source relative to observer | |
vobs | [AU/day] Velocity vector of observer, relative to the solar system barycenter, components in AU/day. | |
lighttime | [day] Light time from object to Earth in days (if known). Or set to 0, and this function will compute it. | |
[out] | out | [AU] Position vector, referred to origin at center of mass of the Earth, corrected for aberration, components in AU. It can be the same vector as one of the inputs. |
References C_AUDAY.
double accum_prec | ( | double | t | ) |
Returns the general precession in longitude (Simon et al. 1994), equivalent to 5028.8200 arcsec/cy at J2000.
t | [cy] Julian centuries since J2000 |
References TWOPI.
short app_planet | ( | double | jd_tt, |
const object * | ss_body, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec, | ||
double * | dis | ||
) |
Computes the apparent place of a solar system body. This is the same as calling place() for the body with NOVAS_TOD as the system, except the different set of return values used.
REFERENCES:
jd_tt | [day] Terretrial Time (TT) based Julian date. | |
ss_body | Pointer to structure containing the body designation for the solar system body. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Apparent right ascension in hours, referred to true equator and equinox of date 'jd_tt'. (It may be NULL if not required) |
[out] | dec | [deg] Apparent declination in degrees, referred to true equator and equinox of date 'jd_tt'. (It may be NULL if not required) |
[out] | dis | [AU] True distance from Earth to the body at 'jd_tt' in AU (can be NULL if not needed). |
References NOVAS_TOD, and radec_planet().
short app_star | ( | double | jd_tt, |
const cat_entry * | star, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec | ||
) |
Computes the apparent place of a star, referenced to dynamical equator at date 'jd_tt', given its catalog mean place, proper motion, parallax, and radial velocity.
Notwithstanding the different set of return values, this is the same as calling place_star() with a NULL observer location and NOVAS_TOD as the system for an object that specifies the star.
REFERENCES:
jd_tt | [day] Terretrial Time (TT) based Julian date. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Apparent right ascension in hours, referred to true equator and equinox of date 'jd_tt' (it may be NULL if not required). |
[out] | dec | [deg] Apparent declination in degrees, referred to true equator and equinox of date 'jd_tt' (it may be NULL if not required). |
References NOVAS_TOD, and radec_star().
short astro_planet | ( | double | jd_tt, |
const object * | ss_body, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec, | ||
double * | dis | ||
) |
Computes the astrometric place of a solar system body, referenced to the ICRS without light deflection or aberration. This is the same as calling place_icrs() for the body, except the different set of return values used.
REFERENCES:
jd_tt | [day] Terretrial Time (TT) based Julian date. | |
ss_body | Pointer to structure containing the body designation for the solar system body. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Astrometric right ascension in hours, referred to the ICRS, without light deflection or aberration. (It may be NULL if not required) |
[out] | dec | [deg] Astrometric declination in degrees, referred to the ICRS, without light deflection or aberration. (It may be NULL if not required) |
[out] | dis | [AU] True distance from Earth to the body at 'jd_tt' in AU (can be NULL if not needed). |
References NOVAS_ICRS, and radec_planet().
short astro_star | ( | double | jd_tt, |
const cat_entry * | star, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec | ||
) |
Computes the astrometric place of a star, referred to the ICRS without light deflection or aberration, at date 'jd_tt', given its catalog mean place, proper motion, parallax, and radial velocity.
Notwithstanding the different set of return values, this is the same as calling place_star() with a NULL observer location and NOVAS_ICRS as the system, or place_icrs() for an object that specifies the star.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Astrometric right ascension in hours, referred to the ICRS, without light deflection or aberration. (It may be NULL if not required) |
[out] | dec | [deg] Astrometric declination in degrees, referred to the ICRS, without light deflection or aberration. (It may be NULL if not required) |
References NOVAS_ICRS, and radec_star().
int bary2obs | ( | const double * | pos, |
const double * | pos_obs, | ||
double * | out, | ||
double * | lighttime | ||
) |
Moves the origin of coordinates from the barycenter of the solar system to the observer (or the geocenter); i.e., this function accounts for parallax (annual+geocentric or just annual).
REFERENCES:
pos | [AU] Position vector, referred to origin at solar system barycenter, components in AU. | |
pos_obs | [AU] Position vector of observer (or the geocenter), with respect to origin at solar system barycenter, components in AU. | |
[out] | out | [AU] Position vector, referred to origin at center of mass of the Earth, components in AU. It may be NULL if not required, or be the same vector as either of the inputs. |
[out] | lighttime | [day] Light time from object to Earth in days. It may be NULL if not required. |
References C_AUDAY.
int cal_date | ( | double | tjd, |
short * | year, | ||
short * | month, | ||
short * | day, | ||
double * | hour | ||
) |
This function will compute a broken down date on the Gregorian calendar for given the Julian date input. Input Julian date can be based on any UT-like time scale (UTC, UT1, TT, etc.) - output time value will have same basis.
REFERENCES:
tjd | [day] Julian date | |
[out] | year | [yr] Gregorian calendar year. It may be NULL if not required. |
[out] | month | [month] Gregorian calendat month [1:12]. It may be NULL if not required. |
[out] | day | [day] Day of the month [1:31]. It may be NULL if not required. |
[out] | hour | [h] Hour of day [0:24]. It may be NULL if not required. |
short cel2ter | ( | double | jd_ut1_high, |
double | jd_ut1_low, | ||
double | ut1_to_tt, | ||
enum novas_earth_rotation_measure | erot, | ||
enum novas_accuracy | accuracy, | ||
enum novas_equatorial_class | class, | ||
double | xp, | ||
double | yp, | ||
const double * | in, | ||
double * | out | ||
) |
Rotates a vector from the celestial to the terrestrial system. Specifically, it transforms a vector in the GCRS, or the dynamcal CIRS or TOD frames to the ITRS (a rotating earth-fixed system) by applying rotations for the GCRS-to-dynamical frame tie, precession, nutation, Earth rotation, and polar motion.
If 'system' is NOVAS_CIRS then method EROT_ERA must be used. Similarly, if 'system' is NOVAS_TOD then method must be EROT_ERA. Otherwise an error 3 is returned.
If both 'xp' and 'yp' are set to 0 no polar motion is included in the transformation.
REFERENCES:
jd_ut1_high | [day] High-order part of UT1 Julian date. | |
jd_ut1_low | [day] Low-order part of UT1 Julian date. | |
ut1_to_tt | [s] TT - UT1 Time difference in seconds | |
erot | EROT_ERA (0) or EROT_GST (1), depending on whether to use GST relative to equinox of date (pre IAU 2006) or ERA relative to the CIO (IAU 2006 standard) as the Earth rotation measure. The main effect of this option is that it specifies the input coordinate system as CIRS or TOD when the input coordinate class is NOVAS_DYNAMICAL_CLASS. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
class | Input coordinate class, NOVAS_REFERENCE_CLASS (0) or NOVAS_DYNAMICAL_CLASS (1). Use the former if the input coordinates are in the GCRS, and the latter if they are CIRS or TOD (the 'erot' parameter selects which dynamical system the input is specified in.) | |
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 | Input position vector, geocentric equatorial rectangular coordinates in the specified input coordinate system (GCRS if 'class' is NOVAS_REFERENCE_CLASS; or else either CIRS if 'erot' is EROT_ERA, or TOD if 'erot' is EROT_GST). | |
[out] | out | ITRS position vector, geocentric equatorial rectangular coordinates (terrestrial system). It can be the same vector as the input. |
References era(), EROT_ERA, EROT_GST, gcrs_to_cirs(), gcrs_to_tod(), NOVAS_DYNAMICAL_CLASS, NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, NOVAS_TRUE_EQUINOX, sidereal_time(), spin(), tt2tdb(), wobble(), and WOBBLE_PEF_TO_ITRS.
short cel_pole | ( | double | jd_tt, |
enum novas_pole_offset_type | type, | ||
double | dpole1, | ||
double | dpole2 | ||
) |
specifies the celestial pole offsets for high-precision applications. Each set of offsets is a correction to the modeled position of the pole for a specific date, derived from observations and published by the IERS.
The variables 'PSI_COR' and 'EPS_COR' are used only in NOVAS function e_tilt().
This function, if used, should be called before any other NOVAS functions for a given date. Values of the pole offsets specified via a call to this function will be used until explicitly changed.
'tjd' is used only if 'type' is POLE_OFFSETS_X_Y (2), to transform dx and dy to the equivalent Δδψ and Δδε values.
If 'type' is POLE_OFFSETS_X_Y (2), dx and dy are unit vector component corrections, but are expressed in milliarcseconds simply by multiplying by 206264806, the number of milliarcseconds in one radian.
NOTES:
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. |
type | POLE_OFFSETS_DPSI_DEPS (1) or POLE_OFFSETS_X_Y (2) |
dpole1 | [mas] Value of celestial pole offset in first coordinate, (Δδψ or dx) in milliarcseconds. |
dpole2 | [mas] Value of celestial pole offset in second coordinate, (Δδε or dy) in milliarcseconds. |
References EPS_COR, POLE_OFFSETS_DPSI_DEPS, POLE_OFFSETS_X_Y, and PSI_COR.
short cio_array | ( | double | jd_tdb, |
long | n_pts, | ||
ra_of_cio * | cio | ||
) |
Given an input TDB Julian date and the number of data points desired, this function returns a set of Julian dates and corresponding values of the GCRS right ascension of the celestial intermediate origin (CIO). The range of dates is centered (at least approximately) on the requested date. The function obtains the data from an external data file.
This function assumes that a CIO locator file (CIO_RA.TXT
or cio_ra.bin
) exists in the default location (configured at build time), or else was specified via set_cio_locator_file()
prior to calling this function.
NOTES:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
n_pts | Number of Julian dates and right ascension values requested (not less than 2 or more than NOVAS_CIO_CACHE_SIZE). | |
[out] | cio | A time series (array) of the right ascension of the Celestial Intermediate Origin (CIO) with respect to the GCRS. |
References DEFAULT_CIO_LOCATOR_FILE, NOVAS_CIO_CACHE_SIZE, and set_cio_locator_file().
short cio_basis | ( | double | jd_tdb, |
double | ra_cio, | ||
enum novas_cio_location_type | loc_type, | ||
enum novas_accuracy | accuracy, | ||
double * | x, | ||
double * | y, | ||
double * | z | ||
) |
Computes the orthonormal basis vectors, with respect to the GCRS (geocentric ICRS), of the celestial intermediate system defined by the celestial intermediate pole (CIP) (in the z direction) and the celestial intermediate origin (CIO) (in the x direction). A TDB Julian date and the right ascension of the CIO at that date is required as input. The right ascension of the CIO can be with respect to either the GCRS origin or the true equinox of date – different algorithms are used in the two cases.
This function effectively constructs the matrix C in eq. (3) of the reference.
REFERENCES:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
ra_cio | [h] Right ascension of the CIO at epoch (hours). | |
loc_type | CIO_VS_GCRS (1) if the cio location is relative to the GCRS or else CIO_VS_EQUINOX (2) if relative to the true equinox of date. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | x | Unit 3-vector toward the CIO, equatorial rectangular coordinates, referred to the GCRS. |
[out] | y | Unit 3-vector toward the y-direction, equatorial rectangular coordinates, referred to the GCRS. |
[out] | z | Unit 3-vector toward north celestial pole (CIP), equatorial rectangular coordinates, referred to the GCRS. |
References CIO_VS_EQUINOX, CIO_VS_GCRS, NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, and tod_to_gcrs().
short cio_location | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
double * | ra_cio, | ||
short * | loc_type | ||
) |
Returns the location of the celestial intermediate origin (CIO) for a given Julian date, as a right ascension with respect to either the GCRS (geocentric ICRS) origin or the true equinox of date. The CIO is always located on the true equator (= intermediate equator) of date.
The user may specify an interpolation file to use via set_cio_locator_file() prior to calling this function. In that case the call will return CIO location relative to GCRS. In the absence of the table, it will calculate the CIO location relative to the true equinox. In either case the type of the location is returned alongside the corresponding CIO location value.
NOTES:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra_cio | [h] Right ascension of the CIO, in hours, or NAN if returning with an error. |
[out] | loc_type | Pointer in which to return the reference system in which right ascension is given, which is either CIO_VS_GCRS (1) if the location was obtained via interpolation of the available data file, or else CIO_VS_EQUINOX (2) if it was calculated locally. It is set to -1 if returning with an error. |
References cio_array(), CIO_VS_EQUINOX, CIO_VS_GCRS, ira_equinox(), novas_debug(), NOVAS_DEBUG_OFF, NOVAS_DEBUG_ON, NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, and NOVAS_TRUE_EQUINOX.
short cio_ra | ( | double | jd_tt, |
enum novas_accuracy | accuracy, | ||
double * | ra_cio | ||
) |
Computes the true right ascension of the celestial intermediate origin (CIO) at a given TT Julian date. This is the negative value for the equation of the origins.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra_cio | [h] Right ascension of the CIO, with respect to the true equinox of date, in hours (+ or -), or NAN when returning with an error code. |
References cio_basis(), cio_location(), NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, tod_to_gcrs(), and tt2tdb().
int cirs_to_gcrs | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
const double * | in, | ||
double * | out | ||
) |
Transforms a rectangular equatorial (x, y, z) vector from the Celestial Intermediate Reference System (CIRS) frame at the given epoch to the Geocentric Celestial Reference System (GCRS).
jd_tdb | [day] Barycentric Dynamical Time (TDB) 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 GCRS coordinate frame. It can be the same vector as the input. |
References cio_basis(), and cio_location().
double d_light | ( | const double * | pos_src, |
const double * | pos_body | ||
) |
Returns the difference in light-time, for a star, between the barycenter of the solar system and the observer (or the geocenter) (Usage A).
Alternatively (Usage B), this function returns the light-time from the observer (or the geocenter) to a point on a light ray that is closest to a specific solar system body. For this purpose, 'pos_src' is the position vector toward observed object, with respect to origin at observer (or the geocenter); 'pos_body' is the position vector of solar system body, with respect to origin at observer (or the geocenter), components in AU; and the returned value is the light time to point on line defined by 'pos' that is closest to solar system body (positive if light passes body before hitting observer, i.e., if 'pos_body' is within 90 degrees of 'pos_src').
NOTES:
pos_src | Position vector towards observed object, with respect to the SSB (Usage A), or relative to the observer / geocenter (Usage B). |
pos_body | [AU] Position of observer relative to SSB (Usage A), or position of intermediate solar-system body with respect to the observer / geocenter (Usage B). |
References C_AUDAY.
int e_tilt | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
double * | mobl, | ||
double * | tobl, | ||
double * | ee, | ||
double * | dpsi, | ||
double * | deps | ||
) |
Computes quantities related to the orientation of the Earth's rotation axis at Julian date 'jd_tdb'.
Values of the celestial pole offsets 'PSI_COR' and 'EPS_COR' are set using function 'cel_pole', if desired. See the prolog of cel_pole() for details.
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | mobl | [deg] Mean obliquity of the ecliptic in degrees. It may be NULL if not required. |
[out] | tobl | [deg] True obliquity of the ecliptic in degrees. It may be NULL if not required. |
[out] | ee | [deg] Equation of the equinoxes in seconds of time. It may be NULL if not required. |
[out] | dpsi | [arcsec] Nutation in longitude in arcseconds. It may be NULL if not required. |
[out] | deps | [arcsec] Nutation in obliquity in arcseconds. It may be NULL if not required. |
References ee_ct(), EPS_COR, mean_obliq(), NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, nutation_angles(), and PSI_COR.
short ecl2equ_vec | ( | double | jd_tt, |
enum novas_equator_type | coord_sys, | ||
enum novas_accuracy | accuracy, | ||
const double * | in, | ||
double * | out | ||
) |
Converts an ecliptic position vector to an equatorial position vector. To convert ecliptic coordinates (mean ecliptic and equinox of J2000.0) to GCRS RA and dec to, 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 type of the coordinates | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
in | Position vector, referred to specified ecliptic and equinox of date. | |
[out] | out | Position vector, referred to specified equator and equinox of date. It can be the same vector as the input. |
References e_tilt(), frame_tie(), J2000_TO_ICRS, mean_obliq(), NOVAS_FULL_ACCURACY, NOVAS_GCRS_EQUATOR, NOVAS_MEAN_EQUATOR, NOVAS_REDUCED_ACCURACY, NOVAS_TRUE_EQUATOR, and tt2tdb().
double ee_ct | ( | double | jd_tt_high, |
double | jd_tt_low, | ||
enum novas_accuracy | accuracy | ||
) |
Computes the "complementary terms" of the equation of the equinoxes. The input Julian date can be split into high and low order parts for improved accuracy. Typically, the split is into integer and fractiona parts. If the precision of a single part is sufficient, you may set the low order part to 0.
The series used in this function was derived from the first reference. This same series was also adopted for use in the IAU's Standards of Fundamental Astronomy (SOFA) software (i.e., subroutine eect00.for and function eect00.c
).
The low-accuracy series used in this function is a simple implementation derived from the first reference, in which terms smaller than 2 microarcseconds have been omitted.
REFERENCES:
jd_tt_high | [day] High-order part of TT based Julian date. |
jd_tt_low | [day] Low-order part of TT based Julian date. |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) |
References accum_prec(), novas_delaunay_args::D, novas_delaunay_args::F, fund_args(), NOVAS_FULL_ACCURACY, NOVAS_MERCURY, NOVAS_NEPTUNE, novas_delaunay_args::Omega, and planet_lon().
short ephemeris | ( | const double * | jd_tdb, |
const object * | body, | ||
enum novas_origin | origin, | ||
enum novas_accuracy | accuracy, | ||
double * | pos, | ||
double * | vel | ||
) |
Retrieves the position and velocity of a solar system body from a fundamental ephemeris.
It is recommended that the input structure 'cel_obj' be created using make_object()
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
body | Pointer to structure containing the designation of the body of interest | |
origin | NOVAS_BARYCENTER (0) or NOVAS_HELIOCENTER (1) | |
accuracy | NOCAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos | [AU] Pointer to structure containing the designation of the body of interest |
[out] | vel | [AU/day] Velocity vector of the body at 'jd_tdb'; equatorial rectangular coordinates in AU/day referred to the ICRS. |
cel_obj->type
is invalid, 10 + the error code from solarsystem(), or 20 + the error code from readeph().References novas_orbital_system::center, ephemeris(), make_planet(), object::name, NOVAS_BARYCENTER, NOVAS_EPHEM_OBJECT, NOVAS_FULL_ACCURACY, NOVAS_HELIOCENTER, novas_orbit_posvel(), NOVAS_ORBITAL_OBJECT, NOVAS_ORIGIN_TYPES, NOVAS_PLANET, NOVAS_SSB, NOVAS_SUN, object::number, object::orbit, readeph(), novas_orbital::system, and object::type.
short equ2ecl | ( | double | jd_tt, |
enum novas_equator_type | coord_sys, | ||
enum novas_accuracy | accuracy, | ||
double | ra, | ||
double | dec, | ||
double * | elon, | ||
double * | elat | ||
) |
Convert right ascension and declination to ecliptic longitude and latitude. To convert GCRS RA and dec to 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) | |
ra | [h] Right ascension in hours, referred to specified equator and equinox of date. | |
dec | [deg] Declination in degrees, referred to specified equator and equinox of date. | |
[out] | elon | [deg] Ecliptic longitude in degrees, referred to specified ecliptic and equinox of date. |
[out] | elat | [deg] Ecliptic latitude in degrees, referred to specified ecliptic and equinox of date. |
References equ2ecl_vec().
short equ2ecl_vec | ( | double | jd_tt, |
enum novas_equator_type | coord_sys, | ||
enum novas_accuracy | accuracy, | ||
const double * | in, | ||
double * | out | ||
) |
Converts an equatorial position vector to an ecliptic position vector. To convert ICRS RA and dec to 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 type of the coordinates. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
in | Position vector, referred to specified equator and equinox of date. | |
[out] | out | Position vector, referred to specified ecliptic and equinox of date. It can be the same vector as the input. If 'coord_sys' is NOVAS_GCRS_EQUATOR(2), the input GCRS coordinates are converted to J2000 ecliptic coordinates. |
References e_tilt(), frame_tie(), ICRS_TO_J2000, mean_obliq(), NOVAS_FULL_ACCURACY, NOVAS_GCRS_EQUATOR, NOVAS_MEAN_EQUATOR, NOVAS_REDUCED_ACCURACY, NOVAS_TRUE_EQUATOR, and tt2tdb().
int equ2gal | ( | double | ra, |
double | dec, | ||
double * | glon, | ||
double * | glat | ||
) |
Converts ICRS right ascension and declination to galactic longitude and latitude.
REFERENCES:
ra | [h] ICRS right ascension in hours. | |
dec | [deg] ICRS declination in degrees. | |
[out] | glon | [deg] Galactic longitude in degrees. |
[out] | glat | [deg] Galactic latitude in degrees. |
int equ2hor | ( | double | jd_ut1, |
double | ut1_to_tt, | ||
enum novas_accuracy | accuracy, | ||
double | xp, | ||
double | yp, | ||
const on_surface * | location, | ||
double | ra, | ||
double | dec, | ||
enum novas_refraction_model | ref_option, | ||
double * | zd, | ||
double * | az, | ||
double * | rar, | ||
double * | decr | ||
) |
Transforms topocentric (TOD) right ascension and declination to zenith distance and azimuth. This method should not be used to convert CIRS apparent coordinates (IAU 2000 standard) – for those you should use cirs_to_itrs() followed by itrs_to_hor() instead.
It uses a method that properly accounts for polar motion, which is significant at the sub-arcsecond level. This function can also adjust coordinates for atmospheric refraction.
NOTES:
REFERENCES:
jd_ut1 | [day] UT1 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 reference pole, in arcseconds. | |
yp | [arcsec] Conventionally-defined y coordinate of celestial intermediate pole with respect to ITRS reference pole, in arcseconds. | |
location | The observer location | |
ra | [h] Topocentric right ascension of object of interest, in hours, referred to true equator and equinox of date. | |
dec | [deg] Topocentric declination of object of interest, in degrees, referred to true equator and equinox of date. | |
ref_option | NOVAS_STANDARD_ATMOSPHERE (1), or NOVAS_WEATHER_AT_LOCATION (2) if to use the weather | |
[out] | zd | [deg] Topocentric zenith distance in degrees (unrefracted). |
[out] | az | [deg] Topocentric azimuth (measured east from north) in degrees. |
[out] | rar | [h] Topocentric right ascension of object of interest, in hours, referred to true equator and equinox of date, affected by refraction if 'ref_option' is non-zero. (It may be NULL if not required) |
[out] | decr | [deg] Topocentric declination of object of interest, in degrees, referred to true equator and equinox of date. (It may be NULL if not required) |
References EROT_GST, on_surface::latitude, on_surface::longitude, NOVAS_DYNAMICAL_CLASS, refract_astro(), and ter2cel().
double era | ( | double | jd_ut1_high, |
double | jd_ut1_low | ||
) |
Returns the value of the Earth Rotation Angle (theta) for a given UT1 Julian date. The expression used is taken from the note to IAU Resolution B1.8 of 2000. The input Julian date cane be split into an into high and low order parts (e.g. integer and fractional parts) for improved accuracy, or else one of the components (e.g. the low part) can be set to zero if no split is desired.
The algorithm used here is equivalent to the canonical theta = 0.7790572732640 + 1.00273781191135448 * t, where t is the time in days from J2000 (t = jd_high + jd_low - JD_J2000), but it avoids many two-PI 'wraps' that decrease precision (adopted from SOFA Fortran routine iau_era00; see also expression at top of page 35 of IERS Conventions (1996)).
REFERENCES:
jd_ut1_high | [day] High-order part of UT1 Julian date. |
jd_ut1_low | [day] Low-order part of UT1 Julian date. |
int frame_tie | ( | const double * | in, |
enum novas_frametie_direction | direction, | ||
double * | out | ||
) |
Transforms a vector from the dynamical reference system to the International Celestial Reference System (ICRS), or vice versa. The dynamical reference system is based on the dynamical mean equator and equinox of J2000.0. The ICRS is based on the space-fixed ICRS axes defined by the radio catalog positions of several hundred extragalactic objects.
For geocentric coordinates, the same transformation is used between the dynamical reference system and the GCRS.
NOTES:
REFERENCES:
in | Position vector, equatorial rectangular coordinates. | |
direction | <0 for for dynamical to ICRS transformation, or else >=0 for ICRS to dynamical transformation. Alternatively you may use the constants J2000_TO_ICRS (-1; or negative) or ICRS_TO_J2000 (0; or positive). | |
[out] | out | Position vector, equatorial rectangular coordinates. It can be the same vector as the input. |
int fund_args | ( | double | t, |
novas_delaunay_args * | a | ||
) |
Compute the fundamental arguments (mean elements) of the Sun and Moon.
REFERENCES:
t | [cy] TDB time in Julian centuries since J2000.0 | |
[out] | a | [rad] Fundamental arguments data to populate (5 doubles) [0:2π] |
References novas_delaunay_args::D, novas_delaunay_args::F, novas_delaunay_args::l, novas_delaunay_args::l1, norm_ang(), and novas_delaunay_args::Omega.
short gcrs2equ | ( | double | jd_tt, |
enum novas_dynamical_type | sys, | ||
enum novas_accuracy | accuracy, | ||
double | rag, | ||
double | decg, | ||
double * | ra, | ||
double * | dec | ||
) |
Converts GCRS right ascension and declination to coordinates with respect to the equator of date (mean or true). For coordinates with respect to the true equator of date, the origin of right ascension can be either the true equinox or the celestial intermediate origin (CIO). This function only supports the CIO-based method.
jd_tt | [day] Terrestrial Time (TT) based Julian date. (Unused if 'coord_sys' is NOVAS_ICRS_EQUATOR) | |
sys | Dynamical equatorial system type | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) (unused if 'coord_sys' is not NOVAS_ICRS [3]) | |
rag | [h] GCRS right ascension in hours. | |
decg | [deg] GCRS declination in degrees. | |
[out] | ra | [h] Right ascension in hours, referred to specified equator and right ascension origin of date. |
[out] | dec | [deg] Declination in degrees, referred to specified equator of date. |
References DEG2RAD, gcrs_to_cirs(), gcrs_to_mod(), gcrs_to_tod(), NOVAS_DYNAMICAL_CIRS, NOVAS_DYNAMICAL_MOD, NOVAS_DYNAMICAL_TOD, tt2tdb(), and vector2radec().
int gcrs_to_cirs | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
const double * | in, | ||
double * | out | ||
) |
Transforms a rectangular equatorial (x, y, z) vector from the Geocentric Celestial Reference System (GCRS) to the Celestial Intermediate Reference System (CIRS) frame at the given epoch
jd_tdb | [day] Barycentric Dynamical Time (TDB) 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 | GCRS Input (x, y, z) position or velocity vector | |
[out] | out | Output position or velocity 3-vector in the True equinox of Date coordinate frame. It can be the same vector as the input. |
References cio_basis(), and cio_location().
int gcrs_to_j2000 | ( | const double * | in, |
double * | out | ||
) |
Change GCRS coordinates to J2000 coordinates. Same as frame_tie() called with ICRS_TO_J2000
in | GCRS input 3-vector | |
[out] | out | J2000 output 3-vector |
References frame_tie(), and ICRS_TO_J2000.
int gcrs_to_mod | ( | double | jd_tdb, |
const double * | in, | ||
double * | out | ||
) |
Transforms a rectangular equatorial (x, y, z) vector from the Geocentric Celestial Reference System (GCRS) to the Mean of Date (MOD) reference frame at the given epoch
jd_tdb | [day] Barycentric Dynamical 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 | |
in | GCRS Input (x, y, z) position or velocity vector | |
[out] | out | Output position or velocity 3-vector in the Mean wquinox of Date coordinate frame. It can be the same vector as the input. |
References frame_tie(), ICRS_TO_J2000, NOVAS_JD_J2000, and precession().
int gcrs_to_tod | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
const double * | in, | ||
double * | out | ||
) |
Transforms a rectangular equatorial (x, y, z) vector from the Geocentric Celestial Reference System (GCRS) to the True of Date (TOD) reference frame at the given epoch
jd_tdb | [day] Barycentric Dynamical 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 | GCRS Input (x, y, z) position or velocity vector | |
[out] | out | Output position or velocity 3-vector in the True equinox of Date coordinate frame. It can be the same vector as the input. |
References frame_tie(), ICRS_TO_J2000, and j2000_to_tod().
short geo_posvel | ( | double | jd_tt, |
double | ut1_to_tt, | ||
enum novas_accuracy | accuracy, | ||
const observer * | obs, | ||
double * | pos, | ||
double * | vel | ||
) |
Computes the geocentric position and velocity of an observer. The final vectors are expressed in the GCRS.
jd_tt | [day] 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) | |
obs | Observer location | |
[out] | pos | [AU] Position 3-vector of observer, with respect to origin at geocenter, referred to GCRS axes, components in AU. (It may be NULL if not required.) |
[out] | vel | [AU/day] Velocity 3-vector of observer, with respect to origin at geocenter, referred to GCRS axes, components in AU/day. (It must be distinct from the pos output vector, and may be NULL if not required) |
References AU_KM, e_tilt(), ephemeris(), EROT_ERA, geo_posvel(), observer::near_earth, NOVAS_AIRBORNE_OBSERVER, NOVAS_BARYCENTER, NOVAS_EARTH_INIT, NOVAS_FULL_ACCURACY, NOVAS_OBSERVER_AT_GEOCENTER, NOVAS_OBSERVER_IN_EARTH_ORBIT, NOVAS_OBSERVER_ON_EARTH, NOVAS_REDUCED_ACCURACY, NOVAS_SOLAR_SYSTEM_OBSERVER, observer::on_surf, in_space::sc_pos, in_space::sc_vel, sidereal_time(), terra(), tod_to_gcrs(), tt2tdb(), and observer::where.
novas_ephem_provider get_ephem_provider | ( | ) |
Returns the user-defined ephemeris accessor function.
novas_planet_provider get_planet_provider | ( | ) |
Returns the custom (low-precision) ephemeris provider function for major planets (and Sun, Moon, SSB...), if any.
novas_planet_provider_hp get_planet_provider_hp | ( | ) |
Returns the custom high-precision ephemeris provider function for major planets (and Sun, Moon, SSB...), if any.
short grav_def | ( | double | jd_tdb, |
enum novas_observer_place | unused, | ||
enum novas_accuracy | accuracy, | ||
const double * | pos_src, | ||
const double * | pos_obs, | ||
double * | out | ||
) |
Computes the total gravitational deflection of light for the observed object 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 NOVAS_FULL_ACCURACY (0), the deflections due to the Sun, Jupiter, Saturn, and Earth are calculated. Otherwise, only the deflection due to the Sun is calculated. In either case, deflection for a given body is ignored if the observer is within ~1500 km of the center of the gravitating body.
NOTES:
REFERENCES:
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date | |
unused | The type of observer frame (no longer used) | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1). In full accuracy mode, it will calculate the deflection due to the Sun, Jupiter, Saturn and Earth. In reduced accuracy mode, only the deflection due to the Sun is calculated. | |
pos_src | [AU] 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] Position vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, corrected for 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_planets(), NOVAS_FULL_ACCURACY, and obs_planets().
int grav_planets | ( | const double * | pos_src, |
const double * | pos_obs, | ||
const novas_planet_bundle * | planets, | ||
double * | out | ||
) |
Computes the total gravitational deflection of light for the observed object due to the specified gravitating bodies in the solar system. This function is valid for an observed body within the solar system as well as for a star.
NOTES:
REFERENCES:
pos_src | [AU] 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] Position vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, corrected for gravitational deflection, components in AU. It can be the same vector as the input, but not the same as pos_obs. |
References C_AUDAY, d_light(), grav_vec(), novas_planet_bundle::mask, NOVAS_PLANETS, NOVAS_RMASS_INIT, novas_planet_bundle::pos, and novas_planet_bundle::vel.
int grav_vec | ( | const double * | pos_src, |
const double * | pos_obs, | ||
const double * | pos_body, | ||
double | rmass, | ||
double * | out | ||
) |
Corrects position vector for the deflection of light in the gravitational field of an arbitrary body. This function valid for an observed body within the solar system as well as for a star.
NOTES:
REFERENCES:
pos_src | [AU] Position 3-vector of observed object, with respect to origin at observer (or the geocenter), components in AU. | |
pos_obs | [AU] Position vector of gravitating body, with respect to origin at solar system barycenter, components in AU. | |
pos_body | [AU] Position 3-vector of gravitating body, with respect to origin at solar system barycenter, components in AU. | |
rmass | [1/Msun] Reciprocal mass of gravitating body in solar mass units, that is, Sun mass / body mass. | |
[out] | out | [AU] Position 3-vector of observed object, with respect to origin at observer (or the geocenter), corrected for gravitational deflection, components in AU. It can the same vector as the input. |
double ira_equinox | ( | double | jd_tdb, |
enum novas_equinox_type | equinox, | ||
enum novas_accuracy | accuracy | ||
) |
Compute the intermediate right ascension of the equinox at the input Julian date, using an analytical expression for the accumulated precession in right ascension. For the true equinox, the result is the equation of the origins.
NOTES:
REFERENCES:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date |
equinox | NOVAS_MEAN_EQUINOX (0) or NOVAS_TRUE_EQUINOX (1; or non-zero) |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) |
References e_tilt(), NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, and NOVAS_TRUE_EQUINOX.
int j2000_to_tod | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
const double * | in, | ||
double * | out | ||
) |
Transforms a rectangular equatorial (x, y, z) vector from J2000 coordinates to the True of Date (TOD) reference frame at the given epoch
jd_tdb | [day] Barycentric Dynamical Time (TDB) 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 | Input (x, y, z) position or velocity vector in rectangular equatorial coordinates at J2000 | |
[out] | out | Output position or velocity 3-vector in the True equinox of Date coordinate frame. It can be the same vector as the input. |
References NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, NUTATE_MEAN_TO_TRUE, nutation(), and precession().
double julian_date | ( | short | year, |
short | month, | ||
short | day, | ||
double | hour | ||
) |
Returns the Julian date for a given Gregorian calendar date. This function makes no checks for a valid input calendar date. Input calendar date must be Gregorian. Input time value can be based on any UT-like time scale (UTC, UT1, TT, etc.) - output Julian date will have the same basis.
REFERENCES:
year | [yr] Gregorian calendar year |
month | [month] Gregorian calendar month [1:12] |
day | [day] Day of month [1:31] |
hour | [hr] Hour of day [0:24] |
short light_time | ( | double | jd_tdb, |
const object * | body, | ||
const double * | pos_obs, | ||
double | tlight0, | ||
enum novas_accuracy | accuracy, | ||
double * | pos_src_obs, | ||
double * | tlight | ||
) |
Computes the geocentric position of a solar system body, as antedated for light-time.
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date | |
body | Pointer to structure containing the designation for the solar system body | |
pos_obs | [AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes, components in AU. | |
tlight0 | [day] First approximation to light-time, in days (can be set to 0.0 if unknown). | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos_src_obs | [AU] Position 3-vector of body, with respect to origin at observer (or the geocenter), referred to ICRS axes, components in AU. It can be the same vector as either of the inputs. |
[out] | tlight | [day] Calculated light time |
References light_time2().
int light_time2 | ( | double | jd_tdb, |
const object * | body, | ||
const double * | pos_obs, | ||
double | tlight0, | ||
enum novas_accuracy | accuracy, | ||
double * | p_src_obs, | ||
double * | v_ssb, | ||
double * | tlight | ||
) |
Computes the geocentric position and velocity of a solar system body, as antedated for light-time. It is effectively the same as the original NOVAS C light_time(), except that this returns the antedated source velocity vector also.
NOTES:
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date | |
body | Pointer to structure containing the designation for the solar system body | |
pos_obs | [AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes, components in AU. | |
tlight0 | [day] First approximation to light-time, in days (can be set to 0.0 if unknown). | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | p_src_obs | [AU] Position 3-vector of body, relative to observer, referred to ICRS axes, components in AU. |
[out] | v_ssb | [AU/day] Velocity 3-vector of body, with respect to the Solar-system barycenter, referred to ICRS axes, components in AU/day. |
[out] | tlight | [day] Calculated light time, or NAN when returning with an error code. |
References bary2obs(), ephemeris(), NOVAS_BARYCENTER, NOVAS_FULL_ACCURACY, and novas_inv_max_iter.
int limb_angle | ( | const double * | pos_src, |
const double * | pos_obs, | ||
double * | limb_ang, | ||
double * | nadir_ang | ||
) |
Determines the angle of an object above or below the Earth's limb (horizon). The geometric limb is computed, assuming the Earth to be an airless sphere (no refraction or oblateness is included). The observer can be on or above the Earth. For an observer on the surface of the Earth, this function returns the approximate unrefracted elevation.
pos_src | [AU] Position 3-vector of observed object, with respect to origin at geocenter, components in AU. | |
pos_obs | [AU] Position 3-vector of observer, with respect to origin at geocenter, components in AU. | |
[out] | limb_ang | [deg] Angle of observed object above (+) or below (-) limb in degrees, or NAN if reurning with an error. It may be NULL if not required. |
[out] | nadir_ang | Nadir angle of observed object as a fraction of apparent radius of limb: lt;1.0 if below the limb; 1.0 on the limb; or >1.0 if above the limb. Returns NAN in case of an error return. It may be NULL if not required. |
References M_PI.
short local_planet | ( | double | jd_tt, |
const object * | ss_body, | ||
double | ut1_to_tt, | ||
const on_surface * | position, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec, | ||
double * | dis | ||
) |
Computes the local apparent place of a solar system body, in the GCRS. This is the same as calling place() for the body for the same observer location and NOVAS_GCRS as the reference system, except the different set of return values used.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
ss_body | Pointer to structure containing the body designation for the solar system body. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
position | Position of the observer | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Local right ascension in hours, referred to the GCRS (it may be NULL if not required). |
[out] | dec | [deg] Local right ascension in hours, referred to the GCRS (it may be NULL if not required). |
[out] | dis | [AU] True distance from Earth to the body at 'jd_tt' in AU (it may be NULL if not required). |
References make_observer(), NOVAS_GCRS, NOVAS_OBSERVER_ON_EARTH, and radec_planet().
short local_star | ( | double | jd_tt, |
double | ut1_to_tt, | ||
const cat_entry * | star, | ||
const on_surface * | position, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec | ||
) |
Computes the local apparent place of a star at date 'jd_tt', in the GCRS, given its catalog mean place, proper motion, parallax, and radial velocity.
Notwithstanding the different set of return values, this is the same as calling place_star() with the same observer location NOVAS_GCRS for an object that specifies the star.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
position | Position of the observer | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Local right ascension in hours, referred to the GCRS (it may be NULL if not required). |
[out] | dec | [deg] Local right ascension in hours, referred to the GCRS (it may be NULL if not required). |
References make_observer(), NOVAS_GCRS, NOVAS_OBSERVER_ON_EARTH, and radec_star().
short make_cat_entry | ( | const char * | star_name, |
const char * | catalog, | ||
long | cat_num, | ||
double | ra, | ||
double | dec, | ||
double | pm_ra, | ||
double | pm_dec, | ||
double | parallax, | ||
double | rad_vel, | ||
cat_entry * | star | ||
) |
Populates the data structure for a 'catalog' source, such as a star.
star_name | Object name (less than SIZE_OF_OBJ_NAME in length). It may be NULL if not relevant. | |
catalog | Catalog identifier (less than SIZE_OF_CAT_NAME in length). E.g. 'HIP' = Hipparcos, 'TY2' = Tycho-2. It may be NULL if not relevant. | |
cat_num | Object number in the catalog. | |
ra | [h] Right ascension of the object (hours). | |
dec | [deg] Declination of the object (degrees). | |
pm_ra | [mas/yr] Proper motion in right ascension (milliarcseconds/year). | |
pm_dec | [mas/yr] Proper motion in declination (milliarcseconds/year). | |
parallax | [mas] Parallax (milliarcseconds). | |
rad_vel | [km/s] Radial velocity (LSR) | |
[out] | star | Pointer to data structure to populate. |
References cat_entry::catalog, cat_entry::dec, cat_entry::parallax, cat_entry::promodec, cat_entry::promora, cat_entry::ra, rad_vel(), cat_entry::radialvelocity, cat_entry::starname, and cat_entry::starnumber.
int make_in_space | ( | const double * | sc_pos, |
const double * | sc_vel, | ||
in_space * | loc | ||
) |
Populates an 'in_space' data structure, for an observer situated on a near-Earth spacecraft, with the provided position and velocity components. Both input vectors are assumed with respect to true equator and equinox of date.
sc_pos | [km] Geocentric (x, y, z) position vector in km. NULL defaults to the origin | |
sc_vel | [km/s] Geocentric (x, y, z) velocity vector in km/s. NULL defaults to zero speed. | |
[out] | loc | Pointer to earth-orbit location data structure to populate. |
References in_space::sc_pos, and in_space::sc_vel.
short make_object | ( | enum novas_object_type | type, |
long | number, | ||
const char * | name, | ||
const cat_entry * | star, | ||
object * | source | ||
) |
Populates an object data structure using the parameters provided. By default (for compatibility with NOVAS C) source names are converted to upper-case internally. You can however enable case-sensitive processing by calling novas_case_sensitive() before.
NOTES:
orbit
field (added in v1.2) with zeroes to remain ABI compatible with versions <1.2, and to avoid the possiblity of segfaulting if used to initialize a legacy object
variable. type | The type of object. NOVAS_PLANET (0), NOVAS_EPHEM_OBJECT (1) or NOVAS_CATALOG_OBJECT (2), or NOVAS_ORBITAL_OBJECT (3). | |
number | The novas ID number (for solar-system bodies only, otherwise ignored) | |
name | The name of the object (case insensitive). It should be shorter than SIZE_OF_OBJ_NAME or else an error will be returned. The name is converted to upper internally unless novas_case_sensitive() was called before to change that. | |
star | Pointer to structure to populate with the catalog data for a celestial object located outside the solar system. Used only if type is NOVAS_CATALOG_OBJECT, otherwise ignored and can be NULL. | |
[out] | source | Pointer to the celestial object data structure to be populated. |
References object::name, NOVAS_CATALOG_OBJECT, NOVAS_OBJECT_TYPES, NOVAS_ORBITAL_OBJECT, NOVAS_PLANET, NOVAS_PLANETS, object::number, object::orbit, object::star, and object::type.
short make_observer | ( | enum novas_observer_place | where, |
const on_surface * | loc_surface, | ||
const in_space * | loc_space, | ||
observer * | obs | ||
) |
Populates an 'observer' data structure given the parameters. The output data structure may be used an the the inputs to NOVAS-C function 'place()'.
where | The location type of the observer | |
loc_surface | Pointer to data structure that defines a location on Earth's surface. Used only if 'where' is NOVAS_OBSERVER_ON_EARTH, otherwise can be NULL. | |
loc_space | Pointer to data structure that defines a near-Earth location in space. Used only if 'where' is NOVAS_OBSERVER_IN_EARTH_ORBIT, otherwise can be NULL. | |
[out] | obs | Pointer to observer data structure to populate. |
References observer::near_earth, NOVAS_AIRBORNE_OBSERVER, NOVAS_OBSERVER_AT_GEOCENTER, NOVAS_OBSERVER_IN_EARTH_ORBIT, NOVAS_OBSERVER_ON_EARTH, NOVAS_SOLAR_SYSTEM_OBSERVER, observer::on_surf, in_space::sc_vel, and observer::where.
int make_observer_at_geocenter | ( | observer * | obs | ) |
Populates an 'observer' data structure for a hypothetical observer located at Earth's geocenter. The output data structure may be used an the the inputs to NOVAS-C function 'place()'.
[out] | obs | Pointer to data structure to populate. |
References make_observer(), and NOVAS_OBSERVER_AT_GEOCENTER.
int make_observer_in_space | ( | 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. Both input vectors are with respect to true equator and equinox of date. The output data structure may be used an the the inputs to NOVAS-C function 'place()'.
sc_pos | [km] Geocentric (x, y, z) position vector in km. | |
sc_vel | [km/s] Geocentric (x, y, z) velocity vector in km/s. | |
[out] | obs | Pointer to the data structure to populate |
References make_in_space(), make_observer(), and NOVAS_OBSERVER_IN_EARTH_ORBIT.
int make_observer_on_surface | ( | double | latitude, |
double | longitude, | ||
double | height, | ||
double | temperature, | ||
double | pressure, | ||
observer * | obs | ||
) |
Populates and 'on_surface' data structure with the specified location defining parameters of the observer. The output data structure may be used an the the inputs to NOVAS-C function 'place()'.
latitude | [deg] Geodetic (ITRS) latitude in degrees; north positive. | |
longitude | [deg] Geodetic (ITRS) longitude in degrees; east positive. | |
height | [m] Altitude over se level of the observer (meters). | |
temperature | [C] Temperature (degrees Celsius). | |
pressure | [mbar] Atmospheric pressure (millibars). | |
[out] | obs | Pointer to the data structure to populate. |
References make_observer(), make_on_surface(), and NOVAS_OBSERVER_ON_EARTH.
int make_on_surface | ( | double | latitude, |
double | longitude, | ||
double | height, | ||
double | temperature, | ||
double | pressure, | ||
on_surface * | loc | ||
) |
Populates an 'on_surface' data structure, for an observer on the surface of the Earth, with the given parameters.
Note, that because this is an original NOVAS C routine, it does not have an argument to set a humidity value (e.g. for radio refraction). As such, the humidity value remains undefined after this call. To set the humidity, set the output structure's field after calling this funcion. Its unit is [%], and so the range is 0.0–100.0.
latitude | [deg] Geodetic (ITRS) latitude in degrees; north positive. | |
longitude | [deg] Geodetic (ITRS) longitude in degrees; east positive. | |
height | [m] Altitude over se level of the observer (meters). | |
temperature | [C] Temperature (degrees Celsius). | |
pressure | [mbar] Atmospheric pressure (millibars). | |
[out] | loc | Pointer to Earth location data structure to populate. |
References on_surface::height, on_surface::latitude, on_surface::longitude, on_surface::pressure, and on_surface::temperature.
int make_planet | ( | enum novas_planet | num, |
object * | planet | ||
) |
Sets a celestial object to be a major planet, or the Sun, Moon, Solar-system Barycenter, etc.
num | Planet ID number (NOVAS convention) | |
[out] | planet | Pointer to structure to populate. |
References make_object(), NOVAS_PLANET, NOVAS_PLANET_NAMES_INIT, and NOVAS_PLANETS.
double mean_obliq | ( | double | jd_tdb | ) |
Computes the mean obliquity of the ecliptic.
REFERENCES:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date |
short mean_star | ( | double | jd_tt, |
double | tra, | ||
double | tdec, | ||
enum novas_accuracy | accuracy, | ||
double * | ira, | ||
double * | idec | ||
) |
Computes the ICRS position of a star, given its True of Date (TOD) apparent place at date 'jd_tt'. Proper motion, parallax and radial velocity are assumed to be zero.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
tra | [h] Apparent (TOD) right ascension in hours, referred to true equator and equinox of date. | |
tdec | [deg] Apparent (TOD) declination in degrees, referred to true equator and equinox of date. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ira | [h] ICRS right ascension in hours, or NAN when returning with an error code. |
[out] | idec | [deg] ICRS declination in degrees, or NAN when returning with an error code. |
References app_star(), CAT_ENTRY_INIT, cat_entry::dec, novas_inv_max_iter, precession(), cat_entry::ra, starvectors(), and vector2radec().
int mod_to_gcrs | ( | double | jd_tdb, |
const double * | in, | ||
double * | out | ||
) |
Transforms a rectangular equatorial (x, y, z) vector from Mean of Date (MOD) reference frame at the given epoch to the Geocentric Celestial Reference System(GCRS)
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date that defines the input epoch. Typically it does not require much precision, and Julian dates in other time measures will be unlikely to affect the result | |
in | Input (x, y, z) position or velocity 3-vector in the Mean equinox of Date coordinate frame. | |
[out] | out | Output GCRS position or velocity vector. It can be the same vector as the input. |
References frame_tie(), J2000_TO_ICRS, NOVAS_JD_J2000, and precession().
double norm_ang | ( | double | angle | ) |
Returns the normalized angle in the [0:2π) range.
angle | [rad] an angle in radians. |
References TWOPI.
void novas_case_sensitive | ( | int | value | ) |
Enables or disables case-sensitive processing of the object name. The effect is not retroactive. The setting will only affect the celestial objects that are defined after the call. Note, that catalog names, set via make_cat_entry() are always case sensitive regardless of this setting.
value | (boolean) TRUE (non-zero) to enable case-sensitive object names, or else FALSE (0) to convert names to upper case only (NOVAS C compatible behavior). |
void novas_debug | ( | enum novas_debug_mode | mode | ) |
Enables or disables reporting errors and traces to the standard error stream.
mode | NOVAS_DEBUG_OFF (0; or <0), NOVAS_DEBUG_ON (1), or NOVAS_DEBUG_EXTRA (2; or >2). |
References NOVAS_DEBUG_EXTRA.
enum novas_debug_mode novas_get_debug_mode | ( | ) |
Returns the current, thread-local, mode for reporting errors encountered (and traces).
int novas_orbit_posvel | ( | double | jd_tdb, |
const novas_orbital * | orbit, | ||
enum novas_accuracy | accuracy, | ||
double * | pos, | ||
double * | vel | ||
) |
Calculates a rectangular equatorial position and velocity vector for the given orbital elements for the specified time of observation.
REFERENCES:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
orbit | Orbital parameters | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1). | |
[out] | pos | [AU] Output ICRS equatorial position vector, or NULL if not required |
[out] | vel | [AU/day] Output ICRS equatorial velocity vector, or NULL if not required |
References novas_orbital::a, novas_orbital::apsis_period, novas_orbital::e, novas_orbital::i, novas_orbital::jd_tdb, novas_orbital::M0, novas_orbital::n, novas_orbital::node_period, novas_inv_max_iter, novas_orbital::omega, novas_orbital::Omega, novas_orbital::system, and TWOPI.
double novas_z2v | ( | double | z | ) |
Converts a redshift value (z = δf / frest) to a radial velocity (i.e. rate) of recession. It is based on the relativistic formula:
1 + z = sqrt((1 + β) / (1 - β))
where β = v / c.
z | the redshift value (δλ / λrest). |
int nutation | ( | double | jd_tdb, |
enum novas_nutation_direction | direction, | ||
enum novas_accuracy | accuracy, | ||
const double * | in, | ||
double * | out | ||
) |
Nutates equatorial rectangular coordinates from mean equator and equinox of epoch to true equator and equinox of epoch. Inverse transformation may be applied by setting flag 'direction'.
This is the old (pre IAU 2006) method of nutation calculation. If you follow the now standard IAU 2000/2006 methodology you will want to use nutation_angles() instead.
REFERENCES:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
direction | NUTATE_MEAN_TO_TRUE (0) or NUTATE_TRUE_TO_MEAN (-1; or non-zero) | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
in | Position 3-vector, geocentric equatorial rectangular coordinates, referred to mean equator and equinox of epoch. | |
[out] | out | Position vector, geocentric equatorial rectangular coordinates, referred to true equator and equinox of epoch. It can be the same as the input position. |
References e_tilt(), and NUTATE_MEAN_TO_TRUE.
int nutation_angles | ( | double | t, |
enum novas_accuracy | accuracy, | ||
double * | dpsi, | ||
double * | deps | ||
) |
Returns the values for nutation in longitude and nutation in obliquity for a given TDB Julian date. The nutation model selected depends upon the input value of 'accuracy'. See notes below for important details.
This function selects the nutation model depending first upon the input value of 'accuracy'. If 'accuracy' is NOVAS_FULL_ACCURACY (0), the IAU 2000A nutation model is used. Otherwise the model set by set_nutation_lp_provider() is used, or else the default nu2000k().
See the prologs of the nutation functions in file 'nutation.c' for details concerning the models.
REFERENCES:
t | [cy] TDB time in Julian centuries since J2000.0 | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | dpsi | [arcsec] Nutation in longitude in arcseconds. |
[out] | deps | [arcsec] Nutation in obliquity in arcseconds. |
References iau2000a(), and NOVAS_FULL_ACCURACY.
int obs_planets | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
const double * | pos_obs, | ||
int | pl_mask, | ||
novas_planet_bundle * | planets | ||
) |
Calculates the positions and velocities for the Solar-system bodies, e.g. for use for graviational deflection calculations. The planet positions are calculated relative to the observer location, while velocities are w.r.t. the SSB. Both positions and velocities are antedated for light travel time, so they accurately reflect the apparent position (and barycentric motion) of the bodies from the observer's perspective.
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1). In full accuracy mode, it will calculate the deflection due to the Sun, Jupiter, Saturn and Earth. In reduced accuracy mode, only the deflection due to the Sun is calculated. | |
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. | |
pl_mask | Bitwise (1 << planet-number) mask indicating which planets to request data for. See enum novas_planet for the enumeration of planet numbers. | |
[out] | planets | Pointer to apparent planet data to populate. have positions and velocities calculated successfully. See enum novas_planet for the enumeration of planet numbers. |
References light_time2(), make_planet(), novas_planet_bundle::mask, novas_debug(), NOVAS_DEBUG_EXTRA, NOVAS_DEBUG_OFF, novas_get_debug_mode(), NOVAS_PLANETS, NOVAS_SUN, novas_planet_bundle::pos, and novas_planet_bundle::vel.
int obs_posvel | ( | double | jd_tdb, |
double | ut1_to_tt, | ||
enum novas_accuracy | accuracy, | ||
const observer * | obs, | ||
const double * | geo_pos, | ||
const double * | geo_vel, | ||
double * | pos, | ||
double * | vel | ||
) |
Calculates the ICRS position and velocity of the observer relative to the Solar System Barycenter (SSB).
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date. | |
ut1_to_tt | [s] TT - UT1 time difference. Used only when 'location->where' is NOVAS_OBSERVER_ON_EARTH (1) or NOVAS_OBSERVER_IN_EARTH_ORBIT (2). | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
obs | The observer location, relative to which the output positions and velocities are to be calculated | |
geo_pos | [AU] ICRS position vector of the geocenter w.r.t. the Solar System Barycenter (SSB). If either geo_pos or geo_vel is NULL, it will be calculated when needed. | |
geo_vel | [AU/day] ICRS velocity vector of the geocenter w.r.t. the Solar System Barycenter (SSB). If either geo_pos or geo_vel is NULL, it will be calculated when needed. | |
[out] | pos | [AU] Position 3-vector of the observer w.r.t. the Solar System Barycenter (SSB). It may be NULL if not required. |
[out] | vel | [AU/day] Velocity 3-vector of the observer w.r.t. the Solar System Barycenter (SSB). It must be distinct from the pos output vector, and may be NULL if not required. |
References ephemeris(), geo_posvel(), observer::near_earth, NOVAS_AIRBORNE_OBSERVER, NOVAS_BARYCENTER, NOVAS_EARTH_INIT, NOVAS_OBSERVER_IN_EARTH_ORBIT, NOVAS_OBSERVER_ON_EARTH, NOVAS_OBSERVER_PLACES, NOVAS_SOLAR_SYSTEM_OBSERVER, in_space::sc_pos, in_space::sc_vel, and observer::where.
short place | ( | double | jd_tt, |
const object * | source, | ||
const observer * | location, | ||
double | ut1_to_tt, | ||
enum novas_reference_system | coord_sys, | ||
enum novas_accuracy | accuracy, | ||
sky_pos * | output | ||
) |
Computes the apparent direction of a celestial object at a specified time and in a specified coordinate system and a specific near-Earth origin.
While coord_sys
defines the celestial pole (i.e. equator) orientation of the coordinate system, location->where
sets the origin of the reference place relative to which positions and velocities are reported.
For all but ICRS coordinate outputs, the calculated positions and velocities include aberration corrections for the moving frame of the observer as well as gravitational deflection due to the Sun and Earth and other major gravitating bodies in the Solar system, provided planet positions are available via a novas_planet_provider function.
In case of a dynamical equatorial system (such as CIRS or TOD) and an Earth-based observer, the polar wobble parameters set via a prior call to cel_pole() together with he ut1_to_tt argument decide whether the resulting 'topocentric' output frame is Pseudo Earth Fixed (PEF; if cel_pole() was not set and DUT1 is 0) or ITRS (actual rotating Earth; if cel_pole() was set and ut1_to_tt includes the DUT1 component).
NOTES:
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
source | Pointer to a celestrial object data structure | |
location | The observer location, relative to which the output positions and velocities are to be calculated | |
ut1_to_tt | [s] TT - UT1 time difference. Used only when 'location->where' is NOVAS_OBSERVER_ON_EARTH (1) or NOVAS_OBSERVER_IN_EARTH_ORBIT (2). | |
coord_sys | The coordinate system that defines the orientation of the celestial pole. If it is NOVAS_ICRS (3), a geometric position and radial velocity is returned. For all other systems, the returned position is the apparent position including aberration and gravitational deflection corrections, and the radial velocity is in the direction the eflected light was emitted from the source. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | output | Data structure to populate with the result. |
References aberration(), bary2obs(), C_AUDAY, d_light(), sky_pos::dec, sky_pos::dis, ephemeris(), gcrs_to_cirs(), gcrs_to_j2000(), gcrs_to_mod(), gcrs_to_tod(), grav_bodies_full_accuracy, grav_bodies_reduced_accuracy, grav_planets(), light_time2(), make_observer_at_geocenter(), make_planet(), NOVAS_BARYCENTER, NOVAS_CATALOG_OBJECT, NOVAS_CIRS, NOVAS_EARTH, NOVAS_FULL_ACCURACY, NOVAS_ICRS, NOVAS_J2000, NOVAS_MOD, NOVAS_REDUCED_ACCURACY, NOVAS_REFERENCE_SYSTEMS, NOVAS_SUN, NOVAS_TOD, obs_planets(), obs_posvel(), proper_motion(), sky_pos::r_hat, sky_pos::ra, rad_vel2(), sky_pos::rv, object::star, starvectors(), tt2tdb(), object::type, and vector2radec().
int place_star | ( | double | jd_tt, |
const cat_entry * | star, | ||
const observer * | obs, | ||
double | ut1_to_tt, | ||
enum novas_reference_system | system, | ||
enum novas_accuracy | accuracy, | ||
sky_pos * | pos | ||
) |
Computes the apparent place of a star, referenced to dynamical equator at date 'jd_tt', given its catalog mean place, proper motion, parallax, and radial velocity. See place()
for more information.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
obs | Observer location (can be NULL if not relevant) | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
system | The type of coordinate reference system in which coordinates are to be returned. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos | The position and radial velocity of of the catalog source in the specified coordinate system and relative to the specified observer location (if applicable) |
References NOVAS_CATALOG_OBJECT, place(), object::star, and object::type.
double planet_lon | ( | double | t, |
enum novas_planet | planet | ||
) |
Returns the planetary longitude, for Mercury through Neptune, w.r.t. mean dynamical ecliptic and equinox of J2000, with high order terms omitted (Simon et al. 1994, 5.8.1-5.8.8).
t | [cy] Julian centuries since J2000 |
planet | Novas planet id, e.g. NOVAS_MARS. |
planet
id is out of range.References NOVAS_EARTH, NOVAS_JUPITER, NOVAS_MARS, NOVAS_MERCURY, NOVAS_NEPTUNE, NOVAS_SATURN, NOVAS_URANUS, NOVAS_VENUS, and TWOPI.
short precession | ( | double | jd_tdb_in, |
const double * | in, | ||
double | jd_tdb_out, | ||
double * | out | ||
) |
Precesses equatorial rectangular coordinates from one epoch to another. Unlike the original NOVAS routine, this routine works for any pairing of the time arguments.
This function calculates precession for the old (pre IAU 2000) methodology. Its main use for NOVAS users is to allow converting older catalog coordinates e.g. to J2000 coordinates, which then can be converted to the now standard ICRS system via frame_tie().
NOTE:
REFERENCES:
jd_tdb_in | [day] Barycentric Dynamic Time (TDB) based Julian date of the input epoch | |
in | Position 3-vector, geocentric equatorial rectangular coordinates, referred to mean dynamical equator and equinox of the initial epoch. | |
jd_tdb_out | [day] Barycentric Dynamic Time (TDB) based Julian date of the output epoch | |
[out] | out | Position 3-vector, geocentric equatorial rectangular coordinates, referred to mean dynamical equator and equinox of the final epoch. It can be the same vector as the input. |
References precession().
int proper_motion | ( | double | jd_tdb_in, |
const double * | pos, | ||
const double * | vel, | ||
double | jd_tdb_out, | ||
double * | out | ||
) |
Applies proper motion, including foreshortening effects, to a star's position.
REFERENCES:
jd_tdb_in | [day] Barycentric Dynamical Time (TDB) based Julian date of the first epoch. | |
pos | [AU] Position vector at first epoch. | |
vel | [AU/day] Velocity vector at first epoch. | |
jd_tdb_out | [day] Barycentric Dynamical Time (TDB) based Julian date of the second epoch. | |
[out] | out | Position vector at second epoch. It can be the same vector as the input. |
int rad_vel | ( | const object * | source, |
const double * | pos_src, | ||
const double * | vel_src, | ||
const double * | vel_obs, | ||
double | d_obs_geo, | ||
double | d_obs_sun, | ||
double | d_src_sun, | ||
double * | rv | ||
) |
Predicts the radial velocity of the observed object as it would be measured by spectroscopic means. Radial velocity is here defined as the radial velocity measure (z) times the speed of light. For major planets (and Sun and Moon), it includes gravitational corrections for light originating at the surface and observed from near Earth or else from a large distance away. For other solar system bodies, it applies to a fictitious emitter at the center of the observed object, assumed massless (no gravitational red shift). The corrections do not in general apply to reflected light. For stars, it includes all effects, such as gravitational redshift, contained in the catalog barycentric radial velocity measure, a scalar derived from spectroscopy. Nearby stars with a known kinematic velocity vector (obtained independently of spectroscopy) can be treated like solar system objects.
Gravitational blueshift corrections for the Solar and Earth potential for observers are included. However, the result does not include a blueshift correction for observers (e.g. spacecraft) orbiting other major Solar-system bodies. You may adjust the amount of gravitational redshift correction applied to the radial velocity via redshift_vrad()
, unredshift_vrad()
and grav_redshift()
if necessary.
All the input arguments are BCRS quantities, expressed with respect to the ICRS axes. 'vel_src' and 'vel_obs' are kinematic velocities - derived from geometry or dynamics, not spectroscopy.
If the object is outside the solar system, the algorithm used will be consistent with the IAU definition of stellar radial velocity, specifically, the barycentric radial velocity measure, which is derived from spectroscopy. In that case, the vector 'vel_src' can be very approximate – or, for distant stars or galaxies, zero – as it will be used only for a small geometric correction that is proportional to proper motion.
Any of the distances (last three input arguments) can be set to zero (0.0) or negative if the corresponding general relativistic gravitational potential term is not to be evaluated. These terms generally are important at the meter/second level only. If 'd_obs_geo' and 'd_obs_sun' are both zero, an average value will be used for the relativistic term for the observer, appropriate for an observer on the surface of the Earth. 'd_src_sun', if given, is used only for solar system objects.
NOTES:
d_obs_geo
and d_obs_sun
were zero. As of SuperNOVAS v1.1, the relatistic corrections for a moving observer will be included in the radial velocity measure always. REFERENCES:
source | Celestial object observed | |
pos_src | [AU|*] Geometric position vector of object with respect to observer. For solar system sources it should be corrected for light-time. For non-solar-system objects, the position vector defines a direction only, with arbitrary magnitude. | |
vel_src | [AU/day] Velocity vector of object with respect to solar system barycenter. | |
vel_obs | [AU/day] Velocity vector of observer with respect to solar system barycenter. | |
d_obs_geo | [AU] Distance from observer to geocenter, or <=0.0 if gravitational blueshifting due to Earth potential around observer can be ignored. | |
d_obs_sun | [AU] Distance from observer to Sun, or <=0.0 if gravitational bluehifting due to Solar potential around observer can be ignored. | |
d_src_sun | [AU] Distance from object to Sun, or <=0.0 if gravitational redshifting due to Solar potential around source can be ignored. | |
[out] | rv | [km/s] The observed radial velocity measure times the speed of light, or NAN if there was an error. |
References rad_vel2().
double rad_vel2 | ( | const object * | source, |
const double * | pos_emit, | ||
const double * | vel_src, | ||
const double * | pos_det, | ||
const double * | vel_obs, | ||
double | d_obs_geo, | ||
double | d_obs_sun, | ||
double | d_src_sun | ||
) |
Predicts the radial velocity of the observed object as it would be measured by spectroscopic means. This is a modified version of the original NOVAS C 3.1 rad_vel(), to account for the different directions in which light is emitted vs in which it detected, e.g. when it is gravitationally deflected.
Radial velocity is here defined as the radial velocity measure (z) times the speed of light. For major planets (and Sun and Moon), it includes gravitational corrections for light originating at the surface and observed from near Earth or else from a large distance away. For other solar system bodies, it applies to a fictitious emitter at the center of the observed object, assumed massless (no gravitational red shift). The corrections do not in general apply to reflected light. For stars, it includes all effects, such as gravitational redshift, contained in the catalog barycentric radial velocity measure, a scalar derived from spectroscopy. Nearby stars with a known kinematic velocity vector (obtained independently of spectroscopy) can be treated like solar system objects.
Gravitational blueshift corrections for the Solar and Earth potential for observers are included. However, the result does not include a blueshift correction for observers (e.g. spacecraft) orbiting other major Solar-system bodies. You may adjust the amount of gravitational redshift correction applied to the radial velocity via redshift_vrad()
, unredshift_vrad()
and grav_redshift()
if necessary.
All the input arguments are BCRS quantities, expressed with respect to the ICRS axes. 'vel_src' and 'vel_obs' are kinematic velocities - derived from geometry or dynamics, not spectroscopy.
If the object is outside the solar system, the algorithm used will be consistent with the IAU definition of stellar radial velocity, specifically, the barycentric radial velocity measure, which is derived from spectroscopy. In that case, the vector 'vel_src' can be very approximate – or, for distant stars or galaxies, zero – as it will be used only for a small geometric and relativistic (time dilation) correction, including the proper motion.
Any of the distances (last three input arguments) can be set to a negative value if the corresponding general relativistic gravitational potential term is not to be evaluated. These terms generally are important only at the meter/second level. If 'd_obs_geo' and 'd_obs_sun' are both zero, an average value will be used for the relativistic term for the observer, appropriate for an observer on the surface of the Earth. 'd_src_sun', if given, is used only for solar system objects.
NOTES:
REFERENCES:
source | Celestial object observed |
pos_emit | [AU|*] position vector of object with respect to observer in the direction that light was emitted from the source. For solar system sources it should be corrected for light-time. For non-solar-system objects, the position vector defines a direction only, with arbitrary magnitude. |
vel_src | [AU/day] Velocity vector of object with respect to solar system barycenter. |
pos_det | [AU|*] apparent position vector of source, as seen by the observer. It may be the same vector as pos_emit , in which case the routine behaves like the original NOVAS_C rad_vel(). |
vel_obs | [AU/day] Velocity vector of observer with respect to solar system barycenter. |
d_obs_geo | [AU] Distance from observer to geocenter, or <=0.0 if gravitational blueshifting due to Earth potential around observer can be ignored. |
d_obs_sun | [AU] Distance from observer to Sun, or <=0.0 if gravitational bluehifting due to Solar potential around observer can be ignored. |
d_src_sun | [AU] Distance from object to Sun, or <=0.0 if gravitational redshifting due to Solar potential around source can be ignored. Additionally, a value <0 will also skip corrections for light originating at the surface of the observed major solar-system body. |
References AU, C, C_AUDAY, cat_entry::dec, GE, GS, NOVAS_CATALOG_OBJECT, NOVAS_EARTH_RADIUS, NOVAS_EPHEM_OBJECT, NOVAS_KM, NOVAS_ORBITAL_OBJECT, NOVAS_PLANET, NOVAS_PLANET_GRAV_Z_INIT, NOVAS_PLANETS, NOVAS_SOLAR_RADIUS, novas_z2v(), object::number, cat_entry::parallax, cat_entry::ra, cat_entry::radialvelocity, object::star, and object::type.
int radec2vector | ( | double | ra, |
double | dec, | ||
double | dist, | ||
double * | pos | ||
) |
Converts equatorial spherical coordinates to a vector (equatorial rectangular coordinates).
ra | [h] Right ascension (hours). | |
dec | [deg] Declination (degrees). | |
dist | [AU] Distance (AU) | |
[out] | pos | [AU] Position 3-vector, equatorial rectangular coordinates (AU). |
int radec_planet | ( | double | jd_tt, |
const object * | ss_body, | ||
const observer * | obs, | ||
double | ut1_to_tt, | ||
enum novas_reference_system | sys, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec, | ||
double * | dis, | ||
double * | rv | ||
) |
Computes the place of a solar system body at the specified time for an observer in the specified coordinate system. This is the same as calling place() with the same arguments, except the different set of return values used.
REFERENCES:
jd_tt | [day] Terretrial Time (TT) based Julian date. | |
ss_body | Pointer to structure containing the body designation for the solar system body. | |
obs | Observer location. It may be NULL if not relevant. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
sys | Coordinate reference system in which to produce output values | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Topocentric apparent right ascension in hours, referred to the true equator and equinox of date, or NAN when returning with an error code. (It may be NULL if not required) |
[out] | dec | [deg] Topocentric apparent declination in degrees referred to the true equator and equinox of date, or NAN when returning with an error code. (It may be NULL if not required) |
[out] | dis | [AU] True distance from Earth to the body at 'jd_tt' in AU, or NAN when returning with an error code. (It may be NULL if not needed). |
[out] | rv | [AU/day] radial velocity relative ot observer, or NAN when returning with an error code. (It may be NULL if not required) |
References sky_pos::dec, sky_pos::dis, NOVAS_EPHEM_OBJECT, NOVAS_ORBITAL_OBJECT, NOVAS_PLANET, place(), sky_pos::ra, sky_pos::rv, SKY_POS_INIT, and object::type.
int radec_star | ( | double | jd_tt, |
const cat_entry * | star, | ||
const observer * | obs, | ||
double | ut1_to_tt, | ||
enum novas_reference_system | sys, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec, | ||
double * | rv | ||
) |
Computes the place of a star at date 'jd_tt', for an observer in the specified coordinate system, given the star's ICRS catalog place, proper motion, parallax, and radial velocity.
Notwithstanding the different set of return values, this is the same as calling place_star() with the same arguments.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
obs | Observer location. It may be NULL if not relevant. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
sys | Coordinate reference system in which to produce output values | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Topocentric right ascension in hours, referred to true equator and equinox of date 'jd_tt' or NAN when returning with an error code. (It may be NULL if not required) |
[out] | dec | [deg] Topocentric declination in degrees, referred to true equator and equinox of date 'jd_tt' or NAN when returning with an error code. (It may be NULL if not required) |
[out] | rv | [AU/day] radial velocity relative ot observer, or NAN when returning with an error code. (It may be NULL if not required) |
References sky_pos::dec, place_star(), sky_pos::ra, sky_pos::rv, and SKY_POS_INIT.
double refract | ( | const on_surface * | location, |
enum novas_refraction_model | option, | ||
double | zd_obs | ||
) |
Computes atmospheric optical refraction for an observed (already refracted!) zenith distance through the atmosphere. In other words this is suitable to convert refracted zenith angles to astrometric (unrefracted) zenith angles. For the reverse, see refract_astro().
The returned value is the approximate refraction for optical wavelengths. This function can be used for planning observations or telescope pointing, but should not be used for precise positioning.
NOTES:
REFERENCES:
location | Pointer to structure containing observer's location. It may also contains weather data (optional) for the observer's location. |
option | NOVAS_STANDARD_ATMOSPHERE (1), or NOVAS_WEATHER_AT_LOCATION (2) if to use the weather values contained in the 'location' data structure. |
zd_obs | [deg] Observed (already refracted!) zenith distance through the atmosphere. |
References on_surface::height, on_surface::latitude, NOVAS_NO_ATMOSPHERE, NOVAS_STANDARD_ATMOSPHERE, NOVAS_WEATHER_AT_LOCATION, on_surface::pressure, and on_surface::temperature.
double refract_astro | ( | const on_surface * | location, |
enum novas_refraction_model | option, | ||
double | zd_astro | ||
) |
Computes atmospheric optical refraction for a source at an astrometric zenith distance (e.g. calculated without accounting for an atmosphere). This is suitable for converting astrometric (unrefracted) zenith angles to observed (refracted) zenith angles. See refract() for the reverse correction.
The returned value is the approximate refraction for optical wavelengths. This function can be used for planning observations or telescope pointing, but should not be used for precise positioning.
REFERENCES:
location | Pointer to structure containing observer's location. It may also contains weather data (optional) for the observer's location. |
option | NOVAS_STANDARD_ATMOSPHERE (1), or NOVAS_WEATHER_AT_LOCATION (2) if to use the weather values contained in the 'location' data structure. |
zd_astro | [deg] Astrometric (unrefracted) zenith distance angle of the source. |
References novas_inv_max_iter, and refract().
int set_cio_locator_file | ( | const char * | filename | ) |
Sets the CIO interpolaton data file to use to interpolate CIO locations vs the GCRS. You can specify either the original CIO_RA.TXT
file included in the distribution (preferred since v1.1), or else a platform-specific binary data file compiled from it via the cio_file
utility (the old way).
filename | Path (preferably absolute path) CIO_RA.TXT or else to the binary cio_ra.bin data. |
int set_ephem_provider | ( | novas_ephem_provider | func | ) |
Sets the function to use for obtaining position / velocity information for minor planets, or sattelites.
func | new function to use for accessing ephemeris data for minor planets or satellites. |
int set_nutation_lp_provider | ( | novas_nutation_provider | func | ) |
Set the function to use for low-precision IAU 2000 nutation calculations instead of the default nu2000k().
func | the new function to use for low-precision IAU 2000 nutation calculations |
int set_planet_provider | ( | novas_planet_provider | func | ) |
Set a custom function to use for regular precision (see NOVAS_REDUCED_ACCURACY) ephemeris calculations instead of the default solarsystem() routine.
func | The function to use for solar system position/velocity calculations. See solarsystem() for further details on what is required of this function. |
int set_planet_provider_hp | ( | novas_planet_provider_hp | func | ) |
Set a custom function to use for high precision (see NOVAS_FULL_ACCURACY) ephemeris calculations instead of the default solarsystem_hp() routine.
func | The function to use for solar system position/velocity calculations. See solarsystem_hp() for further details on what is required of this function. |
short sidereal_time | ( | double | jd_ut1_high, |
double | jd_ut1_low, | ||
double | ut1_to_tt, | ||
enum novas_equinox_type | gst_type, | ||
enum novas_earth_rotation_measure | erot, | ||
enum novas_accuracy | accuracy, | ||
double * | gst | ||
) |
Computes the Greenwich sidereal time, either mean or apparent, at the specified Julian date. The Julian date can be broken into two parts if convenient, but for the highest precision, set 'jd_high' to be the integral part of the Julian date, and set 'jd_low' to be the fractional part.
NOTES:
REFERENCES:
jd_ut1_high | [day] High-order part of UT1 Julian date. | |
jd_ut1_low | [day] Low-order part of UT1 Julian date. (You can leave it at zero if 'jd_high' specified the date with sufficient precision) | |
ut1_to_tt | [s] TT - UT1 Time difference in seconds | |
gst_type | NOVAS_MEAN_EQUINOX (0) or NOVAS_TRUE_EQUINOX (1), depending on whether wanting mean or apparent GST, respectively. | |
erot | EROT_ERA (0) or EROT_GST (1), depending on whether to use GST relative to equinox of date (pre IAU 2006) or ERA relative to the CIO (IAU 2006 standard). | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | gst | [h] Greenwich (mean or apparent) sidereal time, in hours [0:24]. (In case the returned error code is >1 the gst value will be set to NAN.) |
References cio_basis(), cio_location(), e_tilt(), era(), EROT_ERA, EROT_GST, NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, NOVAS_TRUE_EQUINOX, tod_to_gcrs(), and tt2tdb().
int spin | ( | double | angle, |
const double * | in, | ||
double * | out | ||
) |
Transforms a vector from one coordinate system to another with same origin and axes rotated about the z-axis.
REFERENCES:
angle | [deg] Angle of coordinate system rotation, positive counterclockwise when viewed from +z, in degrees. | |
in | Input position vector. | |
[out] | out | Position vector expressed in new coordinate system rotated about z by 'angle'. It can be the same vector as the input. |
References TWOPI.
int starvectors | ( | const cat_entry * | star, |
double * | pos, | ||
double * | vel | ||
) |
Converts angular quantities for stars to vectors.
REFERENCES:
star | Pointer to catalog entry structure containing ICRS catalog | |
[out] | pos | [AU] Position vector, equatorial rectangular coordinates, components in AU. It may be NULL if not required. |
[out] | vel | [AU/day] Velocity vector, equatorial rectangular coordinates, components in AU/Day. It must be distinct from the pos output vector, and may be NULL if not required. |
References AU, C, cat_entry::dec, NOVAS_KM, cat_entry::parallax, cat_entry::promodec, cat_entry::promora, cat_entry::ra, and cat_entry::radialvelocity.
int tdb2tt | ( | double | jd_tdb, |
double * | jd_tt, | ||
double * | secdiff | ||
) |
Computes the Terrestrial Time (TT) or Terrestrial Dynamical Time (TDT) Julian date corresponding to a Barycentric Dynamical Time (TDB) Julian date.
Expression used in this function is a truncated form of a longer and more precise series given in the first reference. The result is good to about 10 microseconds.
REFERENCES:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
[out] | jd_tt | [day] Terrestrial Time (TT) based Julian date. (It may be NULL if not required) |
[out] | secdiff | [s] Difference 'tdb_jd'-'tt_jd', in seconds. (It may be NULL if not required) |
short ter2cel | ( | double | jd_ut1_high, |
double | jd_ut1_low, | ||
double | ut1_to_tt, | ||
enum novas_earth_rotation_measure | erot, | ||
enum novas_accuracy | accuracy, | ||
enum novas_equatorial_class | class, | ||
double | xp, | ||
double | yp, | ||
const double * | in, | ||
double * | out | ||
) |
Rotates a vector from the terrestrial to the celestial system. Specifically, it transforms a vector in the ITRS (rotating earth-fixed system) to the True of Date (TOD), CIRS, or GCRS (a local space-fixed system) by applying rotations for polar motion, Earth rotation (for TOD); and nutation, precession, and the dynamical-to-GCRS frame tie (for GCRS).
If 'system' is NOVAS_CIRS then method EROT_ERA must be used. Similarly, if 'system' is NOVAS_TOD then method must be EROT_ERA. Otherwise an error 3 is returned.
If both 'xp' and 'yp' are set to 0 no polar motion is included in the transformation.
REFERENCES:
jd_ut1_high | [day] High-order part of UT1 Julian date. | |
jd_ut1_low | [day] Low-order part of UT1 Julian date. | |
ut1_to_tt | [s] TT - UT1 Time difference in seconds | |
erot | EROT_ERA (0) or EROT_GST (1), depending on whether to use GST relative to equinox of date (pre IAU 2006) or ERA relative to the CIO (IAU 2006 standard) as the Earth rotation measure. The main effect of this option is that it selects the output coordinate system as CIRS or TOD if the output coordinate class is NOVAS_DYNAMICAL_CLASS. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
class | Output coordinate class NOVAS_REFERENCE_CLASS (0, or any value other than 1) or NOVAS_DYNAMICAL_CLASS (1). Use the former if the output coordinates are to be in the GCRS, and the latter if they are to be in CIRS or TOD (the 'erot' parameter selects which dynamical system to use for the output.) | |
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) in the normal case where 'option' is NOVAS_GCRS (0). | |
[out] | out | Position vector, equatorial rectangular coordinates in the specified output system (GCRS if 'class' is NOVAS_REFERENCE_CLASS; or else either CIRS if 'erot' is EROT_ERA, or TOD if 'erot' is EROT_GST). It may be the same vector as the input. |
References cirs_to_gcrs(), era(), EROT_ERA, EROT_GST, NOVAS_DYNAMICAL_CLASS, NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, NOVAS_TRUE_EQUINOX, sidereal_time(), spin(), tod_to_gcrs(), tt2tdb(), wobble(), and WOBBLE_ITRS_TO_PEF.
int terra | ( | const on_surface * | location, |
double | lst, | ||
double * | pos, | ||
double * | vel | ||
) |
Computes the position and velocity vectors of a terrestrial observer with respect to the center of the Earth.
This function ignores polar motion, unless the observer's longitude and latitude have been corrected for it, and variation in the length of day (angular velocity of earth).
The true equator and equinox of date do not form an inertial system. Therefore, with respect to an inertial system, the very small velocity component (several meters/day) due to the precession and nutation of the Earth's axis is not accounted for here.
REFERENCES:
location | Location of observer in Earth's rotating frame | |
lst | [h] Local apparent sidereal time at reference meridian in hours. | |
[out] | pos | [AU] Position vector of observer with respect to center of Earth, equatorial rectangular coordinates, referred to true equator and equinox of date, components in AU. If reference meridian is Greenwich and 'lst' = 0, 'pos' is effectively referred to equator and Greenwich. (It may be NULL if no position data is required). |
[out] | vel | [AU/day] Velocity vector of observer with respect to center of Earth, equatorial rectangular coordinates, referred to true equator and equinox of date, components in AU/day. (It must be distinct from the pos output vector, and may be NULL if no velocity data is required). |
References ANGVEL, AU_KM, ERAD, on_surface::height, on_surface::latitude, on_surface::longitude, and NOVAS_KM.
int tod_to_gcrs | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
const double * | in, | ||
double * | out | ||
) |
Transforms a rectangular equatorial (x, y, z) vector from True of Date (TOD) reference frame at the given epoch to the Geocentric Celestial Reference System(GCRS)
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date that defines the input 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 | Input (x, y, z) position or velocity 3-vector in the True equinox of Date coordinate frame. | |
[out] | out | Output GCRS position or velocity vector. It can be the same vector as the input. |
References frame_tie(), J2000_TO_ICRS, and tod_to_j2000().
int tod_to_j2000 | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
const double * | in, | ||
double * | out | ||
) |
Transforms a rectangular equatorial (x, y, z) vector from True of Date (TOD) reference frame at the given epoch to the J2000 coordinates.
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date that defines the input 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 | Input (x, y, z) position or velocity 3-vector in the True equinox of Date coordinate frame. | |
[out] | out | Output position or velocity vector in rectangular equatorial coordinates at J2000. It can be the same vector as the input. |
References NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, NUTATE_TRUE_TO_MEAN, nutation(), and precession().
short topo_planet | ( | double | jd_tt, |
const object * | ss_body, | ||
double | ut1_to_tt, | ||
const on_surface * | position, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec, | ||
double * | dis | ||
) |
Computes the topocentric apparent place of a solar system body at the specified time. This is the same as calling place() for the body for the same observer location and NOVAS_TOD as the reference system, except the different set of return values used.
REFERENCES:
jd_tt | [day] Terretrial Time (TT) based Julian date. | |
ss_body | Pointer to structure containing the body designation for the solar system body. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
position | Position of the observer | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Topocentric apparent right ascension in hours, referred to the true equator and equinox of date. (It may be NULL if not required) |
[out] | dec | [deg] Topocentric apparent declination in degrees referred to the true equator and equinox of date. (It may be NULL if not required) |
[out] | dis | [AU] True distance from Earth to the body at 'jd_tt' in AU (may be NULL if not needed). |
References make_observer(), NOVAS_OBSERVER_ON_EARTH, NOVAS_TOD, and radec_planet().
short topo_star | ( | double | jd_tt, |
double | ut1_to_tt, | ||
const cat_entry * | star, | ||
const on_surface * | position, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec | ||
) |
Computes the topocentric (True of Date; TOD) apparent place of a star at date 'jd_tt', given its ICRS catalog place, proper motion, parallax, and radial velocity.
Notwithstanding the different set of return values, this is the same as calling place_star() with the same observer location and NOVAS_TOD for an object that specifies the star.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
position | Position of the observer | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Topocentric right ascension in hours, referred to true equator and equinox of date 'jd_tt'. (It may be NULL if not required) |
[out] | dec | [deg] Topocentric declination in degrees, referred to true equator and equinox of date 'jd_tt'. (It may be NULL if not required) |
References make_observer(), NOVAS_OBSERVER_ON_EARTH, NOVAS_TOD, and radec_star().
short transform_cat | ( | enum novas_transform_type | option, |
double | jd_tt_in, | ||
const cat_entry * | in, | ||
double | jd_tt_out, | ||
const char * | out_id, | ||
cat_entry * | out | ||
) |
Transform a star's catalog quantities for a change the coordinate system and/or the date for which the positions are calculated. Also used to rotate catalog quantities on the dynamical equator and equinox of J2000.0 to the ICRS or vice versa.
'date_incat' and 'date_newcat' may be specified either as a Julian date (e.g., 2433282.5 or NOVAS_JD_B1950) or a fractional Julian year and fraction (e.g., 1950.0). Values less than 10000 are assumed to be years. You can also use the supplied constants NOVAS_JD_J2000 or NOVAS_JD_B1950. The date arguments are ignored for the ICRS frame conversion options.
If 'option' is PROPER_MOTION (1), input data can be in any reference system. If 'option' is PRECESSION (2) or CHANGE_EPOCH (3), input data is assume to be in the dynamical system of 'date_incat' and produces output in the dynamical system of 'date_outcat'. If 'option' is CHANGE_J2000_TO_ICRS (4), the input data should be in the J2000.0 dynamical frame. And if 'option' is CHANGE_ICRS_TO_J2000 (5), the input data must be in the ICRS, and the output will be in the J2000 dynamical frame.
This function cannot be properly used to bring data from old star catalogs into the modern system, because old catalogs were compiled using a set of constants that are incompatible with modern values. In particular, it should not be used for catalogs whose positions and proper motions were derived by assuming a precession constant significantly different from the value implicit in function precession().
option | Type of transformation | |
jd_tt_in | [day|yr] Terrestrial Time (TT) based Julian date, or year, of input catalog data. Not used if option is CHANGE_J2000_TO_ICRS (4) or CHANGE_ICRS_TO_J2000 (5). | |
in | An entry from the input catalog, with units as given in the struct definition | |
jd_tt_out | [day|yr] Terrestrial Time (TT) based Julian date, or year, of output catalog data. Not used if option is CHANGE_J2000_TO_ICRS (4) or CHANGE_ICRS_TO_J2000 (5). | |
out_id | Catalog identifier (0 terminated). It may also be NULL in which case the catalog name is inherited from the input. | |
[out] | out | The transformed catalog entry, with units as given in the struct definition |
References AU_KM, C, cat_entry::catalog, CHANGE_EPOCH, CHANGE_ICRS_TO_J2000, CHANGE_J2000_TO_ICRS, cat_entry::dec, frame_tie(), ICRS_TO_J2000, J2000_TO_ICRS, NOVAS_JD_J2000, NOVAS_KM, cat_entry::parallax, PRECESSION, precession(), cat_entry::promodec, cat_entry::promora, PROPER_MOTION, cat_entry::ra, cat_entry::radialvelocity, cat_entry::starname, and cat_entry::starnumber.
Convert Hipparcos catalog data at epoch J1991.25 to epoch J2000.0, for use within NOVAS. To be used only for Hipparcos or Tycho stars with linear space motion. Both input and output data is in the ICRS.
hipparcos | An entry from the Hipparcos catalog, at epoch J1991.25, with 'ra' in degrees(!) as per Hipparcos catalog units. | |
[out] | hip_2000 | The transformed input entry, at epoch J2000.0, with 'ra' in hours(!) as per the NOVAS convention. |
References cat_entry::catalog, NOVAS_JD_HIP, cat_entry::ra, and transform_cat().
double tt2tdb | ( | double | jd_tt | ) |
Returns the TDB - TT time difference in seconds for a given TT date.
Note, as of version 1.1, it uses the same calculation as the more precise original tdb2tt(). It thus has an acuracy of about 10 μs vs around 30 μs with the simpler formula from the references below.
REFERENCES
jd_tt | [day] Terrestrial Time (TT) based Julian date |
References tdb2tt().
short vector2radec | ( | const double * | pos, |
double * | ra, | ||
double * | dec | ||
) |
Converts an vector in equatorial rectangular coordinates to equatorial spherical coordinates.
REFERENCES:
pos | Position 3-vector, equatorial rectangular coordinates. | |
[out] | ra | [h] Right ascension in hours [0:24] or NAN if the position vector is NULL or a null-vector. It may be NULL if notrequired. |
[out] | dec | [deg] Declination in degrees [-90:90] or NAN if the position vector is NULL or a null-vector. It may be NULL if not required. |
short virtual_planet | ( | double | jd_tt, |
const object * | ss_body, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec, | ||
double * | dis | ||
) |
Computes the virtual place of a solar system body, referenced to the GCRS. This is the same as calling place_gcrs() for the body, except the different set of return values used.
REFERENCES:
jd_tt | [day] Terretrial Time (TT) based Julian date. | |
ss_body | Pointer to structure containing the body designation for the solar system body. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Virtual right ascension in hours, referred to the GCRS (it may be NULL if not required). |
[out] | dec | [deg] Virtual declination in degrees, referred to the GCRS (it may be NULL if not required). |
[out] | dis | [AU] True distance from Earth to the body at 'jd_tt' in AU (can be NULL if not needed). |
References NOVAS_GCRS, and radec_planet().
short virtual_star | ( | double | jd_tt, |
const cat_entry * | star, | ||
enum novas_accuracy | accuracy, | ||
double * | ra, | ||
double * | dec | ||
) |
Computes the virtual place of a star, referenced to GCRS, at date 'jd_tt', given its catalog mean place, proper motion, parallax, and radial velocity.
Notwithstanding the different set of return values, this is the same as calling place_star() with a NULL observer location and NOVAS_GCRS as the system, or place_gcrs() for an object that specifies the star.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Virtual right ascension in hours, referred to the GCRS (it may be NULL if not required). |
[out] | dec | [deg] Virtual declination in degrees, referred to the GCRS (it may be NULL if not required). |
References NOVAS_GCRS, and radec_star().
int wobble | ( | double | jd_tt, |
enum novas_wobble_direction | direction, | ||
double | xp, | ||
double | yp, | ||
const double * | in, | ||
double * | out | ||
) |
Corrects a vector in the ITRS (rotating Earth-fixed system) for polar motion, and also corrects the longitude origin (by a tiny amount) to the Terrestrial Intermediate Origin (TIO). The ITRS vector is thereby transformed to the terrestrial intermediate system, based on the true (rotational) equator and TIO. Because the true equator is the plane orthogonal to the direction of the Celestial Intermediate Pole (CIP), the components of the output vector are referred to z and x axes toward the CIP and TIO, respectively.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
direction | WOBBLE_ITRS_TO_PEF (0) or WOBBLE_PEF_TO_ITRS (1; or nonzero) | |
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 | Input position vector, geocentric equatorial rectangular coordinates, in the original system defined by 'direction' | |
[out] | out | Output Position vector, geocentric equatorial rectangular coordinates, in the final system defined by 'direction'. It can be the same vector as the input. |
References WOBBLE_ITRS_TO_PEF.
double EPS_COR = 0.0 |
Celestial pole offset ε for high-precision applications. It was visible to users in NOVAS C 3.1, hence we continue to expose it also for back compatibility.
int grav_bodies_full_accuracy = DEFAULT_GRAV_BODIES_FULL_ACCURACY |
Current set of gravitating bodies to use for deflection calculations in full accuracy mode. Each bit signifies whether a given body is to be accounted for as a gravitating body that bends light, such as the bit (1 << NOVAS_JUPITER)
indicates whether or not Jupiter is considered as a deflecting body. You should also be sure that you provide ephemeris data for bodies that are designated for the deflection calculation.
int grav_bodies_reduced_accuracy = DEFAULT_GRAV_BODIES_REDUCED_ACCURACY |
Current set of gravitating bodies to use for deflection calculations in reduced accuracy mode. Each bit signifies whether a given body is to be accounted for as a gravitating body that bends light, such as the bit (1 << NOVAS_JUPITER)
indicates whether or not Jupiter is considered as a deflecting body. You should also be sure that you provide ephemeris data for bodies that are designated for the deflection calculation.
int novas_inv_max_iter = 100 |
Maximum number of iterations for convergent inverse calculations. Most iterative inverse functions should normally converge in a handful of iterations. In some pathological cases more iterations may be required. This variable sets an absolute maximum for the number of iterations in order to avoid runaway (zombie) behaviour. If inverse functions faile to converge, they will return a value indicating an error, and errno should be set to ECANCELED.
double PSI_COR = 0.0 |
Celestial pole offset ψ for high-precision applications. It was visible to users in NOVAS C 3.1, hence we continue to expose it also for back compatibility.