PCL
pcl::EphemerisFile Class Reference

Solar system ephemerides from XEPH files. More...

#include <EphemerisFile.h>

Classes

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

Public Member Functions

 EphemerisFile ()=default
 
 EphemerisFile (const EphemerisFile &)=delete
 
 EphemerisFile (const String &filePath)
 
 EphemerisFile (EphemerisFile &&x)
 
virtual ~EphemerisFile () noexcept(false)
 
void Close ()
 
const EphemerisConstantListConstants () const
 
double ConstantValue (const IsoString &name) const
 
TimePoint EndTime () const
 
const StringFilePath () const
 
bool IsConstantAvailable (const IsoString &name) const
 
bool IsObjectAvailable (const IsoString &object, const IsoString &origin=IsoString()) const
 
bool IsOpen () const
 
const EphemerisMetadataMetadata () const
 
int NumberOfHandles () const
 
String ObjectName (const IsoString &object, const IsoString &origin=IsoString()) const
 
EphemerisObjectList Objects () const
 
void Open (const String &filePath)
 
EphemerisFileoperator= (const EphemerisFile &)=delete
 
EphemerisFileoperator= (EphemerisFile &&x)
 
TimePoint StartTime () const
 

Static Public Member Functions

static const EphemerisFileAsteroidEphemerides ()
 
static String CIP_ITRSDataFilePath ()
 
static String DeltaATDataFilePath ()
 
static String DeltaTDataFilePath ()
 
static const EphemerisFileFundamentalEphemerides ()
 
static const EphemerisFileKBOEphemerides ()
 
static const EphemerisFileNutationModel ()
 
static void OverrideAsteroidEphemerides (const String &filePath)
 
static void OverrideCIP_ITRSDataFilePath (const String &filePath)
 
static void OverrideDeltaATDataFilePath (const String &filePath)
 
static void OverrideDeltaTDataFilePath (const String &filePath)
 
static void OverrideFundamentalEphemerides (const String &filePath)
 
static void OverrideKBOEphemerides (const String &filePath)
 
static void OverrideNutationModel (const String &filePath)
 
static void OverrideShortTermAsteroidEphemerides (const String &filePath)
 
static void OverrideShortTermFundamentalEphemerides (const String &filePath)
 
static void OverrideShortTermKBOEphemerides (const String &filePath)
 
static void OverrideShortTermNutationModel (const String &filePath)
 
static void Serialize (const String &filePath, TimePoint startTime, TimePoint endTime, const SerializableEphemerisObjectDataList &data, const EphemerisMetadata &metadata=EphemerisMetadata(), const EphemerisConstantList &constants=EphemerisConstantList())
 
static const EphemerisFileShortTermAsteroidEphemerides ()
 
static const EphemerisFileShortTermFundamentalEphemerides ()
 
static const EphemerisFileShortTermKBOEphemerides ()
 
static const EphemerisFileShortTermNutationModel ()
 

Detailed Description

This class implements ephemerides of solar system bodies computed from data stored in XEPH (Extensible Ephemeris Data format) files. It also implements serialization of ephemeris data in the XEPH file format.

On the PixInsight/PCL platform, the XEPH file format allows for efficient ephemeris calculations through Chebyshev polynomial expansions stored as raw binary data. An XEPH file stores multiple arrays of Chebyshev coefficient series accessible by means of fast indexed file search algorithms and structures, along with auxiliary data and metadata required for ephemeris calculations.

Calculation of state vectors (such as position and velocity vectors) for specific objects is performed through a dedicated client subclass, namely EphemerisFile::Handle. This class implements transparent file seek and read operations, as well as fast, lock-free multithreaded evaluation of Chebyshev polynomials.

XEPH ephemeris files allow for calculation of 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.

An XEPH file storing up-to-date JPL DE/LE ephemeris data is part of all standard PixInsight distributions since 1.8.5 versions released Fall 2018. As of writing this documentation, the standard XEPH file provides the complete JPL DE440/LE440 ephemerides. See the EphemerisFile::FundamentalEphemerides() static member function for detailed information.

Definition at line 626 of file EphemerisFile.h.

Constructor & Destructor Documentation

◆ EphemerisFile() [1/4]

pcl::EphemerisFile::EphemerisFile ( )
default

Default constructor.

Constructs an invalid instance that cannot be used until initialized by calling the Open() member function.

◆ EphemerisFile() [2/4]

pcl::EphemerisFile::EphemerisFile ( const String filePath)
inline

Constructs an &EphemerisFile instance initialized from the specified ephemeris data file in XEPH format.

In the event of errors or invalid data, this constructor will throw the appropriate Error exception.

Definition at line 645 of file EphemerisFile.h.

◆ EphemerisFile() [3/4]

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

Move constructor.

Definition at line 653 of file EphemerisFile.h.

◆ EphemerisFile() [4/4]

pcl::EphemerisFile::EphemerisFile ( const EphemerisFile )
delete

Deleted copy constructor. EphemerisFile instances are unique, hence cannot be copied.

◆ ~EphemerisFile()

virtual pcl::EphemerisFile::~EphemerisFile ( )
virtualnoexcept

Virtual destructor.

Member Function Documentation

◆ AsteroidEphemerides()

static const EphemerisFile& pcl::EphemerisFile::AsteroidEphemerides ( )
static

Returns a reference to the global asteroid ephemerides file currently defined by the running PixInsight platform.

Under normal running conditions, the returned object provides ephemeris data for a set of asteroids with relevant masses. In a standard asteroid ephemeris file, object identifiers are asteroid numbers and object names are asteroid designations; for example:

IdentifierName
1Ceres
2Pallas
3Juno
4Vesta
5Astraea
......
702Alauda
703Noemi
704Interamnia
......

Asteroid ephemeris data are provided relative to the solar system barycenter ("SSB" identifier), with position and velocity coordinates coherent with global fundamental ephemerides.

As of writing this documentation, the standard asteroid ephemeris file provides the complete set of 343 asteroids used for the numerical integration of DE430 ephemerides, with barycentric coordinates coherent with DE440.

The asteroid ephemeris file can be overridden by the caller module. See the OverrideAsteroidEphemerides() member function for more information.

This static member function is thread-safe. It can be called safely from multiple execution threads running concurrently.

◆ CIP_ITRSDataFilePath()

static String pcl::EphemerisFile::CIP_ITRSDataFilePath ( )
static

Returns the path to the global database file of CIP positions referred to the ITRS.

The position of the Celestial Intermediate Pole (CIP) in the International Terrestrial Reference System (ITRS) is necessary to compute polar motion corrections applied to topocentric coordinates of solar system bodies. These corrections are relevant for the topocentric position of the Moon at the milliarcsecond level.

In current versions of PixInsight the CIP_ITRS database is a plain text file generated with values provided by the IERS Rapid Service/Prediction Center. As of writing this documentation, the main online reference is:

http://hpiers.obspm.fr/iers/eop/eopc01/eopc01.iau2000.1846-now

◆ Close()

void pcl::EphemerisFile::Close ( )

Closes the ephemeris file represented by this object and resets all internal structures to a default, uninitialized state.

If a previous file was already opened by this instance, it will be closed and all associated control and file indexing structures will be destroyed and deallocated. If no file is currently open, calling this member has no effect.

Warning
If this object has active ephemeris calculation handles, no action will be taken and an Error exception will be thrown. See EphemerisFile::Handle and NumberOfHandles().

◆ Constants()

const EphemerisConstantList& pcl::EphemerisFile::Constants ( ) const
inline

Returns a reference to the dynamic list of constants stored in the ephemeris data file represented by this object. These constants are name/value pairs, where keys are 8-bit strings and values are double precision floating point numbers.

Typically, for XEPH files that store fundamental JPL DE/LE ephemerides this member function returns the list of constants used by the corresponding DE/LE numerical integration.

The returned list is sorted by constant name in ascending order.

Definition at line 775 of file EphemerisFile.h.

◆ ConstantValue()

double pcl::EphemerisFile::ConstantValue ( const IsoString name) const
inline

Returns the value of an ephemeris constant given by its name.

If no constant is available with the specified name (case-insensitive) in the XEPH file represented by this object, this function throws an Error exception.

Definition at line 797 of file EphemerisFile.h.

References pcl::BinarySearch().

◆ DeltaATDataFilePath()

static String pcl::EphemerisFile::DeltaATDataFilePath ( )
static

Returns the path to the global database file of Delta AT values.

Delta AT is the difference TAI-UTC in seconds. In current versions of PixInsight the Delta AT database is a plain text file generated with values taken from the following online reference:

http://maia.usno.navy.mil/ser7/tai-utc.dat

Delta AT data files are expected to follow a simple format where each text line provides a UTC/DeltaAT pair. The exact format is described at the top of the corresponding file included in the current PixInsight distribution.

The Delta AT database file can be overridden by the caller module. See the OverrideDeltaATDataFilePath() member function for more information.

This static member function is thread-safe. It can be called safely from multiple execution threads running concurrently.

See also
TimePoint::DeltaAT()

◆ DeltaTDataFilePath()

static String pcl::EphemerisFile::DeltaTDataFilePath ( )
static

Returns the path to the global database file of observed Delta T values.

Delta T is the difference TT-UT1 in seconds. In current versions of PixInsight the Delta T database is a plain text file generated with values taken from the following online references:

For the period 1657-1973.0:
http://maia.usno.navy.mil/ser7/historic_deltat.data

For the period 1973 Feb 1 to the beginning of the current year:
http://maia.usno.navy.mil/ser7/deltat.data

For predicted Delta T values until 2025 approximately (as of writing this documentation, October 2018):
http://maia.usno.navy.mil/ser7/deltat.preds

Delta T data files are expected to follow a simple format where each text line provides a TT/DeltaT pair. The exact format is described at the top of the corresponding file included in the current PixInsight distribution.

Outside of the period from 1657 to the current year, the current PixInsight/PCL implementation uses polynomial expressions taken from Five Millennium Canon of Solar Eclipses, by Fred Espenak and Jean Meeus (NASA/TP–2006–214141, Revision 1.0, 2007).

The Delta T database file can be overridden by the caller module. See the OverrideDeltaTDataFilePath() member function for more information.

This static member function is thread-safe. It can be called safely from multiple execution threads running concurrently.

See also
TimePoint::DeltaT()

◆ EndTime()

TimePoint pcl::EphemerisFile::EndTime ( ) const
inline

Ending point of the time span covered by this ephemeris, or the latest time point for which ephemerides can be calculated using this object.

Definition at line 758 of file EphemerisFile.h.

◆ FilePath()

const String& pcl::EphemerisFile::FilePath ( ) const
inline

Returns the path of the ephemeris file represented by this object. Returned file paths are always absolute, full file paths.

Definition at line 739 of file EphemerisFile.h.

◆ FundamentalEphemerides()

static const EphemerisFile& pcl::EphemerisFile::FundamentalEphemerides ( )
static

Returns a reference to the global fundamental ephemerides file currently defined by the running PixInsight platform.

Under normal running conditions, the returned object provides ephemeris data for at least the following objects:

IdentifierName
MeMercury
VeVenus
EMBEarth-Moon barycenter
MaMars' barycenter
JuJupiter's barycenter
SaSaturn's barycenter
UrUranus' barycenter
NeNeptune's barycenter
PlPluto's barycenter
MnMoon's geometric center with respect to Earth's center.
SnSun's geometric center
EaEarth's geometric center

With the only exception of the Moon ("Mn" identifier), ephemeris data for all of the objects above are provided relative to the solar system barycenter ("SSB" identifier).

Additional items may also be available, depending on specific file versions and compilations:

IdentifierName
LbrLunar librations (Euler angles) in radians
NutNutation angles in radians
TT_TDBTT-TDB difference at the geocenter in seconds.

As of writing this documentation, the standard fundamental ephemeris file provides the complete JPL's DE440/LE440 ephemerides, but nutations, librations and time differences are not included.

The fundamental ephemeris file can be overridden by the caller module. See the OverrideFundamentalEphemerides() member function for more information.

This static member function is thread-safe. It can be called safely from multiple execution threads running concurrently.

◆ IsConstantAvailable()

bool pcl::EphemerisFile::IsConstantAvailable ( const IsoString name) const
inline

Returns true iff this instance knows a numerical constant with the specified name.

Definition at line 784 of file EphemerisFile.h.

References pcl::BinarySearch().

◆ IsObjectAvailable()

bool pcl::EphemerisFile::IsObjectAvailable ( const IsoString object,
const IsoString origin = IsoString() 
) const
inline

Returns true iff this instance contains ephemeris data for the specified object, given by its identifier or name, with respect to the specified origin.

The specified object string can be either an object identifier (case-sensitive), or an object name (case-insensitive) encoded as UTF-8. For example, all of 'Ju', 'Jupiter', 'jupiter' and 'JUPITER' refer to Jupiter in a standard XEPH file storing fundamental JPL DE/LE ephemerides. In this example, 'Ju' is a case-sensitive object identifier, thus 'ju' and 'JU' are not valid and this function would return false for both of them.

origin is the identifier of an origin of coordinates. If an empty string is specified (which is the default parameter value), this function will return true if the file contains any ephemeris data for the specified object, irrespective of the origin. Otherwise an exact match of the origin identifier will be required.

See the FundamentalEphemerides() static member function for information about the objects supported by standard XEPH files storing fundamental ephemerides.

Definition at line 846 of file EphemerisFile.h.

◆ IsOpen()

bool pcl::EphemerisFile::IsOpen ( ) const
inline

Returns true iff this object has an open ephemeris file and is ready for ephemeris data retrieval.

Definition at line 730 of file EphemerisFile.h.

References IsOpen().

Referenced by IsOpen().

◆ KBOEphemerides()

static const EphemerisFile& pcl::EphemerisFile::KBOEphemerides ( )
static

Returns a reference to the global Kuiper belt objects (KBOs) ephemerides file currently defined by the running PixInsight platform.

As of writing this documentation (March 2021), the default KBO ephemerides file includes the set of 30 most massive known trans-Neptunian objects used in JPL's DE440 numerical integration:

IdentifierName
19521Chaos
20000Varuna
28978Ixion
423012001 UR163
50000Quaoar
555652002 AW197
556372002 UX25
845222002 TC302
90377Sedna
90482Orcus
905682004 GV9
120347Salacia
136108Haumea
136199Eris
136472Makemake
1454522005 RN43
174567Varda
2089962003 AZ84
225088Gonggong
2309652004 XA192
2783612007 JJ43
3072612002 MS4
4555022003 UZ413
5236392010 RE64
5283812008 ST291
2004 XR190
2006 QH181
2010 FX86
2010 KZ39
2010 RF43

KBO ephemeris data are provided relative to the solar system barycenter ("SSB" identifier), with position and velocity coordinates coherent with global fundamental ephemerides. These ephemerides have been generated by numerical integration with starting state vectors provided by official NASA/JPL asteroid databases.

The KBO ephemeris file can be overridden by the caller module. See the OverrideKBOEphemerides() member function for more information.

This static member function is thread-safe. It can be called safely from multiple execution threads running concurrently.

◆ Metadata()

const EphemerisMetadata& pcl::EphemerisFile::Metadata ( ) const
inline

Returns a reference to the (immutable) metadata items available in the ephemeris data file loaded by this object.

Definition at line 884 of file EphemerisFile.h.

◆ NumberOfHandles()

int pcl::EphemerisFile::NumberOfHandles ( ) const
inline

Returns the number of handles currently active for this ephemeris file.

If an EphemerisFile object has active handles, also known as child handles, destroying it will most likely lead to a crash, since any activity performed by a child handle will make reference to a nonexistent parent object.

This function is thread-safe.

See also
pcl::EphemerisFile::Handle

Definition at line 901 of file EphemerisFile.h.

◆ NutationModel()

static const EphemerisFile& pcl::EphemerisFile::NutationModel ( )
static

Returns a reference to the global nutation model ephemeris file currently defined by the running PixInsight platform.

Under normal running conditions, the returned object provides Chebyshev polynomial expansions for the current IAU nutation model. As of writing this documentation, the standard nutation model file provides the IAU 2006/2000A_R nutation model (MHB2000 luni-solar and planetary nutation with adjustments to match the IAU 2006 precession).

The returned object should provide at least one object with the "IAUNut" identifier, which can be used to approximate the implemented nutation theory with a child EphemerisFile::Handle object.

The nutation model ephemeris file can be overridden by the caller module. See the OverrideNutationModel() member function for more information.

This static member function is thread-safe. It can be called safely from multiple execution threads running concurrently.

◆ ObjectName()

String pcl::EphemerisFile::ObjectName ( const IsoString object,
const IsoString origin = IsoString() 
) const
inline

Returns the name of the specified object, with respect to the specified origin.

Both object and origin must be object identifiers. If origin is an empty string (which is the default parameter value), this function will return the name of the first object found with the specified identifier, irrespective of the origin. Otherwise an exact match of the origin identifier will be required.

If no object with the required conditions is available in this ephemeris file, this function returns an empty string.

Definition at line 872 of file EphemerisFile.h.

◆ Objects()

EphemerisObjectList pcl::EphemerisFile::Objects ( ) const
inline

Returns a dynamic list of EphemerisObject instances describing all of the objects available in this file for ephemeris calculations.

The returned list is sorted by object and origin identifiers (in that order of precedence) in ascending order.

Definition at line 813 of file EphemerisFile.h.

◆ Open()

void pcl::EphemerisFile::Open ( const String filePath)

Initializes this object to give access to the specified ephemeris data file in XEPH format.

This member function opens an existing file at the specified filePath, loads and parses its XML header, and loads the file indexes ready for fast access to ephemeris data. The file will remain open until this object is destroyed, or until a new call to this function is made.

If a previous file was already opened by this instance, it will be closed and all associated control and file indexing structures will be destroyed and deallocated, before accessing the new file.

Warning
If this object has active ephemeris calculation handles, no action will be taken and an Error exception will be thrown. See EphemerisFile::Handle and NumberOfHandles().

◆ operator=() [1/2]

EphemerisFile& pcl::EphemerisFile::operator= ( const EphemerisFile )
delete

Deleted copy assignment operator. EphemerisFile instances are unique, hence cannot be copied.

◆ operator=() [2/2]

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

Move assignment operator. Returns a reference to this object.

Definition at line 665 of file EphemerisFile.h.

◆ OverrideAsteroidEphemerides()

static void pcl::EphemerisFile::OverrideAsteroidEphemerides ( const String filePath)
static

Override the default asteroid ephemerides file.

The specified filePath must be a valid path to an existing file in XEPH format, which must provide asteroid ephemerides coherent with the global fundamental ephemerides being used. See the AsteroidEphemerides() member function for a more comprehensive description.

After calling this member function, all asteroid ephemerides will be calculated using the specified XEPH file, which will be installed automatically upon the first call to AsteroidEphemerides().

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ OverrideCIP_ITRSDataFilePath()

static void pcl::EphemerisFile::OverrideCIP_ITRSDataFilePath ( const String filePath)
static

Override the global database of CIP positions in the ITRS.

The specified filePath must be a valid path to an existing file in a plain text database format compatible with the standard CIP_ITRS database included in PixInsight distributions. See the CIP_ITRSDataFilePath() member function for more details.

After calling this member function, CIP coordinates in the ITRS will be calculated by interpolation from the data provided by the specified file, which will be loaded and parsed automatically upon the first call (explicit or implicit) to Position::CIP_ITRS(). However, for performance and modularization reasons, once a CIP_ITRS database has been loaded there is no way to change it, so calling this function again will have no effect.

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ OverrideDeltaATDataFilePath()

static void pcl::EphemerisFile::OverrideDeltaATDataFilePath ( const String filePath)
static

Override the global Delta AT database.

The specified filePath must be a valid path to an existing file in a plain text database format compatible with the standard Delta AT database included in PixInsight distributions. See the DeltaATDataFilePath() member function for more details.

After calling this member function, Delta AT values will be calculated by interpolation from the data provided by the specified file, which will be loaded and parsed automatically upon the first call (explicit or implicit) to TimePoint::DeltaAT(). However, for performance and modularization reasons, once a Delta AT database has been loaded there is no way to change it, so calling this function again will have no effect.

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ OverrideDeltaTDataFilePath()

static void pcl::EphemerisFile::OverrideDeltaTDataFilePath ( const String filePath)
static

Override the global Delta T database.

The specified filePath must be a valid path to an existing file in a plain text database format compatible with the standard Delta T database included in PixInsight distributions. See the DeltaTDataFilePath() member function for more details.

After calling this member function, Delta T values will be calculated by interpolation from the data provided by the specified file, which will be loaded and parsed automatically upon the first call (explicit or implicit) to TimePoint::DeltaT(). However, for performance and modularization reasons, once a Delta T database has been loaded there is no way to change it, so calling this function again will have no effect.

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ OverrideFundamentalEphemerides()

static void pcl::EphemerisFile::OverrideFundamentalEphemerides ( const String filePath)
static

Override the default fundamental ephemerides file.

The specified filePath must be a valid path to an existing file in XEPH format, which must provide fundamental solar system ephemerides. See the FundamentalEphemerides() member function for a comprehensive description.

After calling this member function, all fundamental ephemerides will be calculated using the specified XEPH file, which will be installed automatically upon the first call to FundamentalEphemerides().

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ OverrideKBOEphemerides()

static void pcl::EphemerisFile::OverrideKBOEphemerides ( const String filePath)
static

Override the default Kuiper belt objects (KBOs) ephemerides file.

The specified filePath must be a valid path to an existing file in XEPH format, which must provide KBO ephemerides coherent with the global fundamental ephemerides being used. See the KBOEphemerides() member function for a more comprehensive description.

After calling this member function, all KBO ephemerides will be calculated using the specified XEPH file, which will be installed automatically upon the first call to KBOEphemerides().

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ OverrideNutationModel()

static void pcl::EphemerisFile::OverrideNutationModel ( const String filePath)
static

Override the global nutation model ephemeris file.

The specified filePath must be a valid path to an existing file in XEPH format, which must provide estimates of nutation angles compliant with the current IAU nutation model. See the NutationModel() member function for more details.

After calling this member function, the nutation model used for ephemeris calculations will be provided by the specified XEPH file, which will be installed automatically upon the first call to NutationModel().

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ OverrideShortTermAsteroidEphemerides()

static void pcl::EphemerisFile::OverrideShortTermAsteroidEphemerides ( const String filePath)
static

Override the default short-term asteroid ephemerides file.

The specified filePath must be a valid path to an existing file in XEPH format, which must provide asteroid ephemerides coherent with the global fundamental ephemerides being used. See the AsteroidEphemerides() member function for a more comprehensive description.

After calling this member function, all short-term asteroid ephemerides will be calculated using the specified XEPH file, which will be installed automatically upon the first call to ShortTermAsteroidEphemerides().

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ OverrideShortTermFundamentalEphemerides()

static void pcl::EphemerisFile::OverrideShortTermFundamentalEphemerides ( const String filePath)
static

Override the default short-term fundamental ephemerides file.

The specified filePath must be a valid path to an existing file in XEPH format, which must provide fundamental solar system ephemerides. See the FundamentalEphemerides() member function for a comprehensive description.

After calling this member function, all short-term fundamental ephemerides will be calculated using the specified XEPH file, which will be installed automatically upon the first call to ShortTermFundamentalEphemerides().

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ OverrideShortTermKBOEphemerides()

static void pcl::EphemerisFile::OverrideShortTermKBOEphemerides ( const String filePath)
static

Override the default short-term Kuiper belt objects (KBOs) ephemerides file.

The specified filePath must be a valid path to an existing file in XEPH format, which must provide KBO ephemerides coherent with the global fundamental ephemerides being used. See the KBOEphemerides() member function for a more comprehensive description.

After calling this member function, short-term KBO ephemerides will be calculated using the specified XEPH file, which will be installed automatically upon the first call to ShortTermKBOEphemerides().

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ OverrideShortTermNutationModel()

static void pcl::EphemerisFile::OverrideShortTermNutationModel ( const String filePath)
static

Override the global short-term nutation model ephemeris file.

The specified filePath must be a valid path to an existing file in XEPH format, which must provide estimates of nutation angles compliant with the current IAU nutation model. See the NutationModel() member function for more details.

After calling this member function, the nutation model used for short-term ephemeris calculations will be provided by the specified XEPH file, which will be installed automatically upon the first call to ShortTermNutationModel().

This function is useful to build standalone applications that don't depend on a running PixInsight core application, which is necessary to retrieve the default file names and directories from global settings.

◆ Serialize()

static void pcl::EphemerisFile::Serialize ( const String filePath,
TimePoint  startTime,
TimePoint  endTime,
const SerializableEphemerisObjectDataList data,
const EphemerisMetadata metadata = EphemerisMetadata(),
const EphemerisConstantList constants = EphemerisConstantList() 
)
static

Generates a file to store solar system ephemeris data in XEPH format.

Parameters
filePathPath to the file that will be generated in the local filesystem. The file name should carry the '.xeph' suffix.
startTimeLower bound of the entire time span for which the ephemeris data being serialized is valid.
endTimeUpper bound of the entire time span for which the ephemeris data being serialized is valid.
dataReference to an array of per-object ephemeris data, including Chebyshev polynomial expansions for the ephemeris function (such as position) and, optionally, for its first derivative (such as velocity), as well as object identifiers and object names.
metadataReference to an EphemerisMetadata structure with optional metadata information that will be included in the generated XEPH file.
constantsReference to an array of name/value pairs used to represent a set of numerical constants relevant to the ephemerides being serialized If an empty array is specified (as the default parameter value), no numerical constant will be included in the generated XEPH file.

The entire time span covered by an ephemeris file, from startTime to endTime, is usually subdivided into many small chunks or subspans, each of them with a relatively short polynomial expansion. The duration of each subspan is defined in a way such that the motion of the object for which positions are being calculated is sufficiently smooth to be fitted within the time subspan by truncated Chebyshev polynomials with relatively few coefficients (typically in the range of 15 to 30 coefficients) to achieve a prescribed accuracy. The faster and more perturbed the object's motion is, the more and shorter subspans are necessary to fit an accurate representation of the object's orbit.

In the event of invalid, incongruent or malformed data, or if an I/O error occurs, this function will throw an Error exception.

Warning
If a file already exists at the specified path, its previous contents will be lost after calling this function.

◆ ShortTermAsteroidEphemerides()

static const EphemerisFile& pcl::EphemerisFile::ShortTermAsteroidEphemerides ( )
static

Returns a reference to the global short-term asteroid ephemerides file currently defined by the running PixInsight platform.

See the AsteroidEphemerides() static member function for information on asteroid ephemerides and their status in current versions of PixInsight.

Under normal running conditions, the returned object should be a shortened version (that is, covering a shorter time span) of the standard asteroid ephemerides file. As of writing this documentation, the standard short-term asteroid ephemeris file covers the period from 1950 January 1.0 to 2100 January 32.0.

The short-term asteroid ephemeris file can be overridden by the caller module. See the OverrideShortTermAsteroidEphemerides() member function for more information.

This static member function is thread-safe. It can be called safely from multiple execution threads running concurrently.

◆ ShortTermFundamentalEphemerides()

static const EphemerisFile& pcl::EphemerisFile::ShortTermFundamentalEphemerides ( )
static

Returns a reference to the global short-term fundamental ephemerides file currently defined by the running PixInsight platform.

See the FundamentalEphemerides() static member function for information on fundamental ephemerides and their status in current versions of PixInsight.

Under normal running conditions, the returned object should be a shortened version (that is, covering a shorter time span) of the standard fundamental ephemerides file. As of writing this documentation, the standard short-term fundamental ephemeris file provides DE440/LE440 ephemerides for the period from 1850 January 1.0 to 2150 December 32.0.

The short-term fundamental ephemeris file can be overridden by the caller module. See the OverrideShortTermFundamentalEphemerides() member function for more information.

This static member function is thread-safe. It can be called safely from multiple execution threads running concurrently.

◆ ShortTermKBOEphemerides()

static const EphemerisFile& pcl::EphemerisFile::ShortTermKBOEphemerides ( )
static

Returns a reference to the global short-term Kuiper belt objects (KBOs) ephemerides file currently defined by the running PixInsight platform.

See the KBOEphemerides() static member function for information on asteroid ephemerides and their status in current versions of PixInsight.

Under normal running conditions, the returned object should be a shortened version (that is, covering a shorter time span) of the standard KBO ephemerides file. As of writing this documentation, the standard short-term asteroid ephemeris file covers the period from 1950 January 1.0 to 2100 January 32.0.

The short-term KBO ephemeris file can be overridden by the caller module. See the OverrideShortTermKBOEphemerides() member function for more information.

This static member function is thread-safe. It can be called safely from multiple execution threads running concurrently.

◆ ShortTermNutationModel()

static const EphemerisFile& pcl::EphemerisFile::ShortTermNutationModel ( )
static

Returns a reference to the global short-term nutation model ephemeris file currently defined by the running PixInsight platform.

See the NutationModel() static member function for information on nutation model ephemerides and their status in current versions of PixInsight.

Under normal running conditions, the returned object should be a shortened version (that is, covering a shorter time span) of the standard nutation model ephemerides file. As of writing this documentation, the standard short-term nutation model file provides the IAU 2006/2000A_R nutation model for the period from 1850 January 1.0 to 2150 December 32.0.

The returned object should provide at least one object with the "IAUNut" identifier, which can be used to approximate the implemented nutation theory with a child EphemerisFile::Handle object.

The short-term nutation model ephemeris file can be overridden by the caller module. See the OverrideShortTermNutationModel() member function for more information.

This static member function is thread-safe. It can be called safely from multiple execution threads running concurrently.

◆ StartTime()

TimePoint pcl::EphemerisFile::StartTime ( ) const
inline

Starting point of the time span covered by this ephemeris, or the earliest time point for which ephemerides can be calculated using this object.

Definition at line 749 of file EphemerisFile.h.


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