![]() |
SuperNOVAS C++ API v1.6
High-precision C/C++ astrometry library
|
Apparent places, defined on the local sky of an observer. More...
Classes | |
| class | supernovas::Apparent |
| Apparent position on sky as seen by an observer at a specific time of observation. More... | |
Macros | |
| #define | DEFAULT_GRAV_BODIES_FULL_ACCURACY ( DEFAULT_GRAV_BODIES_REDUCED_ACCURACY | (1 << NOVAS_JUPITER) | (1 << NOVAS_SATURN) ) |
| Bitwise mask defining a default set of gravitating bodies to use for deflection calculations in full accuracy mode. | |
| #define | DEFAULT_GRAV_BODIES_REDUCED_ACCURACY ( (1 << NOVAS_SUN) ) |
| Bitwise mask defining a sefault set of gravitating bodies to use for deflection calculations in reduced accuracy mode. | |
Functions | |
| Angle | supernovas::Source::angle_to (const Source &source, const Frame &frame) const |
| Returns the angular separation of this source from another source, for the given observer location and time of observation. | |
| Apparent | supernovas::Source::apparent_in (const Frame &frame) const |
| Returns the apparent position of a source (if possible), or else an invalid position. | |
| Apparent | supernovas::Frame::apparent_moon_elp2000 (double limit_term=0.0) const |
| Returns the Moon's apparent position using the ELP/MPP02 model by Chapront & Francou (2003) down to the specified limiting term amplitude. | |
| Apparent | supernovas::Planet::approx_apparent_in (const Frame &frame) const |
| Return an approximate apparent location on sky for a major planet, Sun, Moon, Earth-Moon Barycenter (EMB) – typically to arcmin level accuracy – using Keplerian orbital elements. | |
| Angle | supernovas::Source::moon_angle (const Frame &frame) const |
| Returns the angular separation of this source from the Moon, for the given observer location and time of observation. | |
| Angle | supernovas::Source::sun_angle (const Frame &frame) const |
| Returns the angular separation of this source from the Sun, for the given observer location and time of observation. | |
| Apparent | supernovas::Horizontal::to_apparent (const Frame &frame, const ScalarVelocity &rv=ScalarVelocity::stationary(), const Coordinate &distance=Coordinate::at_Gpc()) const |
| Converts these horizontal coordinates to an apparent place on the sky for an observer on or near Earth's surface. | |
| Apparent | supernovas::Horizontal::to_apparent (const Frame &frame, double rv=0.0, double distance=Unit::Gpc) const |
| Converts these horizontal coordinates to an apparent place on the sky. | |
Variables | |
| int | grav_bodies_full_accuracy |
| Current set of gravitating bodies to use for deflection calculations in full accuracy mode. | |
| int | grav_bodies_reduced_accuracy |
| Current set of gravitating bodies to use for deflection calculations in reduced accuracy mode. | |
Apparent places, defined on the local sky of an observer.
They are mainly a direction on the (e.g. R.A./Dec) on the sky and a spectroscopic radial velocity measure, from the observer's point of view. Unlike geometric locations, apparent positions are corrected for aberration for the observer's relative movement, and for gravitational deflection around the major gravitating solar-system bodies as light transits the Solar-system from the source to the observer.
| #define DEFAULT_GRAV_BODIES_FULL_ACCURACY ( DEFAULT_GRAV_BODIES_REDUCED_ACCURACY | (1 << NOVAS_JUPITER) | (1 << NOVAS_SATURN) ) |
Bitwise mask defining a default set of gravitating bodies to use for deflection calculations in full accuracy mode.
| #define DEFAULT_GRAV_BODIES_REDUCED_ACCURACY ( (1 << NOVAS_SUN) ) |
Bitwise mask defining a sefault set of gravitating bodies to use for deflection calculations in reduced accuracy mode.
(only apply gravitational deflection for the Sun.)
Returns the angular separation of this source from another source, for the given observer location and time of observation.
| Astronomical object of interest | the other source. |
| Observing frames | observing frame (observer location and time of observation) |
References Source(), supernovas::Frame::_novas_frame(), _object, supernovas::Unit::deg, and novas_object_sep().
Returns the apparent position of a source (if possible), or else an invalid position.
After the return, you should probably check for validity, e.g.:
There are multiple reasons why we might not be able to calculate valid apparent positions, such as:
The apparent position of a source is where it appears to the observer on the celestial sphere. As such it is mainly a direction on sky, which is corrected for light-travel time (i.e. where the source was at the time light originated from the Solar-system body, or the differential light-travel time between the Solar-system barycenter and the observer location for sidereal sources).
Unlike geometric positions, the apparent location is also corrected for the observer's motion (aberration), as well as gravitational deflection around the major Solar-system bodies.
| Observing frames | observer frame, which defines the observer location and the time of observation, as well as the accuracy requirement. |
References supernovas::Frame::_novas_frame(), _object, supernovas::Apparent::from_tod_sky_pos(), novas_sky_pos(), NOVAS_TOD, and supernovas::Apparent::undefined().
| Apparent supernovas::Frame::apparent_moon_elp2000 | ( | double | limit_term = 0.0 | ) | const |
Returns the Moon's apparent position using the ELP/MPP02 model by Chapront & Francou (2003) down to the specified limiting term amplitude.
NOTES:
REFERENCES:
| limit_term | (optional) [arcsec|km] Sum only the harmonic terms with amplitudes larger than this limit. |
References _novas_frame(), supernovas::Apparent::from_tod_sky_pos(), supernovas::Validating::is_valid(), novas_moon_elp_sky_pos_fp(), NOVAS_TOD, and supernovas::Apparent::undefined().
Return an approximate apparent location on sky for a major planet, Sun, Moon, Earth-Moon Barycenter (EMB) – typically to arcmin level accuracy – using Keplerian orbital elements.
The returned position is antedated for light-travel time (for Solar-System bodies). It also applies an appropriate aberration correction (but not gravitational deflection).
The orbitals can provide planet positions to arcmin-level precision for the rocky inner planets, and to a fraction of a degree precision for the gas and ice giants and Pluto. The accuracies for Uranus, Neptune, and Pluto are significantly improved (to the arcmin level) if used in the time range of 1800 AD to 2050 AD. For a more detailed summary of the typical accuracies, see either of the top two references below.
For accurate positions, you should use planetary ephemerides (such as the JPL ephemerides via the CALCEPH or CSPICE plugins) and novas_sky_pos() instead.
While this function is generally similar to creating an orbital object with an orbit initialized with novas_make_planet_orbit() or novas_make_moon_orbit(), and then calling novas_sky_pos(), there are a few important differences to note:
REFERENCES:
| Observing frames | The observing frame for which to calculate the apparent positions. |
References supernovas::Frame::_novas_frame(), supernovas::Apparent::from_tod_sky_pos(), supernovas::Validating::is_valid(), novas_approx_sky_pos(), novas_id(), and NOVAS_TOD.
Returns the angular separation of this source from the Moon, for the given observer location and time of observation.
| Observing frames | observing frame (observer location and time of observation) |
References supernovas::Frame::_novas_frame(), _object, supernovas::Unit::deg, and novas_moon_angle().
Returns the angular separation of this source from the Sun, for the given observer location and time of observation.
| Observing frames | observing frame (observer location and time of observation) |
References supernovas::Frame::_novas_frame(), _object, supernovas::Unit::deg, and novas_sun_angle().
| Apparent supernovas::Horizontal::to_apparent | ( | const Frame & | frame, |
| const ScalarVelocity & | rv = ScalarVelocity::stationary(), | ||
| const Coordinate & | distance = Coordinate::at_Gpc() ) const |
Converts these horizontal coordinates to an apparent place on the sky for an observer on or near Earth's surface.
Typically you should call this on unrefracted (astrometric) horizontal coordinates. If starting with observed (refracted) coordinates you should call to_unrefracted() first, before calling this function. The returned value may be invalid, e.g. if the frame is not that for an Earth-bound observer with a horizon. As such, it is best practice to check on the valifity of the returned value before using it, e.g. as:
| Observing frames | an Earth-based observing frame, defining the time of observation and the observer location, above (or slightly below) Earth's surface. |
| rv | (optional) observed radial velocity, if any (default: 0.0). |
| distance | (optional) apparent distance at which the observed light originated. |
References supernovas::Coordinate::m(), supernovas::ScalarVelocity::m_per_s(), and to_apparent().
| Apparent supernovas::Horizontal::to_apparent | ( | const Frame & | frame, |
| double | rv = 0.0, | ||
| double | distance = Unit::Gpc ) const |
Converts these horizontal coordinates to an apparent place on the sky.
Typically you should call this on unrefracted (astrometric) horizontal coordinates. If starting with observed (refracted) coordinates you should call to_unrefracted() first, before calling this function. The returned value may be invalid, e.g. if the frame is not that for an Earth-bound observer with a horizon. As such, it is best practice to check on the valifity of the returned value before using it, e.g. as:
| Observing frames | an Earth-based observing frame, defining the time of observation and the observer location, above (or slightly below) Earth's surface. |
| rv | [m/s] (optional) observed radial velocity, if any (default: 0.0). |
| distance | [m] (optional) apparent distance at which the observed light originated. |
References supernovas::Frame::_novas_frame(), supernovas::Unit::au, supernovas::Unit::day, novas_sky_pos::dec, novas_sky_pos::dis, supernovas::Apparent::from_tod_sky_pos(), supernovas::Validating::is_valid(), supernovas::Spherical::latitude(), supernovas::Spherical::longitude(), novas_hor_to_app(), NOVAS_TOD, novas_sky_pos::r_hat, novas_sky_pos::ra, radec2vector(), novas_sky_pos::rv, and supernovas::Apparent::undefined().
Referenced by to_apparent().
|
extern |
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.
Referenced by grav_def(), grav_undef(), novas_change_observer(), and place().
|
extern |
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.
Referenced by grav_def(), grav_undef(), novas_change_observer(), and place().