![]() |
SuperNOVAS v1.3
The NOVAS C library, made better
|
Functions | |
int | aberration (const double *pos, const double *vobs, double lighttime, double *out) |
int | bary2obs (const double *pos, const double *pos_obs, double *out, double *restrict lighttime) |
short | light_time (double jd_tdb, const object *restrict body, const double *pos_obs, double tlight0, enum novas_accuracy accuracy, double *pos_src_obs, double *restrict tlight) |
int | light_time2 (double jd_tdb, const object *restrict body, const double *restrict pos_obs, double tlight0, enum novas_accuracy accuracy, double *p_src_obs, double *restrict v_ssb, double *restrict tlight) |
int | make_airborne_observer (const on_surface *location, const double *vel, observer *obs) |
int | make_in_space (const double *sc_pos, const double *sc_vel, in_space *loc) |
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 *restrict 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 *restrict obs) |
int | make_on_surface (double latitude, double longitude, double height, double temperature, double pressure, on_surface *restrict loc) |
int | make_solar_system_observer (const double *sc_pos, const double *sc_vel, observer *obs) |
int | novas_e2h_offset (double dra, double ddec, double pa, double *restrict daz, double *restrict del) |
double | novas_epa (double ha, double dec, double lat) |
int | novas_h2e_offset (double daz, double del, double pa, double *restrict dra, double *restrict ddec) |
double | novas_hpa (double az, double el, double lat) |
int | novas_los_to_xyz (const double *los, double lon, double lat, double *xyz) |
int | novas_xyz_to_los (const double *xyz, double lon, double lat, double *los) |
int | novas_xyz_to_uvw (const double *xyz, double ha, double dec, double *uvw) |
int | obs_planets (double jd_tdb, enum novas_accuracy accuracy, const double *restrict pos_obs, int pl_mask, novas_planet_bundle *restrict planets) |
int | obs_posvel (double jd_tdb, double ut1_to_tt, enum novas_accuracy accuracy, const observer *restrict obs, const double *restrict geo_pos, const double *restrict geo_vel, double *restrict pos, double *restrict vel) |
Function that define an astronomical observer location or are related to observer location.
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 novas_vlen().
int bary2obs | ( | const double * | pos, |
const double * | pos_obs, | ||
double * | out, | ||
double *restrict | 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 novas_vlen().
short light_time | ( | double | jd_tdb, |
const object *restrict | body, | ||
const double * | pos_obs, | ||
double | tlight0, | ||
enum novas_accuracy | accuracy, | ||
double * | pos_src_obs, | ||
double *restrict | tlight | ||
) |
Computes the geocentric position of a solar system body, as antedated for light-time.
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 *restrict | body, | ||
const double *restrict | pos_obs, | ||
double | tlight0, | ||
enum novas_accuracy | accuracy, | ||
double * | p_src_obs, | ||
double *restrict | v_ssb, | ||
double *restrict | tlight | ||
) |
Computes the geocentric position and velocity of a solar system body, as antedated for light-time. It is effectively the same as the original NOVAS C light_time(), except that this returns the antedated source velocity vector also.
NOTES:
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 make_airborne_observer | ( | const on_surface * | location, |
const double * | vel, | ||
observer * | obs | ||
) |
Populates an 'observer' data structure for an observer moving relative to the surface of Earth, such as an airborne observer. Airborne observers have an earth fixed momentary location, defined by longitude, latitude, and altitude, the same was as for a stationary observer on Earth, but are moving relative to the surface, such as in an aircraft or balloon observatory.
location | Current longitude, latitude and altitude, and local weather (temperature and pressure) | |
vel | [km/s] Surface velocity. | |
[out] | obs | Pointer to data structure to populate. |
References IN_SPACE_INIT, make_observer(), NOVAS_AIRBORNE_OBSERVER, and in_space::sc_vel.
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_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 *restrict | 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 *restrict | 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 *restrict | 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.
NOTES
humidity
) that was not yet part of the on_surface
structure in v1.0. As such, linking SuperNOVAS v1.1 or later with application code compiled for SuperNOVAS v1.0 can result in memory corruption or segmentation fault when this function is called. To be safe, make sure your application has been (re)compiled against SuperNOVAS v1.1 or later. 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. |
int make_solar_system_observer | ( | const double * | sc_pos, |
const double * | sc_vel, | ||
observer * | obs | ||
) |
Populates an 'observer' data structure, for an observer situated on a near-Earth spacecraft, with the specified geocentric position and velocity vectors. Solar-system observers are similar to observers in Earth-orbit but their momentary position and velocity is defined relative to the Solar System Barycenter, instead of the geocenter.
sc_pos | [AU] Solar-system barycentric (x, y, z) position vector in ICRS. | |
sc_vel | [AU/day] Solar-system barycentric (x, y, z) velocity vector in ICRS. | |
[out] | obs | Pointer to the data structure to populate |
References make_in_space(), make_observer(), and NOVAS_SOLAR_SYSTEM_OBSERVER.
int novas_e2h_offset | ( | double | dra, |
double | ddec, | ||
double | pa, | ||
double *restrict | daz, | ||
double *restrict | del | ||
) |
Converts coordinate offsets, from the local equatorial system to local horizontal offsets. Converting between local flat projections and spherical coordinates usually requires a WCS projection.
REFERENCES:
dra | [arcsec] Projected ffset position in the apparent true-of-date R.A. direction. E.g. The projected offset between two RA coordinates at a same reference declination, is δRA = (RA2 - RA1) * cos(Dec0) | |
ddec | [arcsec] Projected offset position in the apparent true-of-date declination direction. | |
pa | [deg] Parallactic Angle | |
[out] | daz | [arcsec] Output offset position in the local azimuth direction. It can be a pointer to one of the input coordinates, or NULL if not required. |
[out] | del | [arcsec] Output offset position in the local elevation direction. It can be a pointer to one of the input coordinates, or NULL if not required. |
References novas_h2e_offset().
double novas_epa | ( | double | ha, |
double | dec, | ||
double | lat | ||
) |
Returns the equatorial Parallactic Angle (PA) calculated for an R.A./Dec location of the sky at a given sidereal time. The PA is the angle between the local horizontal coordinate directions and the local true-of-date equatorial coordinate directions, at the given location and time. The polar wobble is not included in the calculation.
The Parallactic Angle is sometimes referrred to as the Vertical Position Angle (VPA). Both define the same quantity.
ha | [h] Hour angle (LST - RA) i.e., the difference between the Local (apparent) Sidereal Time and the apparent (true-of-date) Right Ascension of observed source. |
dec | [deg] Apparent (true-of-date) declination of observed source |
lat | [deg] Geodetic latitude of observer |
int novas_h2e_offset | ( | double | daz, |
double | del, | ||
double | pa, | ||
double *restrict | dra, | ||
double *restrict | ddec | ||
) |
Converts coordinate offsets, from the local horizontal system to local equatorial offsets. Converting between local flat projections and spherical coordinates usually requires a WCS projection.
REFERENCES:
daz | [arcsec] Projected offset position in the azimuth direction. The projected offset between two azimuth positions at the same reference elevation is δAz = (Az2 - Az1) * cos(El0). | |
del | [arcsec] projected offset position in the elevation direction | |
pa | [deg] Parallactic Angle | |
[out] | dra | [arcsec] Output offset position in the local true-of-date R.A. direction. It can be a pointer to one of the input coordinates, or NULL if not required. |
[out] | ddec | [arcsec] Output offset position in the local true-of-date declination direction. It can be a pointer to one of the input coordinates, or NULL if not required. |
double novas_hpa | ( | double | az, |
double | el, | ||
double | lat | ||
) |
Returns the horizontal Parallactic Angle (PA) calculated for a gorizontal Az/El location of the sky. The PA is the angle between the local horizontal coordinate directions and the local true-of-date equatorial coordinate directions at the given location. The polar wobble is not included in the calculation.
The Parallactic Angle is sometimes referrred to as the Vertical Position Angle (VPA). Both define the same quantity.
az | [deg] Azimuth angle |
el | [deg] Elevation angle |
lat | [deg] Geodetic latitude of observer |
int novas_los_to_xyz | ( | const double * | los, |
double | lon, | ||
double | lat, | ||
double * | xyz | ||
) |
Converts a 3D line-of-sight vector (δφ, δθ δr) to a rectangular equatorial (δx, δy, δz) vector.
los | [arb.u.] Line-of-sight 3-vector (δφ, δθ δr). | |
lon | [deg] Line-of-sight longitude. | |
lat | [deg] Line-of-sight latitude. | |
[out] | xyz | [arb.u.] Output rectangular equatorial 3-vector (δx, δy, δz), in the same units as the input. It may be the same vector as the input. |
int novas_xyz_to_los | ( | const double * | xyz, |
double | lon, | ||
double | lat, | ||
double * | los | ||
) |
Converts a 3D rectangular equatorial (δx, δy, δz) vector to a polar (δφ, δθ δr) vector along a line-of-sight.
xyz | [arb.u.] Rectangular equatorial 3-vector (δx, δy, δz). | |
lon | [deg] Line-of-sight longitude. | |
lat | [deg] Line-of-sight latitude. | |
[out] | los | [arb.u.] Output line-of-sight 3-vector (δφ, δθ δr), in the same units as the input. It may be the same vector as the input. |
int novas_xyz_to_uvw | ( | const double * | xyz, |
double | ha, | ||
double | dec, | ||
double * | uvw | ||
) |
Converts rectangular telescope x,y,z (absolute or relative) coordinates (in ITRS) to equatorial u,v,w projected coordinates for a specified line of sight.
x,y,z are Cartesian coordinates w.r.t the Greenwich meridian. The directions are x: long=0, lat=0; y: long=90, lat=0; z: lat=90.
u,v,w are Cartesian coordinates (u,v) along the local equatorial R.A. and declination directions as seen from a direction on the sky (w).
xyz | [arb.u.] Absolute or relative x,y,z coordinates (double[3]). | |
ha | [h] Hourangle (LST - RA) i.e., the difference between the Local (apparent) Sidereal Time and the apparent (true-of-date) Right Ascension of observed source. | |
dec | [deg] Apparent (true-of-date) declination of source | |
[out] | uvw | [arb.u.] Converted u,v,w coordinates (double[3]) in same units as xyz. It may be the same vector as the input. |
References novas_xyz_to_los().
int obs_planets | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
const double *restrict | pos_obs, | ||
int | pl_mask, | ||
novas_planet_bundle *restrict | 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_debug(), NOVAS_DEBUG_EXTRA, NOVAS_DEBUG_OFF, novas_get_debug_mode(), NOVAS_PLANETS, and NOVAS_SUN.
int obs_posvel | ( | double | jd_tdb, |
double | ut1_to_tt, | ||
enum novas_accuracy | accuracy, | ||
const observer *restrict | obs, | ||
const double *restrict | geo_pos, | ||
const double *restrict | geo_vel, | ||
double *restrict | pos, | ||
double *restrict | vel | ||
) |
Calculates the ICRS position and velocity of the observer relative to the Solar System Barycenter (SSB).
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(), NOVAS_AIRBORNE_OBSERVER, NOVAS_BARYCENTER, NOVAS_EARTH_INIT, NOVAS_OBSERVER_IN_EARTH_ORBIT, NOVAS_OBSERVER_ON_EARTH, NOVAS_OBSERVER_PLACES, and NOVAS_SOLAR_SYSTEM_OBSERVER.