![]() |
SuperNOVAS C++ API v1.6
High-precision C/C++ astrometry library
|
Functions | |
| int | novas_equals_cat_entry (const cat_entry *a, const cat_entry *b) |
| Checks if two catalog entries define the same parameters, within the typical tolerances associated to these. | |
| int | novas_equals_frame (const novas_frame *a, const novas_frame *b) |
| Checks if two observing frames are essentially the same, within typical tolerances. | |
| int | novas_equals_near_earth (const in_space *a, const in_space *b) |
| Checks if two near-Earth locations and motions match within 1 mm and 1 mm/s, respectively. | |
| int | novas_equals_object (const object *a, const object *b) |
| Checks if two astronomical targets are the same within typical tolerances. | |
| int | novas_equals_observer (const observer *a, const observer *b) |
| Checks if two observers are essentially the same within the tolerances associated to their defining components. | |
| int | novas_equals_on_surface (const on_surface *a, const on_surface *b) |
| Checks if two geodetic locations, and the weather parameters defined for them, match. | |
| int | novas_equals_orbital (const novas_orbital *a, const novas_orbital *b) |
| Checks if two Keplerian orbitals match, within typical tolerances. | |
| int | novas_equals_orbsys (const novas_orbital_system *a, const novas_orbital_system *b) |
| Checks if two orbital systems match within typical tolerances. | |
| int | novas_equals_planet_bundle (const novas_planet_bundle *a, const novas_planet_bundle *b) |
| Checks if two bundles containing Solar-system baricentric planet positions and velocities are effectively the same. | |
| int | novas_equals_sky_pos (const sky_pos *a, const sky_pos *b) |
| Checks if two apparent positions on the observer's sky are the same within typical tolerances. | |
| int | novas_equals_ssb_posvel (const in_space *a, const in_space *b) |
| Checks if two Solar-system locations and motions match within 1 m and ~1 mm/s, respectively. | |
| int | novas_equals_timespec (const novas_timespec *a, const novas_timespec *b) |
| Checks if two time specifications are the same within 10-7 days (~100 μs). | |
| int | novas_equals_vector (const double *a, const double *b, double tol) |
| Checks if two 3D vectors are effectively the same, within the specified absolute tolerance. | |
Set of functions to check for effective equality between SuperNOVAS data sturctures within typical tolerances.
Checks if two catalog entries define the same parameters, within the typical tolerances associated to these.
For two catalog entries to be equals, they must have matching
NOTES:
| a | one of the catalog entries |
| b | another catalog entry |
References novas_cat_entry::catalog, novas_cat_entry::dec, novas_cat_entry::parallax, novas_cat_entry::promodec, novas_cat_entry::promora, novas_cat_entry::ra, novas_cat_entry::radialvelocity, SIZE_OF_CAT_NAME, SIZE_OF_OBJ_NAME, novas_cat_entry::starname, and novas_cat_entry::starnumber.
Referenced by supernovas::CatalogEntry::equals(), and novas_equals_object().
| int novas_equals_frame | ( | const novas_frame * | a, |
| const novas_frame * | b ) |
Checks if two observing frames are essentially the same, within typical tolerances.
Two frames are considered equal, if they are both:
all within the typical tolerances associated to these. For non-geodetic observers two frames are considered equal also if both have all NAN (undefined) polar offsets – since polar offsets are not generally required for non-geodetic observing frames.
NOTES:
| a | an observing frame |
| b | another observing frame |
References novas_frame::accuracy, novas_frame::dx, novas_frame::dy, novas_equals_observer(), novas_equals_planet_bundle(), novas_equals_timespec(), novas_frame::observer, novas_frame::planets, and novas_frame::time.
Referenced by supernovas::Frame::equals().
Checks if two near-Earth locations and motions match within 1 mm and 1 mm/s, respectively.
Such near-Earth data structures are used by airborne observers or observers in Earth orbit.
NOTES:
| a | one near-Earth position/velocity data structure |
| b | another near-Earth position/velocity data structure |
References novas_equals_vector(), novas_in_space::sc_pos, and novas_in_space::sc_vel.
Referenced by novas_equals_observer().
Checks if two astronomical targets are the same within typical tolerances.
Two targets are considered equals only if their types, names (case sensitive), numerical IDs, and defining parameters match within typical tolerances.
NOTES:
| a | an astronomical target |
| b | another astronomical target |
References novas_object::name, NOVAS_CATALOG_OBJECT, NOVAS_EPHEM_OBJECT, novas_equals_cat_entry(), novas_equals_orbital(), NOVAS_ORBITAL_OBJECT, NOVAS_PLANET, novas_object::number, novas_object::orbit, SIZE_OF_OBJ_NAME, novas_object::star, and novas_object::type.
Referenced by supernovas::Source::equals().
Checks if two observers are essentially the same within the tolerances associated to their defining components.
NOTES:
| a | an observer data structure |
| b | another observer data structure |
References novas_observer::near_earth, NOVAS_AIRBORNE_OBSERVER, novas_equals_near_earth(), novas_equals_on_surface(), novas_equals_ssb_posvel(), NOVAS_OBSERVER_AT_GEOCENTER, NOVAS_OBSERVER_IN_EARTH_ORBIT, NOVAS_OBSERVER_ON_EARTH, NOVAS_SOLAR_SYSTEM_OBSERVER, novas_observer::on_surf, and novas_observer::where.
Referenced by supernovas::Observer::equals(), and novas_equals_frame().
| int novas_equals_on_surface | ( | const on_surface * | a, |
| const on_surface * | b ) |
Checks if two geodetic locations, and the weather parameters defined for them, match.
For two on_surface structures to be equal, they must have matching:
NOTES:
| a | one geodetic location (with weather) |
| b | another geodetic location (with weather) |
References novas_on_surface::height, novas_on_surface::humidity, novas_on_surface::latitude, novas_on_surface::longitude, novas_on_surface::pressure, and novas_on_surface::temperature.
Referenced by novas_equals_observer().
| int novas_equals_orbital | ( | const novas_orbital * | a, |
| const novas_orbital * | b ) |
Checks if two Keplerian orbitals match, within typical tolerances.
For two orbitals to be considered equal, they must have matching orbital systems, have the same reference dates to within ~10 ms (see novas_time_equals()), and:
NOTES:
| a | one set of Keplerian orbital parameters |
| b | another set of Keplerian orbital parameters |
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_AU, novas_equals_orbsys(), novas_orbital::Omega, novas_orbital::omega, and novas_orbital::system.
Referenced by supernovas::Orbital::equals(), and novas_equals_object().
| int novas_equals_orbsys | ( | const novas_orbital_system * | a, |
| const novas_orbital_system * | b ) |
Checks if two orbital systems match within typical tolerances.
Two orbital system are considered equal, if they are defined with respect to the same reference plane, around the same major planet (or Solar-system position), and have the obliquity and ascending node (if obliquity is non-zero) to within 1 μas.
NOTES:
| a | one of the orbital systems |
| b | another orbital system |
References novas_orbital_system::center, novas_orbital_system::obl, novas_orbital_system::Omega, novas_orbital_system::plane, and novas_orbital_system::type.
Referenced by supernovas::OrbitalSystem::equals(), and novas_equals_orbital().
| int novas_equals_planet_bundle | ( | const novas_planet_bundle * | a, |
| const novas_planet_bundle * | b ) |
Checks if two bundles containing Solar-system baricentric planet positions and velocities are effectively the same.
The two bundles are considered equals if:
NOTES:
| a | a planet position / velocity bundle |
| b | another planet position / velocity bundle |
References novas_planet_bundle::mask, NOVAS_AU, novas_equals_vector(), NOVAS_PLANETS, novas_planet_bundle::pos, and novas_planet_bundle::vel.
Referenced by novas_equals_frame().
Checks if two apparent positions on the observer's sky are the same within typical tolerances.
Two sky positions are considered to be equal if their:
NOTES:
| a | a position on the observer's sky |
| b | another sky position |
References novas_sky_pos::dec, novas_sky_pos::dis, novas_equals_vector(), novas_sky_pos::r_hat, novas_sky_pos::ra, and novas_sky_pos::rv.
Checks if two Solar-system locations and motions match within 1 m and ~1 mm/s, respectively.
Such near-Earth data structures are used by observers defined relative to the Solar-System Barycenter (SSB).
NOTES:
| a | one Solar-system position/velocity data structure |
| b | another Solar-system position/velocity data structure |
References NOVAS_AU, novas_equals_vector(), novas_in_space::sc_pos, and novas_in_space::sc_vel.
Referenced by novas_equals_observer().
| int novas_equals_timespec | ( | const novas_timespec * | a, |
| const novas_timespec * | b ) |
Checks if two time specifications are the same within 10-7 days (~100 μs).
NOTES:
| a | one of the time specs |
| b | the other time spec |
References novas_timespec::dut1, novas_timespec::fjd_tt, novas_timespec::ijd_tt, novas_timespec::tt2tdb, and novas_timespec::ut1_to_tt.
Referenced by novas_equals_frame().
| int novas_equals_vector | ( | const double * | a, |
| const double * | b, | ||
| double | tol ) |
Checks if two 3D vectors are effectively the same, within the specified absolute tolerance.
Two vectors are equal if the distance between them is less than or equal to the magnitude of the specified tolerance.
NOTES:
| a | one of the vectors |
| b | the other vector |
| tol | absolute tolerance (sign is ignored). |
Referenced by supernovas::Vector::equals(), novas_equals_near_earth(), novas_equals_planet_bundle(), novas_equals_sky_pos(), and novas_equals_ssb_posvel().