SuperNOVAS C++ API v1.6
High-precision C/C++ astrometry library
Loading...
Searching...
No Matches
Apparent equatorial positions on sky

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.

Detailed Description

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.

Macro Definition Documentation

◆ DEFAULT_GRAV_BODIES_FULL_ACCURACY

#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.

Since
1.1
Author
Attila Kovacs
See also
grav_bodies_full_accuracy, grav_planets(), grav_undo_planets()

◆ DEFAULT_GRAV_BODIES_REDUCED_ACCURACY

#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.)

Since
1.1
Author
Attila Kovacs
See also
grav_bodies_reduced_accuracy, grav_planets(), grav_undo_planets()

Function Documentation

◆ angle_to()

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.

Parameters
Astronomical object of interestthe other source.
Observing framesobserving frame (observer location and time of observation)
Returns
the distance between this source and the specified other source.
See also
sun_angle(), moon_angle()

References Source(), supernovas::Frame::_novas_frame(), _object, supernovas::Unit::deg, and novas_object_sep().

◆ apparent_in()

Apparent supernovas::Source::apparent_in ( const Frame & frame) const

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.:

Apparent app = my_source.apparent(frame);
if(!app.is_valid()) {
// We could not obtain apparent coordinates for some reason.
...
}
Apparent position on sky as seen by an observer at a specific time of observation.
Definition supernovas.h:2259
bool is_valid() const
Returns the previously set 'valid' stae of the implementing instance.
Definition supernovas.h:247

There are multiple reasons why we might not be able to calculate valid apparent positions, such as:

  • the frame itself may be invalid.
  • the system parameter may be outside of the enum range
  • for Solar system sources:
    • SuperNOVAS may not have a planet or ephemeris provider function configured for the given accuracy.
    • the planet or ephemeris provider function does not have data for the source, or planet at the orbital center (for orbital sources), for the requested time of observation.

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.

Parameters
Observing framesobserver frame, which defines the observer location and the time of observation, as well as the accuracy requirement.
Returns
the apparent position of the source, which may be invalid if the position cannot be determined.
See also
Planet::approx_apparent_in(), geometric_in()

References supernovas::Frame::_novas_frame(), _object, supernovas::Apparent::from_tod_sky_pos(), novas_sky_pos(), NOVAS_TOD, and supernovas::Apparent::undefined().

◆ apparent_moon_elp2000()

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:

  • The initial implementation (in v1.6) truncates the full series, keeping only terms with amplitudes larger than 1 mas (around 3400 harmonic terms in total), resulting in a limiting accuracy below 1 km level (and less than 100 meter error typically for 1900 – 2100).

REFERENCES:

Parameters
limit_term(optional) [arcsec|km] Sum only the harmonic terms with amplitudes larger than this limit.
Returns
The apparent place of the moon on the observer's sky.
See also
geometric_moon_elp2000()

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().

◆ approx_apparent_in()

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.

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:

  1. This function calculates Earth and Moon positions about the Keplerian orbital position of the Earth-Moon Barycenter (EMB). In contrast, novas_make_planet_orbit() does not provide orbitals for the Earth directly, and make_moot_orbit() references the Moon's orbital to the Earth position returned by the currently configured planet calculator function (see set_planet_provider()).
  2. This function ignores gravitational deflection. It makes little sense to bother about corrections that are orders of magnitude below the accuracy of the orbital positions obtained.

REFERENCES:

  1. E.M. Standish and J.G. Williams 1992.
  2. https://ssd.jpl.nasa.gov/planets/approx_pos.html
  3. Chapront, J. et al., 2002, A&A 387, 700–709
  4. Chapront-Touze, M., and Chapront, J. 1983, Astronomy and Astrophysics (ISSN 0004-6361), vol. 124, no. 1, July 1983, p. 50-62.
Parameters
Observing framesThe observing frame for which to calculate the apparent positions.
Returns
approximate, orbital model based, apparent position for the given planet. It may be invalid either if this planet or the frame argument is invalid.
See also
approx_geometric_in()

References supernovas::Frame::_novas_frame(), supernovas::Apparent::from_tod_sky_pos(), supernovas::Validating::is_valid(), novas_approx_sky_pos(), novas_id(), and NOVAS_TOD.

◆ moon_angle()

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.

Parameters
Observing framesobserving frame (observer location and time of observation)
Returns
the Moon's distance from the source.
See also
sun_angle(), angle_to()

References supernovas::Frame::_novas_frame(), _object, supernovas::Unit::deg, and novas_moon_angle().

◆ sun_angle()

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.

Parameters
Observing framesobserving frame (observer location and time of observation)
Returns
the Sun's distance from the source.
See also
moon_angle(), angle_to()

References supernovas::Frame::_novas_frame(), _object, supernovas::Unit::deg, and novas_sun_angle().

◆ to_apparent() [1/2]

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:

Horizontal hor = ...;
Apparent app = hor.to_apparent(frame);
if(!app) {
// Oops, could not calculate valid apparent positions.
return;
}
Apparent 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.
Definition Horizontal.cpp:248
Parameters
Observing framesan 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.
Returns
the apparent equatorial place corresponding to these astrometric horizontal coordinates on the sky, or Apparent::udefined() if the observing frame is not Earth surface based.
See also
to_unrefracted(), Apparent::to_horizontal()

References supernovas::Coordinate::m(), supernovas::ScalarVelocity::m_per_s(), and to_apparent().

◆ to_apparent() [2/2]

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:

Horizontal hor = ...;
Apparent app = hor.to_apparent(frame);
if(!app) {
// Oops, could not calculate valid apparent positions.
return;
}
Parameters
Observing framesan 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.
Returns
the apparent equatorial place corresponding to these astrometric horizontal coordinates on the sky, or Apparent::undefined() if the observing frame is not Earth surface based.
See also
to_unrefracted(), Apparent::to_horizontal()

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().

Variable Documentation

◆ grav_bodies_full_accuracy

int grav_bodies_full_accuracy
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.

Since
1.1
See also
grav_def(), grav_planets(), DEFAULT_GRAV_BODIES_FULL_ACCURACY
set_ephem_provider_hp()

Referenced by grav_def(), grav_undef(), novas_change_observer(), and place().

◆ grav_bodies_reduced_accuracy

int grav_bodies_reduced_accuracy
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.

Since
1.1
See also
grav_def(), grav_planets(), DEFAULT_GRAV_BODIES_REDUCED_ACCURACY
set_ephem_provider()

Referenced by grav_def(), grav_undef(), novas_change_observer(), and place().