USGS

Isis 3.0 Developer's Reference (API)

Home

Isis::SpicePosition Class Reference
[Spice, Instruments, and Cameras]

Obtain SPICE position information for a body. More...

#include <SpicePosition.h>

List of all members.

Public Types

enum  Source { Spice, Memcache, HermiteCache, PolyFunction }
 

This enum defines indicates the status of the object.

More...
enum  PartialType { WRT_X, WRT_Y, WRT_Z }
enum  OverrideType { NoOverrides, ScaleOnly, BaseAndScale }

Public Member Functions

 SpicePosition (int targetCode, int observerCode)
 Construct an empty SpicePosition class using valid body codes.
virtual ~SpicePosition ()
 Destructor.
void SetTimeBias (double timeBias)
 Apply a time bias when invoking SetEphemerisTime method.
void SetAberrationCorrection (const std::string &correction)
 Set the aberration correction (light time).
const std::vector< double > & SetEphemerisTime (double et)
 Return J2000 coordinate at given time.
double EphemerisTime () const
 Return the current ephemeris time.
const std::vector< double > & GetCenterCoordinate ()
 Compute and return the coordinate at the center time.
const std::vector< double > & Coordinate ()
 Return the current J2000 position.
const std::vector< double > & Velocity ()
 Return the current J2000 velocity.
bool HasVelocity ()
 Return the flag indicating whether the velocity exists.
void LoadCache (double startTime, double endTime, int size)
 Cache J2000 position over a time range.
void LoadCache (double time)
 Cache J2000 position for a time.
void LoadCache (Table &table)
 Cache J2000 positions using a table file.
Table LineCache (const std::string &tableName)
 Return a table with J2000 to reference positions.
Table LoadHermiteCache (const std::string &tableName)
 Cache J2000 position over existing cached time range using polynomials stored as Hermite cubic spline knots.
void ReloadCache ()
 Cache J2000 positions over existing cached time range using polynomials.
void ReloadCache (Table &table)
 Cache J2000 position over existing cached time range using table.
Table Cache (const std::string &tableName)
 Return a table with J2000 positions.
bool IsCached () const
 Is this position cached.
void SetPolynomial ()
 Set the coefficients of a polynomial fit to each of the components (X, Y, Z) of the position vector for the time period covered by the cache, component = c0 + c1*t + c2*t**2 + .
void SetPolynomial (const std::vector< double > &XC, const std::vector< double > &YC, const std::vector< double > &ZC)
 Set the coefficients of a polynomial (parabola) fit to each of the three coordinates of the position vector for the time period covered by the cache, coord = c0 + c1*t + c2*t**2 + .
void GetPolynomial (std::vector< double > &XC, std::vector< double > &YC, std::vector< double > &ZC)
 Return the coefficients of a polynomial fit to each of the three coordinates of the position for the time period covered by the cache, angle = c0 + c1*t + c2*t**2 + .
void SetPolynomialDegree (int degree)
 Set the polynomial degree.
Source GetSource ()
 Return the source of the position.
void ComputeBaseTime ()
 Compute the base time using cached times.
double GetBaseTime ()
 Return the base time for the position.
void SetOverrideBaseTime (double baseTime, double timeScale)
 Set an override base time to be used with observations on scanners to allow all images in an observation to use the same base time and polynomials for the positions.
double GetTimeScale ()
 Return the time scale for the position.
double DPolynomial (const int coeffIndex)
 Evaluate the derivative of the fit polynomial (parabola) defined by the given coefficients with respect to the coefficient at the given index, at the current time.
std::vector< double > CoordinatePartial (SpicePosition::PartialType partialVar, int coeffIndex)
 Set the coefficients of a polynomial fit to each of the three coordinates of the position vector for the time period covered by the cache,.
std::vector< double > VelocityPartial (SpicePosition::PartialType partialVar, int coeffIndex)
 Compute the derivative of the velocity with respect to the specified variable.
void Memcache2HermiteCache (double tolerance)
 This method reduces the cache for position, time and velocity to the minimum number of values needed to interpolate the J2000 coordinates using a Hermite spline, given a tolerance of deviation from the NAIF values.
std::vector< double > Extrapolate (double timeEt)
 Extrapolate position for a given time assuming a constant velocity.

Protected Member Functions

void SetEphemerisTimeMemcache ()
 This is a protected method that is called by SetEphemerisTime() when Source type is Memcache.
void SetEphemerisTimeHermiteCache ()
 This is a protected method that is called by SetEphemerisTime() when Source type is HermiteCache.
void SetEphemerisTimeSpice ()
 This is a protected method that is called by SetEphemerisTime() when Source type is Spice.
void SetEphemerisTimePolyFunction ()
 This is a protected method that is called by SetEphemerisTime() when Source type is PolyFunction.
std::vector< int > HermiteIndices (double tol, std::vector< int > indexList)
 This method is called by Memcache2HermiteCache() to determine which indices from the orginal cache should be saved in the reduced cache.

Detailed Description

Obtain SPICE position information for a body.

This class will obtain the J2000 body position between a target and observer body. For example, a spacecraft and Mars or the Sun and Mars. It is essentially a C++ wrapper to the NAIF spkez_c routine. Therefore, appropriate NAIF kernels are expected to be loaded prior to using this class. The position can be returned with or without one way light time corrections between the two bodies. See NAIF required reading for more information regarding this subject at ftp://naif.jpl.nasa.gov/pub/naif/toolkit_docs/C/ascii/individual_docs/spk.req

An important functionality of this class is the ability to cache the positions so they do not have to be constantly read from the NAIF kernels. Onced the data is cached the NAIF kernels can be unloaded too.

Author:
2005-12-01 Jeff Anderson

Member Enumeration Documentation

Enumerator:
NoOverrides 
ScaleOnly 
BaseAndScale 
Enumerator:
WRT_X 
WRT_Y 
WRT_Z 

This enum defines indicates the status of the object.

Enumerator:
Spice 

Object is reading directly from the kernels.

Memcache 

Object is reading from cached table.

HermiteCache 

Object is reading from splined table.

PolyFunction 

Object is calculated from nth degree polynomial.


Constructor & Destructor Documentation

Isis::SpicePosition::SpicePosition ( int  targetCode,
int  observerCode 
)

Construct an empty SpicePosition class using valid body codes.

See required reading ftp://naif.jpl.nasa.gov/pub/naif/toolkit_docs/C/ascii/individual_docs/naif_ids.req

Parameters:
targetCode Valid naif body name.
observerCode Valid naif body name.

References NoOverrides.

Isis::SpicePosition::~SpicePosition (  )  [virtual]

Destructor.

Free the memory allocated by this SpicePosition instance.


Member Function Documentation

Table Isis::SpicePosition::Cache ( const std::string &  tableName  ) 

Return a table with J2000 positions.

Return a table containg the cached coordinates with the given name. The table will have four or seven columns, J2000 x,y,z (optionally vx,vy,vx) and the ephemeris time.

Parameters:
tableName Name of the table to create and return
void Isis::SpicePosition::ComputeBaseTime (  ) 

Compute the base time using cached times.

Referenced by Memcache2HermiteCache(), SetEphemerisTimeHermiteCache(), and SetPolynomial().

const std::vector<double>& Isis::SpicePosition::Coordinate (  )  [inline]
std::vector< double > Isis::SpicePosition::CoordinatePartial ( SpicePosition::PartialType  partialVar,
int  coeffIndex 
)

Set the coefficients of a polynomial fit to each of the three coordinates of the position vector for the time period covered by the cache,.

coordinate = c0 + c1*t + c2*t**2 + ... cn*t**n, where t = (time - p_basetime) / p_timeScale.

Parameters:
partialVar Designated variable of the partial derivative
Returns:
Derivative of j2000 vector calculated with polynomial with respect to partialVar

References DPolynomial().

Referenced by Isis::RadarGroundMap::GetdXYdPosition(), and Isis::CameraGroundMap::GetdXYdPosition().

double Isis::SpicePosition::DPolynomial ( const int  coeffIndex  ) 

Evaluate the derivative of the fit polynomial (parabola) defined by the given coefficients with respect to the coefficient at the given index, at the current time.

Parameters:
coeffIndex Index of coefficient to differentiate with respect to
Returns:
The derivative evaluated at the current time

Referenced by CoordinatePartial().

double Isis::SpicePosition::EphemerisTime (  )  const [inline]

Return the current ephemeris time.

std::vector< double > Isis::SpicePosition::Extrapolate ( double  timeEt  ) 

Extrapolate position for a given time assuming a constant velocity.

The position and velocity at the current time will be used to extrapolate position at the input time. If velocity does not exist, the value at the current time will be output. The caller must make sure to call SetEphemerisTime to set the time to the time to be used for the extrapolation.

Parameters:
[in] timeEt The time of the position to be extrapolated
[out] An extrapolated position at the input time
double Isis::SpicePosition::GetBaseTime (  )  [inline]

Return the base time for the position.

Referenced by Isis::BundleAdjust::Solve(), and Isis::BundleAdjust::SolveCholesky().

const std::vector< double > & Isis::SpicePosition::GetCenterCoordinate (  ) 

Compute and return the coordinate at the center time.

void Isis::SpicePosition::GetPolynomial ( std::vector< double > &  XC,
std::vector< double > &  YC,
std::vector< double > &  ZC 
)

Return the coefficients of a polynomial fit to each of the three coordinates of the position for the time period covered by the cache, angle = c0 + c1*t + c2*t**2 + .

.. + cn*t**n, where t = (time - p_basetime) / p_timeScale.

Parameters:
[out] XC Coefficients of fit to first coordinate of position
[out] YC Coefficients of fit to second coordinate of position
[out] ZC Coefficients of fit to third coordinate of position

Referenced by Isis::BundleAdjust::Solve(), and Isis::BundleAdjust::SolveCholesky().

Source Isis::SpicePosition::GetSource (  )  [inline]

Return the source of the position.

Referenced by Isis::Spice::CreateCache().

double Isis::SpicePosition::GetTimeScale (  )  [inline]

Return the time scale for the position.

Referenced by Isis::BundleAdjust::Solve(), and Isis::BundleAdjust::SolveCholesky().

bool Isis::SpicePosition::HasVelocity (  )  [inline]

Return the flag indicating whether the velocity exists.

std::vector< int > Isis::SpicePosition::HermiteIndices ( double  tolerance,
std::vector< int >  indexList 
) [protected]

This method is called by Memcache2HermiteCache() to determine which indices from the orginal cache should be saved in the reduced cache.

It is a recursive method that starts with an index list of 3 elements (first, center and last index values) and adds values to this list if the tolerance is not met.

Parameters:
tolerance Maximum error allowed between NAIF kernel coordinate values and values interpolated by the Hermite spline.
indexList Vector containing the list of indices to be kept in the cache. This list grows as the method is recursively called
Returns:
vector/<double> Vector containing final list of indices that will be kept for the Hermite cache.

Referenced by Memcache2HermiteCache().

bool Isis::SpicePosition::IsCached (  )  const [inline]

Is this position cached.

Referenced by Isis::Spice::CreateCache().

Table Isis::SpicePosition::LineCache ( const std::string &  tableName  ) 

Return a table with J2000 to reference positions.

Return a table containing the cached positions with the given name. The table will have seven columns, positionX, positionY, positionZ, angular velocity X, angular velocity Y, angular velocity Z, and time of J2000 position.

Parameters:
tableName Name of the table to create and return
void Isis::SpicePosition::LoadCache ( Table &  table  ) 

Cache J2000 positions using a table file.

This method will load an internal cache with coordinates from an ISIS table file. The table must have 4 columns, or 7 (if velocity) is included, and at least one row. The 4 columns contain the following information, body position x,y,z in J2000 and the ephemeris time of that position. If there are multiple rows it is assumed you can interpolate position at times in between the rows.

Parameters:
table An ISIS table blob containing valid J2000 coordinate/time values
void Isis::SpicePosition::LoadCache ( double  time  ) 

Cache J2000 position for a time.

This method will load an internal cache with coordinates for a single time (e.g. useful for framing cameras). This prevents the NAIF kernels from being read over-and-over again and slowing a application down due to I/O performance. Once the cache has been loaded then the kernels can be unloaded from the NAIF system. This calls the LoadCache(stime,etime,size) method using the time as both the starting and ending time with a size of 1.

Parameters:
time single ephemeris time in seconds to cache
void Isis::SpicePosition::LoadCache ( double  startTime,
double  endTime,
int  size 
)

Cache J2000 position over a time range.

This method will load an internal cache with coordinates over a time range. This prevents the NAIF kernels from being read over-and-over again and slowing a application down due to I/O performance. Once the cache has been loaded then the kernels can be unloaded from the NAIF system.

Parameters:
startTime Starting ephemeris time in seconds for the cache
endTime Ending ephemeris time in seconds for the cache
size Number of coordinates/positions to keep in the cache

Referenced by Isis::Spice::CreateCache(), and ReloadCache().

Table Isis::SpicePosition::LoadHermiteCache ( const std::string &  tableName  ) 

Cache J2000 position over existing cached time range using polynomials stored as Hermite cubic spline knots.

This method will reload an internal cache with positions formed from a cubic Hermite spline over a time range. The method assumes a polynomial function has been fit to the coordinates of the positions and calculates the spline from the polynomial function.

void Isis::SpicePosition::Memcache2HermiteCache ( double  tolerance  ) 

This method reduces the cache for position, time and velocity to the minimum number of values needed to interpolate the J2000 coordinates using a Hermite spline, given a tolerance of deviation from the NAIF values.

Parameters:
tolerance Maximum error allowed between NAIF kernel coordinate values and values interpolated by the Hermite spline.

References _FILEINFO_, ComputeBaseTime(), HermiteCache, HermiteIndices(), Memcache, Isis::iException::Message(), Isis::iException::Programmer, and ScaleOnly.

Referenced by Isis::Spice::CreateCache().

void Isis::SpicePosition::ReloadCache ( Table &  table  ) 

Cache J2000 position over existing cached time range using table.

This method will reload an internal cache with positions formed from coordinates in a table

Parameters:
table An ISIS table blob containing valid J2000 coordinate/time values.

References LoadCache().

void Isis::SpicePosition::ReloadCache (  ) 

Cache J2000 positions over existing cached time range using polynomials.

This method will reload an internal cache with positions calculated from functions fit to the coordinates of the position over a time range.

void Isis::SpicePosition::SetAberrationCorrection ( const std::string &  correction  ) 

Set the aberration correction (light time).

See NAIF required reading for more information on this correction at ftp://naif.jpl.nasa.gov/pub/naif/toolkit_docs/C/ascii/individual_docs/spk.req

Parameters:
correction This value must be one of the following: "NONE", "LT", "LT+S", where LT is a correction for planetary aberration (light time) and S a correction for stellar aberration. If never called the default is "LT+S".

Referenced by Isis::Spice::CreateCache(), Isis::LoHighCamera::LoHighCamera(), Isis::LoMediumCamera::LoMediumCamera(), Isis::Mariner10Camera::Mariner10Camera(), and Isis::Spice::SetTime().

const std::vector< double > & Isis::SpicePosition::SetEphemerisTime ( double  et  ) 

Return J2000 coordinate at given time.

This method returns the J2000 coordinates (x,y,z) of the body at a given et in seconds. The coordinates are obtained from either a valid NAIF spk kernel, or alternatively from an internal cache loaded from an ISIS Table object. In the first case, the SPK kernel must contain positions for the body code specified in the constructor at the given time and it must be loaded using the SpiceKernel class.

Parameters:
et ephemeris time in seconds

Referenced by SetPolynomial(), and Isis::Spice::SetTime().

void Isis::SpicePosition::SetEphemerisTimeHermiteCache (  )  [protected]

This is a protected method that is called by SetEphemerisTime() when Source type is HermiteCache.

It calculates J2000 coordinates (x,y,z) of the body that correspond to a given et in seconds. These coordinates are obtained by using a Hermite spline to interpolate values from an internal reduced cache loaded from an ISIS Table object.

See also:
SetEphemerisTime()

References _FILEINFO_, Isis::NumericalApproximation::AddCubicHermiteDeriv(), Isis::NumericalApproximation::AddData(), ComputeBaseTime(), Isis::NumericalApproximation::CubicHermite, Isis::NumericalApproximation::Evaluate(), Isis::NumericalApproximation::EvaluateCubicHermiteFirstDeriv(), Isis::NumericalApproximation::Extrapolate, Isis::iException::Io, Isis::iException::Message(), and ScaleOnly.

void Isis::SpicePosition::SetEphemerisTimeMemcache (  )  [protected]

This is a protected method that is called by SetEphemerisTime() when Source type is Memcache.

It calculates J2000 coordinates (x,y,z) of the body that correspond to a given et in seconds. These coordinates are obtained from an internal cache loaded from an ISIS Table object.

See also:
SetEphemerisTime()
void Isis::SpicePosition::SetEphemerisTimePolyFunction (  )  [protected]

This is a protected method that is called by SetEphemerisTime() when Source type is PolyFunction.

It calculates J2000 coordinates (x,y,z) of the body that correspond to a given et in seconds. These coordinates are obtained by using an nth degree polynomial function fit to each coordinate of the position vector.

See also:
SetEphemerisTime()

References Isis::BasisFunction::Evaluate(), Isis::BasisFunction::SetCoefficients(), WRT_X, WRT_Y, and WRT_Z.

void Isis::SpicePosition::SetEphemerisTimeSpice (  )  [protected]

This is a protected method that is called by SetEphemerisTime() when Source type is Spice.

It calculates J2000 coordinates (x,y,z) of the body that correspond to a given et in seconds. The coordinates are obtained from a valid NAIF spk kernel. The SPK kernel must contain positions for the body code specified in the constructor at the given time and it must be loaded using the SpiceKernel class.

See also:
SetEphemerisTime()
void Isis::SpicePosition::SetOverrideBaseTime ( double  baseTime,
double  timeScale 
)

Set an override base time to be used with observations on scanners to allow all images in an observation to use the same base time and polynomials for the positions.

Parameters:
[in] baseTime The baseTime to use and override the computed base time
History:
2011-04-08 Debbie A. Cook - Corrected p_override set to BaseAndScale

Referenced by Isis::BundleAdjust::Solve(), and Isis::BundleAdjust::SolveCholesky().

void Isis::SpicePosition::SetPolynomial ( const std::vector< double > &  XC,
const std::vector< double > &  YC,
const std::vector< double > &  ZC 
)

Set the coefficients of a polynomial (parabola) fit to each of the three coordinates of the position vector for the time period covered by the cache, coord = c0 + c1*t + c2*t**2 + .

.. + cn*t**n, where t = (time - p_baseTime) / p_timeScale.

Parameters:
[in] XC Coefficients of fit to X coordinate
[in] YC Coefficients of fit to Y coordinate
[in] ZC Coefficients of fit to Z coordinate

References ComputeBaseTime(), PolyFunction, Isis::BasisFunction::SetCoefficients(), and SetEphemerisTime().

void Isis::SpicePosition::SetPolynomial (  ) 

Set the coefficients of a polynomial fit to each of the components (X, Y, Z) of the position vector for the time period covered by the cache, component = c0 + c1*t + c2*t**2 + .

.. + cn*t**n, where t = (time - p_baseTime) / p_timeScale.

Referenced by Isis::BundleAdjust::Solve(), and Isis::BundleAdjust::SolveCholesky().

void Isis::SpicePosition::SetPolynomialDegree ( int  degree  ) 

Set the polynomial degree.

Set the degree of the polynomials to be fit to the three position coordinates for the time period covered by the cache, coordinate = c0 + c1*t + c2*t**2 + .

.. + cn*t**n, where t = (time - p_baseTime) / p_timeScale, and n = p_degree.

Parameters:
[in] degree Degree of the polynomial to be fit
void Isis::SpicePosition::SetTimeBias ( double  timeBias  ) 

Apply a time bias when invoking SetEphemerisTime method.

The bias is used only when reading from NAIF kernels. It is added to the ephermeris time passed into SetEphemerisTime and then the body position is read from the NAIF kernels and returned. When the cache is loaded from a table the bias is ignored as it is assumed to have already been applied. If this method is never called the default bias is 0.0 seconds.

Parameters:
timeBias time bias in seconds
const std::vector< double > & Isis::SpicePosition::Velocity (  ) 

Return the current J2000 velocity.

Return the velocity vector if available.

Returns:
The velocity vector evaluated at the current time

Referenced by Isis::Spice::InstrumentVelocity(), Isis::MiniRF::MiniRF(), and Isis::RadarGroundMap::RadarGroundMap().

std::vector< double > Isis::SpicePosition::VelocityPartial ( SpicePosition::PartialType  partialVar,
int  coeffIndex 
)

Compute the derivative of the velocity with respect to the specified variable.

The velocity is the derivative of the coordinate with respect to time.

coordinate = C0 + C1*t + C2*t**2 + ... +Cn*t**n , where t = (time - p_basetime) / p_timeScale. velocity = C1 + 2*C2*t + ... + n*Cn*t**(n-1) partial(velocity) with respect to C0 = 0. partial(velocity) with respect to C1 = 1/p_timeScale. partial(velocity) with respect to C2 = 2*t/p_timeScale partial(velocity) with respect to CN = n*t**(n-1)/p_timeScale

Parameters:
partialVar Designated variable of the partial derivative
Returns:
Derivative of j2000 velocity vector calculated with respect to partialVar

Referenced by Isis::RadarGroundMap::GetdXYdPosition().


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