![]() |
SuperNOVAS C++ API v1.6
High-precision C/C++ astrometry library
|
An observing frame, defined by an observer location and precise time of observation. More...
#include <supernovas.h>
Public Member Functions | |
| Frame (const Frame &frame) | |
| Custom copy contructor, that points to a copy of the observer. | |
| Frame (const Observer &obs, const Time &time, enum novas_accuracy accuracy=NOVAS_FULL_ACCURACY) | |
| Instantiates a new observer frame, given the observer location and time of observation, and optionally the required accuracy. | |
| const novas_frame * | _novas_frame () const |
| Returns the pointer to the underlying NOVAS C novas_frame data structure of this observing frame. | |
| enum novas_accuracy | accuracy () const |
| Returns the accuracy type of this bserving frame. | |
| Apparent | 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. | |
| double | clock_skew (enum novas_timescale=NOVAS_TT) const |
| Returns the instantaneous incremental rate at which the observer's clock (i.e. | |
| const EOP | eop () const |
| Returns the Earth Orientation Parameters contained in this frame. | |
| Geometric | geometric (const Position &p, const Velocity &v, enum novas_reference_system system=NOVAS_TOD) const |
| Returns new geometric coordinates, relative to the observer in this frame, in the equatorial coordinate reference system of choice. | |
| Geometric | geometric_moon_elp2000 (double limit_term=0.0) const |
| Returns the Moon's geometric position using the ELP/MPP02 model by Chapront & Francou (2003), down to the specified limiting term amplitude. | |
| double | jd (enum novas_timescale=NOVAS_TT) const |
| Returns the precise Julian Date of this observing frame, in the specific timescale of choice. | |
| const Observer & | observer () const |
| Returns the observer location (and motion) of this observing frame. | |
| Position | observer_ssb_position () const |
| Returns the observer's ICRS position relative to the Solar System Barycenter (SSB). | |
| Velocity | observer_ssb_velocity () const |
| Returns the observer's ICRS velocity relative to the Solar System Barycenter (SSB). | |
| Frame & | operator= (const Frame &frame) |
| Custom copy-assignment operator, which updates the observer pointer to a copy of the original frame's observer. | |
| const Time & | time () const |
| Returns the astrometric time of observation of this observing frame. | |
| std::string | to_string () const |
| Returns a human-readable string representation of this observing frame. | |
| Public Member Functions inherited from supernovas::Validating | |
| bool | is_valid () const |
| Returns the previously set 'valid' stae of the implementing instance. | |
| operator bool () const | |
| Objects that implement Validating can be used in conditionals directly, without explicitly calling is_valid(). | |
Static Public Member Functions | |
| static Frame | reduced_accuracy (const Observer &obs, const Time &time) |
| Returns a reduced accuracy observing frame for the specified observer at the specified time. | |
| static const Frame & | undefined () |
| Returns a reference to a statically defined standard invalid observing frame. | |
Additional Inherited Members | |
| Protected Member Functions inherited from supernovas::Validating | |
| Validating () | |
| dummy constructor; | |
| Protected Attributes inherited from supernovas::Validating | |
| bool | _valid = false |
| the state variable. | |
An observing frame, defined by an observer location and precise time of observation.
Frames can be created with full (default) and reduced accuracy, supporting calculations with mas, or μas precisions typically. However note that full accuracy frames require SuperNOVAS to be configured with an appropriate high-precision planet ephemeris provider (see e.g. novas_use_calceph() or novas_use_cspice()), or else the resulting full-accuracy frame will be invalid.
Reduced accuracy frames may also be invalid if the low precision planet ephemeris provider ( which, by default, calculates approximate positions for the Earth and Sun only) cannot provide position for the Earth, Sun, the observer.
Therefore, one is strongly advised to check the validity of an observing frame after instantiation (using the is_valid() method), or else use the static Frame::create() function to return an optional.
| supernovas::Frame::Frame | ( | const Observer & | obs, |
| const Time & | time, | ||
| enum novas_accuracy | accuracy = NOVAS_FULL_ACCURACY ) |
Instantiates a new observer frame, given the observer location and time of observation, and optionally the required accuracy.
After the new frame is returned you should check that it's valid:
The returned new frame may be invalid for multiple reasons, such as:
In either case, you can obtain more information on why things went wrong, when they do, by enabling debug mode is enabled via novas_debug() prior to constructing a Frame.
| obs | observer location |
| time | time of observation |
| accuracy | (optional) NOVAS_FULL_ACCURACY (default) or NOVAS_REDUCED_ACCURACY. |
References supernovas::Observer::_novas_observer(), supernovas::Validating::_valid, accuracy(), supernovas::Observer::is_geodetic(), supernovas::Validating::is_valid(), novas_make_frame(), and time().
| supernovas::Frame::Frame | ( | const Frame & | frame | ) |
Custom copy contructor, that points to a copy of the observer.
| Observing frames | the frame to be copied. |
References supernovas::Validating::Validating(), and time().
| const novas_frame * supernovas::Frame::_novas_frame | ( | ) | const |
Returns the pointer to the underlying NOVAS C novas_frame data structure of this observing frame.
Referenced by supernovas::Source::angle_to(), supernovas::Source::apparent_in(), apparent_moon_elp2000(), supernovas::Planet::approx_apparent_in(), supernovas::Planet::approx_geometric_in(), supernovas::Source::equatorial_track(), supernovas::Source::geometric_in(), geometric_moon_elp2000(), supernovas::Source::horizontal_track(), supernovas::Source::moon_angle(), observer_ssb_position(), observer_ssb_velocity(), supernovas::Source::rises_above(), supernovas::Source::sets_below(), supernovas::SolarSystemSource::solar_illumination(), supernovas::Source::sun_angle(), supernovas::Horizontal::to_apparent(), and supernovas::Source::transits_in().
| enum novas_accuracy supernovas::Frame::accuracy | ( | ) | const |
Returns the accuracy type of this bserving frame.
Referenced by Frame().
| double supernovas::Frame::clock_skew | ( | enum novas_timescale | timescale = NOVAS_TT | ) | const |
Returns the instantaneous incremental rate at which the observer's clock (i.e.
proper time τ) ticks faster than a clock in the specified timescale in this observing frame. I.e., it returns D, which is defined by:
dτobs / dttimescale = (1 + D)
The instantaneous difference in clock rate includes tiny diurnal or orbital variationd for Earth-bound observers as the they cycle through the tidal potential around the geocenter (mainly due to the Sun and Moon). For a closer match to Earth-based timescales (TCG, TT, TAI, GPS, or UTC) you may want to exclude the periodic tidal effects and calculate the averaged observer clock rate over the geocentric cycle (see Eqs. 10.6 and 10.8 of the IERS Conventions 2010), which is provided by novas_mean_clock_skew() instead.
For reduced accuracy frames, the result will be approximate, because the gravitational effect of the Sun and Earth alone may be accounted for.
NOTES:
REFERENCES:
| timescale | Reference timescale for the comparison. All timescales except NOVAS_UT1 are supported. (UT1 advances at an irregular rate). |
References novas_clock_skew().
| const EOP supernovas::Frame::eop | ( | ) | const |
Returns the Earth Orientation Parameters contained in this frame.
The returned parameters are valid only if this frame was defined for a geodetic observer location such as an observer at a fixed site on Earth or an airborne observer location.
References eop(), supernovas::Unit::mas, and novas_time_leap().
Referenced by eop(), supernovas::Source::rises_above(), supernovas::Source::sets_below(), and supernovas::Source::transits_in().
| Geometric supernovas::Frame::geometric | ( | const Position & | p, |
| const Velocity & | v, | ||
| enum novas_reference_system | system = NOVAS_TOD ) const |
Returns new geometric coordinates, relative to the observer in this frame, in the equatorial coordinate reference system of choice.
| p | equatorial position vector, with respect to the observer |
| v | equatorial velocity vector, with respect to the observer |
| system | equatorial coordinate reference_system, in which position and velocity vectors are defined |
References supernovas::Validating::is_valid().
| double supernovas::Frame::jd | ( | enum novas_timescale | timescale = NOVAS_TT | ) | const |
Returns the precise Julian Date of this observing frame, in the specific timescale of choice.
It is but a shorthand for time().jd(timescale).
| timescale | (optional) the timescale in which to return the result (default: TT). |
Referenced by supernovas::Planet::approx_geometric_in(), and supernovas::Source::equatorial_track().
| const Observer & supernovas::Frame::observer | ( | ) | const |
Returns the observer location (and motion) of this observing frame.
Referenced by supernovas::Source::horizontal_track().
| Position supernovas::Frame::observer_ssb_position | ( | ) | const |
Returns the observer's ICRS position relative to the Solar System Barycenter (SSB).
References _novas_frame(), supernovas::Unit::AU, supernovas::Validating::is_valid(), and supernovas::Position::undefined().
| Velocity supernovas::Frame::observer_ssb_velocity | ( | ) | const |
Returns the observer's ICRS velocity relative to the Solar System Barycenter (SSB).
References _novas_frame(), supernovas::Unit::AU, supernovas::Unit::day, supernovas::Validating::is_valid(), and supernovas::Velocity::undefined().
| Frame & supernovas::Frame::operator= | ( | const Frame & | frame | ) |
Custom copy-assignment operator, which updates the observer pointer to a copy of the original frame's observer.
| Observing frames | the frame to be copied. |
References supernovas::Validating::_valid, and supernovas::Observer::copy().
Returns a reduced accuracy observing frame for the specified observer at the specified time.
Reduced accuracy frames provide 1 mas accuracy typically, and do not require a planet provider to be configured. As such, they offer a simplest way for obtaining astrometric positions for catalog and orbital sources at the 1 mas level.
Note, that the returned frame may be invalid, if the this observer or the time argument themselves are invalid.
| obs | Observer location |
| time | Astrometric time of observation |
References supernovas::Validating::is_valid(), NOVAS_REDUCED_ACCURACY, and time().
Referenced by supernovas::Observer::reduced_accuracy_frame_at().
| const Time & supernovas::Frame::time | ( | ) | const |
Returns the astrometric time of observation of this observing frame.
Referenced by Frame(), Frame(), and reduced_accuracy().
| std::string supernovas::Frame::to_string | ( | ) | const |
Returns a human-readable string representation of this observing frame.
References to_string().
Referenced by to_string().
|
static |
Returns a reference to a statically defined standard invalid observing frame.
This invalid frame may be used inside any object that is invalid itself.