SuperNOVAS v1.2
The NOVAS C library, made better
|
Macros | |
#define | DA0 (-0.01460 * ARCSEC) |
Frame bias term da0 | |
#define | ETA0 (-0.0068192 * ARCSEC) |
Frame bias term η0 | |
#define | XI0 (-0.0166170 * ARCSEC) |
Frame bias term ξ0 | |
Functions | |
int | novas_app_to_geom (const novas_frame *frame, enum novas_reference_system sys, double ra, double dec, double dist, double *geom_icrs) |
int | novas_app_to_hor (const novas_frame *frame, enum novas_reference_system sys, double ra, double dec, RefractionModel ref_model, double *az, double *el) |
int | novas_change_observer (const novas_frame *orig, const observer *obs, novas_frame *out) |
int | novas_geom_posvel (const object *source, const novas_frame *frame, enum novas_reference_system sys, double *pos, double *vel) |
int | novas_geom_to_app (const novas_frame *frame, const double *pos, enum novas_reference_system sys, sky_pos *out) |
int | novas_hor_to_app (const novas_frame *frame, double az, double el, RefractionModel ref_model, enum novas_reference_system sys, double *ra, double *dec) |
int | novas_invert_transform (const novas_transform *transform, novas_transform *inverse) |
int | novas_make_frame (enum novas_accuracy accuracy, const observer *obs, const novas_timespec *time, double dx, double dy, novas_frame *frame) |
int | novas_make_transform (const novas_frame *frame, enum novas_reference_system from_system, enum novas_reference_system to_system, novas_transform *transform) |
int | novas_sky_pos (const object *object, const novas_frame *frame, enum novas_reference_system sys, sky_pos *out) |
int | novas_transform_sky_pos (const sky_pos *in, const novas_transform *transform, sky_pos *out) |
int | novas_transform_vector (const double *in, const novas_transform *transform, double *out) |
SuperNOVAS routines for higher-level and efficient repeat coordinate transformations using observer frames. Observer frames represent an observer location at a specific astronomical time (instant), which can be re-used again and again to calculate or transform positions of celestial sources in a a range of astronomical coordinate systems.
int novas_app_to_geom | ( | const novas_frame * | frame, |
enum novas_reference_system | sys, | ||
double | ra, | ||
double | dec, | ||
double | dist, | ||
double * | geom_icrs | ||
) |
Converts an observed apparent sky position of a source to an ICRS geometric position, by undoing the gravitational deflection and aberration corrections.
frame | The observer frame, defining the location and time of observation | |
sys | The reference system in which the observed position is specified. | |
ra | [h] Observed ICRS right-ascension of the source | |
dec | [deg] Observed ICRS declination of the source | |
dist | [AU] Observed distance from observer. A value of <=0 will translate to 1015 AU (around 5 Gpc). | |
[out] | geom_icrs | [AU] The corresponding geometric position for the source, in ICRS. |
References novas_frame::gcrs_to_cirs, grav_undo_planets(), novas_frame::icrs_to_j2000, NOVAS_CIRS, NOVAS_J2000, NOVAS_MOD, NOVAS_REFERENCE_SYSTEMS, NOVAS_TOD, novas_frame::nutation, novas_frame::obs_pos, novas_frame::planets, novas_frame::precession, and radec2vector().
int novas_app_to_hor | ( | const novas_frame * | frame, |
enum novas_reference_system | sys, | ||
double | ra, | ||
double | dec, | ||
RefractionModel | ref_model, | ||
double * | az, | ||
double * | el | ||
) |
Converts an observed apparent position vector in the specified coordinate system to local horizontal coordinates in the specified observer frame. The observer must be located on the surface of Earth, or else the call will return with an error. The caller may optionally supply a refraction model of choice to calculate an appropriate elevation angle that includes a refraction correction for Earth's atmosphere. If no such model is provided the calculated elevation will be the astrometric elevation without a refraction correction.
frame | Observer frame, defining the time and place of observation (on Earth). | |
sys | Astronomical coordinate system in which the observed position is given. | |
ra | [h] Observed apparent right ascension (R.A.) coordinate | |
dec | [deg] Observed apparent declination coordinate | |
ref_model | An appropriate refraction model, or NULL to calculate unrefracted elevation. Depending on the refraction model, you might want to make sure that the weather parameters were set when the observing frame was defined. | |
[out] | az | [deg] Calculated azimuth angle. It may be NULL if not required. |
[out] | el | [deg] Calculated elevation angle. It may be NULL if not required. |
References novas_frame::era, novas_timespec::fjd_tt, novas_frame::gcrs_to_cirs, novas_frame::gst, novas_timespec::ijd_tt, itrs_to_hor(), NOVAS_AIRBORNE_OBSERVER, NOVAS_CIRS, NOVAS_GCRS, NOVAS_ICRS, NOVAS_J2000, NOVAS_MOD, NOVAS_OBSERVER_ON_EARTH, NOVAS_REFRACT_ASTROMETRIC, NOVAS_TOD, novas_frame::nutation, novas_frame::observer, observer::on_surf, novas_frame::precession, radec2vector(), spin(), novas_frame::time, and observer::where.
int novas_change_observer | ( | const novas_frame * | orig, |
const observer * | obs, | ||
novas_frame * | out | ||
) |
Change the observer location for an observing frame.
orig | Pointer to original observing frame | |
obs | New observer location | |
[out] | out | Observing frame to populate with a original frame data and new observer location. It can be the same as the input. |
References novas_frame::accuracy, grav_bodies_full_accuracy, grav_bodies_reduced_accuracy, NOVAS_FULL_ACCURACY, novas_get_time(), NOVAS_TDB, obs_planets(), novas_frame::obs_pos, novas_frame::observer, novas_frame::planets, novas_frame::state, and novas_frame::time.
int novas_geom_posvel | ( | const object * | source, |
const novas_frame * | frame, | ||
enum novas_reference_system | sys, | ||
double * | pos, | ||
double * | vel | ||
) |
Calculates the geometric position and velocity vectors, relative to the observer, for a source in the given observing frame, in the specified coordinate system of choice. The geometric position includes proper motion, and for solar-system bodies it is antedated for light travel time, so it effectively represents the geometric position as seen by the observer. However, the geometric does not include aberration correction, nor gravitational deflection.
If you want apparent positions, which account for aberration and gravitational deflection, use novas_skypos() instead.
You can also use novas_transform_vector() to convert the output position and velocity vectors to a dfferent coordinate system of choice afterwards if you want the results expressed in more than one coordinate system.
It implements the same geometric transformations as place()
but at a reduced computational cost. See place()
for references.
source | Pointer to a celestial source data structure that is observed | |
frame | Observer frame, defining the location and time of observation | |
sys | The coordinate system in which to return positions and velocities. | |
[out] | pos | [AU] Calculated geometric position vector of the source relative to the observer location, in the designated coordinate system. It may be NULL if not required. |
[out] | vel | [AU/day] The calculated velocity vector of the source relative to the observer in the designated coordinate system. It must be distinct from the pos output vector, and may be NULL if not required. |
References novas_frame::accuracy, bary2obs(), d_light(), light_time2(), novas_planet_bundle::mask, NOVAS_CATALOG_OBJECT, NOVAS_FULL_ACCURACY, novas_get_time(), NOVAS_JD_J2000, NOVAS_PLANET, NOVAS_REDUCED_ACCURACY, NOVAS_TDB, object::number, novas_frame::obs_pos, novas_frame::planets, novas_planet_bundle::pos, proper_motion(), object::star, starvectors(), novas_frame::time, object::type, and novas_planet_bundle::vel.
int novas_geom_to_app | ( | const novas_frame * | frame, |
const double * | pos, | ||
enum novas_reference_system | sys, | ||
sky_pos * | out | ||
) |
Converts an geometric position in ICRS to an apparent position on sky, by applying appropriate corrections for aberration and gravitational deflection for the observer's frame. Unlike place()
the output reports the distance calculated from the parallax for sidereal sources. The radial velocity of the output is set to NAN (if needed use novas_sky_pos() instead).
frame | The observer frame, defining the location and time of observation | |
pos | [AU] Geometric position of source in ICRS coordinates | |
sys | The coordinate system in which to return the apparent sky location | |
[out] | out | Pointer to the data structure which is populated with the calculated apparent location in the designated coordinate system. It may be the same pounter as the input position. |
References novas_frame::accuracy, sky_pos::dec, sky_pos::dis, grav_planets(), NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, novas_frame::obs_pos, novas_frame::planets, sky_pos::r_hat, sky_pos::ra, sky_pos::rv, and vector2radec().
int novas_hor_to_app | ( | const novas_frame * | frame, |
double | az, | ||
double | el, | ||
RefractionModel | ref_model, | ||
enum novas_reference_system | sys, | ||
double * | ra, | ||
double * | dec | ||
) |
Converts an observed azimuth and elevation coordinate to right ascension (R.A.) and declination coordinates expressed in the coordinate system of choice. The observer must be located on the surface of Earth, or else the call will return with an error. The caller may optionally supply a refraction model of choice to calculate an appropriate elevation angle that includes a refraction correction for Earth's atmosphere. If no such model is provided, the provided elevation value will be assumed to be an astrometric elevation without a refraction correction.
frame | Observer frame, defining the time and place of observation (on Earth). | |
az | [deg] Observed azimuth angle. It may be NULL if not required. | |
el | [deg] Observed elevation angle. It may be NULL if not required. | |
ref_model | An appropriate refraction model, or NULL to assume unrefracted elevation. Depending on the refraction model, you might want to make sure that the weather parameters were set when the observing frame was defined. | |
sys | Astronomical coordinate system in which the output is R.A. and declination values are to be calculated. | |
[out] | ra | [h] Calculated apparent right ascension (R.A.) coordinate |
[out] | dec | [deg] Calculated apparent declination coordinate |
References novas_frame::era, novas_timespec::fjd_tt, novas_frame::gcrs_to_cirs, novas_frame::gst, hor_to_itrs(), novas_timespec::ijd_tt, NOVAS_AIRBORNE_OBSERVER, NOVAS_CIRS, NOVAS_GCRS, NOVAS_ICRS, NOVAS_J2000, NOVAS_MOD, NOVAS_OBSERVER_ON_EARTH, NOVAS_REFRACT_OBSERVED, NOVAS_TOD, novas_frame::nutation, novas_frame::observer, observer::on_surf, novas_frame::precession, spin(), novas_frame::time, vector2radec(), and observer::where.
int novas_invert_transform | ( | const novas_transform * | transform, |
novas_transform * | inverse | ||
) |
Inverts a novas coordinate transformation matrix.
transform | Pointer to a coordinate transformation matrix. | |
[out] | inverse | Pointer to a coordinate transformation matrix to populate with the inverse transform. It may be the same as the input. |
References novas_transform::matrix.
int novas_make_frame | ( | enum novas_accuracy | accuracy, |
const observer * | obs, | ||
const novas_timespec * | time, | ||
double | dx, | ||
double | dy, | ||
novas_frame * | frame | ||
) |
Sets up a observing frame for a specific observer location, time of observation, and accuracy requirement. The frame is initialized using the currently configured planet ephemeris provider function (see set_planet_provider() and set_planet_provider_hp()), and in case of reduced accuracy mode, the currently configured IAU nutation model provider (see set_nutation_lp_provider()).
Note, that to construct full accuracy frames, you will need a high-precision ephemeris provider for the major planets (not just the default Earth/Sun), as without it, gravitational bending around massive plannets cannot be accounted for, and therefore μas accuracy cannot be ensured, in general. Attempting to construct a high-accuracy frame without a high-precision ephemeris provider for the major planets will result in an error in the 10–40 range from the required ephemeris() call.
accuracy | Accuracy requirement, NOVAS_FULL_ACCURACY (0) for the utmost precision or NOVAS_REDUCED_ACCURACY (1) if ~1 mas accuracy is sufficient. | |
obs | Observer location | |
time | Time of observation | |
dx | [mas] Earth orientation parameter, polar offset in x. | |
dy | [mas] Earth orientation parameter, polar offset in y. | |
[out] | frame | Pointer to the observing frame to configure. |
References novas_frame::accuracy, CAT_ENTRY_INIT, novas_frame::deps0, novas_frame::dpsi0, novas_frame::dx, novas_frame::dy, e_tilt(), novas_frame::earth_pos, novas_frame::earth_vel, novas_frame::ee, ephemeris(), novas_frame::era, era(), EROT_GST, novas_timespec::fjd_tt, novas_frame::gst, novas_timespec::ijd_tt, novas_frame::mobl, NOVAS_BARYCENTER, novas_change_observer(), NOVAS_EARTH, novas_get_split_time(), NOVAS_JD_J2000, NOVAS_OBSERVER_PLACES, NOVAS_PLANET, NOVAS_REDUCED_ACCURACY, NOVAS_SUN, NOVAS_TRUE_EQUINOX, NOVAS_UT1, nutation_angles(), sidereal_time(), novas_frame::state, novas_frame::sun_pos, novas_frame::sun_vel, novas_frame::time, novas_frame::tobl, tt2tdb(), novas_timespec::ut1_to_tt, and observer::where.
int novas_make_transform | ( | const novas_frame * | frame, |
enum novas_reference_system | from_system, | ||
enum novas_reference_system | to_system, | ||
novas_transform * | transform | ||
) |
Calculates a transformation matrix that can be used to convert positions and velocities from one coordinate reference system to another.
frame | Observer frame, defining the location and time of observation | |
from_system | Original coordinate reference system | |
to_system | New coordinate reference system | |
[out] | transform | Pointer to the transform data structure to populate. |
References novas_transform::frame, novas_transform::from_system, novas_frame::gcrs_to_cirs, novas_frame::icrs_to_j2000, novas_matrix::M, novas_transform::matrix, NOVAS_CIRS, NOVAS_GCRS, NOVAS_ICRS, NOVAS_J2000, NOVAS_MOD, NOVAS_REFERENCE_SYSTEMS, NOVAS_TOD, novas_frame::nutation, novas_frame::precession, and novas_transform::to_system.
int novas_sky_pos | ( | const object * | object, |
const novas_frame * | frame, | ||
enum novas_reference_system | sys, | ||
sky_pos * | out | ||
) |
Calculates an apparent location on sky for the source. The position takes into account the proper motion (for sidereal soure), or is antedated for light-travel time (for Solar-System bodies). It also applies an appropriate aberration correction and gravitational deflection of the light.
To calculate corresponding local horizontal coordinates, you can pass the output RA/Dec coordinates to novas_app_to_hor(). Or to calculate apparent coordinates in other systems, you may pass the result to novas_transform_sy_pos() after.
And if you want geometric positions instead (not corrected for aberration or gravitational deflection), you may want to use novas_geom_posvel() instead.
The approximate 'inverse' of this function is novas_app_to_geom().
This function implements the same aberration and gravitational deflection corrections as place()
, but at reduced computational cost. See place()
for references. Unlike place()
, however, the output always reports the distance calculated from the parallax for sidereal sources. Note also, that while place()
does not apply aberration and gravitational deflection corrections when sys
is NOVAS_ICRS (3), this routine will apply those corrections consistently for all coordinate systems (and you can use novas_geom_posvel() instead to get positions without aberration or deflection in any system).
object | Pointer to a celestial object data structure that is observed | |
frame | The observer frame, defining the location and time of observation | |
sys | The coordinate system in which to return the apparent sky location | |
[out] | out | Pointer to the data structure which is populated with the calculated apparent location in the designated coordinate system. |
References novas_frame::accuracy, sky_pos::dis, novas_frame::earth_pos, grav_planets(), NOVAS_CATALOG_OBJECT, NOVAS_FULL_ACCURACY, novas_geom_posvel(), novas_geom_to_app(), NOVAS_ICRS, NOVAS_REDUCED_ACCURACY, novas_frame::obs_pos, novas_frame::obs_vel, novas_frame::planets, rad_vel2(), sky_pos::rv, novas_frame::sun_pos, and object::type.
int novas_transform_sky_pos | ( | const sky_pos * | in, |
const novas_transform * | transform, | ||
sky_pos * | out | ||
) |
Transforms a position or velocity 3-vector from one coordinate reference system to another.
in | Input apparent position on sky in the original coordinate reference system | |
transform | Pointer to a coordinate transformation matrix | |
[out] | out | Output apparent position on sky in the new coordinate reference system. It may be the same as the input. |
References sky_pos::dec, novas_transform::matrix, sky_pos::r_hat, sky_pos::ra, and vector2radec().
int novas_transform_vector | ( | const double * | in, |
const novas_transform * | transform, | ||
double * | out | ||
) |
Transforms a position or velocity 3-vector from one coordinate reference system to another.
in | Input 3-vector in the original coordinate reference system | |
transform | Pointer to a coordinate transformation matrix | |
[out] | out | Output 3-vector in the new coordinate reference system. It may be the same as the input. |
References novas_transform::matrix.