PCL
pcl::EphemerisFile::Handle Class Reference

Calculation of ephemerides from data stored in XEPH files. More...

#include <EphemerisFile.h>

Public Member Functions

 Handle ()=default
 
 Handle (const EphemerisFile &parent, const IsoString &object, const IsoString &origin=IsoString())
 
 Handle (const Handle &x)
 
 Handle (Handle &&x)
 
virtual ~Handle ()
 
const Optional< double > & A1 () const
 
const Optional< double > & A2 () const
 
const Optional< double > & A3 () const
 
const Optional< double > & B_V () const
 
void ComputeFirstDerivative (Vector &v, TimePoint t)
 
void ComputeSecondDerivative (Vector &a, TimePoint t)
 
void ComputeState (Vector &p, TimePoint t)
 
void ComputeState (Vector &p, Vector &v, TimePoint t)
 
const Optional< double > & D () const
 
const Optional< double > & DT () const
 
TimePoint EndTime (int i=0) const
 
Vector FirstDerivative (TimePoint t)
 
const Optional< double > & G () const
 
const Optional< double > & H () const
 
bool HasDerivative () const
 
const Optional< double > & I_R () const
 
const Optional< double > & K1 () const
 
const Optional< double > & K2 () const
 
const Optional< double > & M1 () const
 
const Optional< double > & M2 () const
 
const IsoStringObjectId () const
 
const StringObjectName () const
 
Handleoperator= (const Handle &x)
 
Handleoperator= (Handle &&x)
 
const IsoStringOriginId () const
 
const EphemerisFileParentFile () const
 
const Optional< double > & PC () const
 
Vector SecondDerivative (TimePoint t)
 
TimePoint StartTime (int i=0) const
 
Vector StateVector (TimePoint t)
 
MultiVector StateVectors (TimePoint t)
 
const Optional< double > & U_B () const
 

Detailed Description

This subclass provides access to ephemeris data for a specific object available in an ephemeris file. It can perform basic ephemeris calculations, including state vectors and its first derivatives (such as position and velocity), and performs all the low-level file seek and read operations transparently and efficiently.

Data stored in an XEPH ephemeris file generates rectangular coordinates referred to the axes of the International Celestial Reference System (ICRS/J2000.0). Positions are given in au and velocities in au/day for all solar system objects, except planetocentric coordinates of natural satellites, including the geocentric Moon, for which positions and velocities are given in kilometers and km/day, respectively. Angles (nutations and librations) are given in radians, and time differences (such as TT-TDB) in seconds.

For performance reasons, this class is not thread-safe. However, as far as a given instance is never used concurrently from two or more threads, multiple instances of this class can be used from different running threads, including instances constructed to calculate ephemerides for the same object. This allows for the implementation of performance-intensive, lock-free multithreaded ephemeris calculation tasks.

Definition at line 1888 of file EphemerisFile.h.

Constructor & Destructor Documentation

◆ Handle() [1/4]

pcl::EphemerisFile::Handle::Handle ( )
default

Default constructor. Constructs an invalid instance that cannot be used unless assigned with a valid EphemerisFile::Handle object.

◆ Handle() [2/4]

pcl::EphemerisFile::Handle::Handle ( const EphemerisFile parent,
const IsoString object,
const IsoString origin = IsoString() 
)
inline

Constructs a new Handle object.

Parameters
parentReference to an open EphemerisFile object providing access to an ephemeris data file in XEPH format.
objectThe identifier or name of the object for which this instance will compute ephemerides.
originThe identifier of the origin of coordinates. If an empty string is specified (which is the default parameter value), this object will be created for the first occurrence of ephemeris data available for object in the parent file, irrespective of the origin. In such case, if parent is a standard fundamental ephemerides file, the origin is the solar system barycenter (identified as "SSB").

For a detailed description of object names and identifiers, see EphemerisFile::IsObjectAvailable().

If the specified parent ephemeris file is not open, or if no ephemeris data are available for the specified object and origin in the parent ephemeris file, this constructor will throw an Error exception.

Definition at line 1924 of file EphemerisFile.h.

◆ Handle() [3/4]

pcl::EphemerisFile::Handle::Handle ( const Handle x)
inline

Copy constructor.

Definition at line 1936 of file EphemerisFile.h.

◆ Handle() [4/4]

pcl::EphemerisFile::Handle::Handle ( Handle &&  x)
inline

Move constructor.

Definition at line 1952 of file EphemerisFile.h.

◆ ~Handle()

virtual pcl::EphemerisFile::Handle::~Handle ( )
inlinevirtual

Virtual destructor.

Definition at line 1965 of file EphemerisFile.h.

Member Function Documentation

◆ A1()

const Optional<double>& pcl::EphemerisFile::Handle::A1 ( ) const
inline

Returns the radial component of the comet non-gravitational acceleration, in au/day^2.

Definition at line 2459 of file EphemerisFile.h.

◆ A2()

const Optional<double>& pcl::EphemerisFile::Handle::A2 ( ) const
inline

Returns the transverse component of the comet non-gravitational acceleration, in au/day^2.

Definition at line 2468 of file EphemerisFile.h.

◆ A3()

const Optional<double>& pcl::EphemerisFile::Handle::A3 ( ) const
inline

Returns the normal component of the comet non-gravitational acceleration, in au/day^2.

Definition at line 2477 of file EphemerisFile.h.

◆ B_V()

const Optional<double>& pcl::EphemerisFile::Handle::B_V ( ) const
inline

Returns the color index B-V, in magnitudes.

Definition at line 2425 of file EphemerisFile.h.

◆ ComputeFirstDerivative()

void pcl::EphemerisFile::Handle::ComputeFirstDerivative ( Vector v,
TimePoint  t 
)
inline

Computes the first derivative of the state vector for the specified time point t.

Parameters
[out]vReference to a vector where the components of the computed first derivative will be stored.
tThe requested time point in the TDB time scale (rigorously, this is Teph as defined by JPL, but is equivalent to TDB or TT for all practical purposes).

Rectangular coordinates for velocity are provided in au/day, except for the geocentric Moon, for which velocity is provided in km/day.

Angle variations are provided in radians/day.

Time difference variations are provided in seconds/day.

The reference system is ICRS/J2000.0.

If the parent ephemeris file provides Chebyshev expansions of state vector derivatives for the object being calculated, the components of v will be calculated directly from these expansions. Otherwise the components of v will be approximated by numerical differentiation of the Chebyshev expansions for state vectors.

If t is either an invalid (uninitialized) TimePoint instance, or a time point outside the time span available from the parent ephemeris file, this member function throws an Error exception.

Definition at line 2104 of file EphemerisFile.h.

◆ ComputeSecondDerivative()

void pcl::EphemerisFile::Handle::ComputeSecondDerivative ( Vector a,
TimePoint  t 
)
inline

Computes the second derivative of the state vector for the specified time point t.

Parameters
[out]aReference to a vector where the components of the computed second derivative will be stored.
tThe requested time point in the TDB time scale (rigorously, this is Teph as defined by JPL, but is equivalent to TDB or TT for all practical purposes).

Rectangular coordinates for acceleration are provided in au/day^2, except for the geocentric Moon, for which acceleration is provided in km/day^2.

Angular second derivatives are provided in radians/day^2.

Time second derivatives are provided in seconds/day^2.

The reference system is ICRS/J2000.0.

The components of the second derivative will be approximated by numerical differentiation of the Chebyshev expansions for the first derivative. The latter can be either provided directly by the parent ephemeris file, or also approximated by numerical differentiation.

If t is either an invalid (uninitialized) TimePoint instance, or a time point outside the time span available from the parent ephemeris file, this member function throws an Error exception.

Definition at line 2148 of file EphemerisFile.h.

◆ ComputeState() [1/2]

void pcl::EphemerisFile::Handle::ComputeState ( Vector p,
TimePoint  t 
)
inline

Computes a state vector for the specified time point t.

Parameters
[out]pReference to a vector where the components of the computed state will be stored.
tThe requested time point in the TDB time scale (rigorously, this is Teph as defined by JPL, but is equivalent to TDB or TT for all practical purposes).

Rectangular coordinates for position are provided in au, except for the geocentric Moon, for which positions are provided in km.

Angles are provided in radians.

TT-TDB differences are provided in seconds.

The reference system is ICRS/J2000.0.

If t is either an invalid (uninitialized) TimePoint instance, or a time point outside the time span available from the parent ephemeris file, this member function throws an Error exception.

Definition at line 2027 of file EphemerisFile.h.

◆ ComputeState() [2/2]

void pcl::EphemerisFile::Handle::ComputeState ( Vector p,
Vector v,
TimePoint  t 
)
inline

Computes a state vector and its first derivative for the specified time point t.

Parameters
[out]pReference to a vector where the components of the computed state will be stored.
[out]vReference to a vector where the components of the computed first derivative will be stored.
tThe requested time point in the TDB time scale (rigorously, this is Teph as defined by JPL, but is equivalent to TDB or TT for all practical purposes).

Rectangular coordinates for position and velocity are provided in au and au/day, respectively, except for the geocentric Moon, for which position and velocity are provided in km and km/day, respectively.

Angles and their variations are provided in radians and radians/day.

Time differences (such as TT-TDB) and their variations are provided in seconds and seconds/day.

The reference system is ICRS/J2000.0.

If the parent ephemeris file provides Chebyshev expansions of state vector derivatives for the object being calculated, the components of v will be calculated directly from these expansions. Otherwise the components of v will be approximated by numerical differentiation of the Chebyshev expansions for state vectors.

If t is either an invalid (uninitialized) TimePoint instance, or a time point outside the time span available from the parent ephemeris file, this member function throws an Error exception.

Definition at line 2068 of file EphemerisFile.h.

◆ D()

const Optional<double>& pcl::EphemerisFile::Handle::D ( ) const
inline

Returns the object's diameter in km. When available, this is normally an IRAS diameter for an asteroid.

Definition at line 2450 of file EphemerisFile.h.

◆ DT()

const Optional<double>& pcl::EphemerisFile::Handle::DT ( ) const
inline

Returns the perihelion time offset parameter of the comet non-gravitational acceleration, in days.

Definition at line 2486 of file EphemerisFile.h.

◆ EndTime()

TimePoint pcl::EphemerisFile::Handle::EndTime ( int  i = 0) const
inline

Returns the upper bound of the time span for which this instance can calculate ephemerides using the ephemeris data already available. If a time point outside this span is requested, new file seek and read operations must be performed.

Parameters
iExpansion index: 0 to select function values (such as position), 1 to select function derivatives (such as velocity).

Definition at line 2315 of file EphemerisFile.h.

References pcl::Range().

◆ FirstDerivative()

Vector pcl::EphemerisFile::Handle::FirstDerivative ( TimePoint  t)
inline

Computes the first derivative of a state vector for the specified time point t.

Parameters
tThe requested time point in the TDB time scale (rigorously, this is Teph as defined by JPL, but is equivalent to TDB or TT for all practical purposes).

Returns the computed first derivative vector.

See ComputeFirstDerivative( Vector&, TimePoint ) for complete information.

Definition at line 2221 of file EphemerisFile.h.

◆ G()

const Optional<double>& pcl::EphemerisFile::Handle::G ( ) const
inline

Returns the asteroid magnitude slope parameter. See H() for references.

Definition at line 2357 of file EphemerisFile.h.

◆ H()

const Optional<double>& pcl::EphemerisFile::Handle::H ( ) const
inline

Returns the asteroid absolute magnitude. H is the visual magnitude of the object as seen at 1 au of the Earth, 1 au from the Sun, and with a phase angle of 0 degrees.

References

E. Bowell et al., Asteroids II, R. P. Binzel et al. (eds.), The University of Arizona Press, Tucson, 1989, pp. 549-554.

Urban, Sean E., Kenneth Seidelmann, P., ed. (2013), The Explanatory Supplement to the Astronomical Almanac 3rd Edition, Section 10.4.3.

Definition at line 2348 of file EphemerisFile.h.

◆ HasDerivative()

bool pcl::EphemerisFile::Handle::HasDerivative ( ) const
inline

Returns true iff the parent ephemeris file provides Chebyshev expansions of state vector derivatives for the object for which this handle computes ephemerides.

When no expansions for derivatives are available, derivatives are approximated by numerical differentiation of the Chebyshev expansions for state vectors.

Definition at line 2329 of file EphemerisFile.h.

◆ I_R()

const Optional<double>& pcl::EphemerisFile::Handle::I_R ( ) const
inline

Returns the color index I-R, in magnitudes.

Definition at line 2441 of file EphemerisFile.h.

◆ K1()

const Optional<double>& pcl::EphemerisFile::Handle::K1 ( ) const
inline

Returns the comet total magnitude slope parameter. See M1() for information on the calculation of comet apparent magnitudes.

Definition at line 2389 of file EphemerisFile.h.

◆ K2()

const Optional<double>& pcl::EphemerisFile::Handle::K2 ( ) const
inline

Returns the comet nuclear magnitude slope parameter. See M1() for information on the calculation of comet apparent magnitudes.

Definition at line 2408 of file EphemerisFile.h.

◆ M1()

const Optional<double>& pcl::EphemerisFile::Handle::M1 ( ) const
inline

Returns the comet total absolute magnitude. M1 is the visual absolute magnitude of the comet's combined nucleus and coma.

For the calculation of apparent comet magnitudes we apply the following equations:

Tmag = M1 + 5*log(d) + K1*log(r)
Nmag = M2 + 5*log(d) + K2*log(r) + PC*beta

where Tmag and Nmag are, respectively, the total (nucleus+coma) and nuclear apparent magnitudes. In these equations, M1 and M2 are the comet's total and nuclear absolute magnitude paranmeters, K1 and K2 are the total and nuclear magnitude slope parameters, PC is the nuclear magnitude phase coefficient, d is the comet's distance to Earth, r is its distance from the Sun, and beta is the phase angle.

Definition at line 2380 of file EphemerisFile.h.

◆ M2()

const Optional<double>& pcl::EphemerisFile::Handle::M2 ( ) const
inline

Returns the comet nuclear absolute magnitude. M2 is the visual absolute magnitude of the comet's nucleus. See M1() for information on the calculation of comet apparent magnitudes.

Definition at line 2399 of file EphemerisFile.h.

◆ ObjectId()

const IsoString& pcl::EphemerisFile::Handle::ObjectId ( ) const
inline

Returns the unique identifier of the object for which this handle can compute ephemerides.

Definition at line 2261 of file EphemerisFile.h.

◆ ObjectName()

const String& pcl::EphemerisFile::Handle::ObjectName ( ) const
inline

Returns the name of the object for which this handle can compute ephemerides. The returned value can be an empty string if no information about the object name is available in the parent ephemeris file.

Definition at line 2285 of file EphemerisFile.h.

◆ operator=() [1/2]

Handle& pcl::EphemerisFile::Handle::operator= ( const Handle x)
inline

Copy assignment operator. Returns a reference to this object.

Definition at line 1974 of file EphemerisFile.h.

◆ operator=() [2/2]

Handle& pcl::EphemerisFile::Handle::operator= ( Handle &&  x)
inline

Move assignment operator. Returns a reference to this object.

Definition at line 1991 of file EphemerisFile.h.

◆ OriginId()

const IsoString& pcl::EphemerisFile::Handle::OriginId ( ) const
inline

Returns the unique identifier of the origin of the coordinates computed by this object. For example, if this member function returns "Ea" or "Sn", the computed coordinates are geocentric or heliocentric, respectively. If this function returns "SSB", the origin of coordinates is the solar system barycenter. This happens, for example, with all standard XEPH files providing fundamental ephemerides.

Definition at line 2274 of file EphemerisFile.h.

◆ ParentFile()

const EphemerisFile& pcl::EphemerisFile::Handle::ParentFile ( ) const
inline

Returns a reference to the parent ephemeris file object. The parent object was specified upon construction.

Definition at line 2252 of file EphemerisFile.h.

◆ PC()

const Optional<double>& pcl::EphemerisFile::Handle::PC ( ) const
inline

Returns the comet nuclear magnitude phase coefficient. See M1() for information on the calculation of comet apparent magnitudes.

Definition at line 2417 of file EphemerisFile.h.

◆ SecondDerivative()

Vector pcl::EphemerisFile::Handle::SecondDerivative ( TimePoint  t)
inline

Computes the second derivative of a state vector for the specified time point t.

Parameters
tThe requested time point in the TDB time scale (rigorously, this is Teph as defined by JPL, but is equivalent to TDB or TT for all practical purposes).

Returns the computed second derivative vector.

See ComputeSecondDerivative( Vector&, TimePoint ) for complete information.

Definition at line 2241 of file EphemerisFile.h.

◆ StartTime()

TimePoint pcl::EphemerisFile::Handle::StartTime ( int  i = 0) const
inline

Returns the lower bound of the time span for which this instance can calculate ephemerides using the ephemeris data already available. If a time point outside this span is requested, new file seek and read operations must be performed.

Parameters
iExpansion index: 0 to select function values (such as position), 1 to select function derivatives (such as velocity).

Definition at line 2300 of file EphemerisFile.h.

References pcl::Range().

◆ StateVector()

Vector pcl::EphemerisFile::Handle::StateVector ( TimePoint  t)
inline

Computes a state vector for the specified time point t.

Parameters
tThe requested time point in the TDB time scale (rigorously, this is Teph as defined by JPL, but is equivalent to TDB or TT for all practical purposes).

Returns the computed state function value.

See ComputeState( Vector&, TimePoint ) for complete information.

Definition at line 2180 of file EphemerisFile.h.

◆ StateVectors()

MultiVector pcl::EphemerisFile::Handle::StateVectors ( TimePoint  t)
inline

Computes a state vector and its first derivative for the specified time point t.

Parameters
tThe requested time point in the TDB time scale (rigorously, this is Teph as defined by JPL, but is equivalent to TDB or TT for all practical purposes).

Returns a two-component multivector, where the first component vector is the state function's value and the second is its first derivative.

See ComputeState( Vector&, Vector&, TimePoint ) for complete information.

Definition at line 2201 of file EphemerisFile.h.

◆ U_B()

const Optional<double>& pcl::EphemerisFile::Handle::U_B ( ) const
inline

Returns the color index U-B, in magnitudes.

Definition at line 2433 of file EphemerisFile.h.


The documentation for this class was generated from the following file: