USGS

Isis 3.0 Object Programmers' Reference

Home

Isis::SpicePosition Class Reference

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

#include <SpicePosition.h>

Inheritance diagram for Isis::SpicePosition:
Inheritance graph
Collaboration diagram for Isis::SpicePosition:
Collaboration graph

Public Types

enum  Source {
  Spice, Memcache, HermiteCache, PolyFunction,
  PolyFunctionOverHermiteConstant
}
 This enum 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.
 
double GetTimeBias () const
 Returns the value of the time bias added to ET.
 
virtual void SetAberrationCorrection (const QString &correction)
 Set the aberration correction (light time)
 
virtual QString GetAberrationCorrection () const
 Returns current state of stellar aberration correction.
 
double GetLightTime () const
 Return the light time coorection value.
 
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 QString &tableName)
 Return a table with J2000 to reference positions.
 
Table LoadHermiteCache (const QString &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 QString &tableName)
 Return a table with J2000 positions.
 
bool IsCached () const
 Is this position cached.
 
void SetPolynomial (const Source type=PolyFunction)
 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, const Source type=PolyFunction)
 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.
 
std::vector< double > HermiteCoordinate ()
 This method returns the Hermite coordinate for the current time for PolyFunctionOverHermiteConstant functions.
 

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.
 
virtual 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.
 
void SetEphemerisTimePolyFunctionOverHermiteConstant ()
 This is a protected method that is called by SetEphemerisTime() when Source type is PolyFunctionOverHermiteConstant.
 
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.
 
 SpicePosition (int targetCode, int observerCode, bool swapObserverTarget)
 Constructor for swapping observer/target order.
 
int getObserverCode () const
 Returns observer code.
 
int getTargetCode () const
 Returns target code.
 
double getAdjustedEphemerisTime () const
 Returns adjusted ephemeris time.
 
void computeStateVector (double et, int target, int observer, const QString &refFrame, const QString &abcorr, double state[6], bool &hasVelocity, double &lightTime) const
 Computes the state vector of the target w.r.t observer.
 
void setStateVector (const double state[6], const bool &hasVelocity)
 Sets the state of target relative to observer.
 
void setLightTime (const double &lightTime)
 Inheritors can set the light time if indicated.
 

Private Member Functions

void init (int targetCode, int observerCode, const bool &swapObserverTarget=false)
 Internal initialization of the object support observer/target swap.
 
void ClearCache ()
 Removes the entire cache from memory.
 
void LoadTimeCache ()
 Load the time cache.
 
void CacheLabel (Table &table)
 Add labels to a SpicePosition table.
 
double ComputeVelocityInTime (PartialType var)
 Compute the velocity with respect to time instead of scaled time.
 

Private Attributes

int p_targetCode
 target body code
 
int p_observerCode
 observer body code
 
double p_timeBias
 iTime bias when reading kernels
 
QString p_aberrationCorrection
 Light time correction to apply.
 
double p_et
 Current ephemeris time.
 
std::vector< double > p_coordinate
 J2000 position at time et.
 
std::vector< double > p_velocity
 J2000 velocity at time et.
 
NumericalApproximationp_xhermite
 Hermite spline for x coordinate if Source == HermiteCache.
 
NumericalApproximationp_yhermite
 Hermite spline for y coordinate if Source == HermiteCache.
 
NumericalApproximationp_zhermite
 Hermite spline for z coordinate if Source == HermiteCache.
 
Source p_source
 Enumerated value for the location of the SPK information used.
 
std::vector< double > p_cacheTime
 iTime for corresponding position
 
std::vector< std::vector
< double > > 
p_cache
 Cached positions.
 
std::vector< std::vector
< double > > 
p_cacheVelocity
 Cached velocities.
 
std::vector< double > p_coefficients [3]
 Coefficients of polynomials fit to 3 coordinates.
 
double p_baseTime
 Base time used in fit equations.
 
double p_timeScale
 Time scale used in fit equations.
 
bool p_degreeApplied
 Flag indicating whether or not a polynomial.
 
int p_degree
 Degree of polynomial function fit to the coordinates of the position.
 
double p_fullCacheStartTime
 Original start time of the complete cache after spiceinit.
 
double p_fullCacheEndTime
 Original end time of the complete cache after spiceinit.
 
double p_fullCacheSize
 Orignial size of the complete cache after spiceinit.
 
bool p_hasVelocity
 Flag to indicate velocity is available.
 
OverrideType p_override
 Time base and scale override options;.
 
double p_overrideBaseTime
 Value set by caller to override computed base time.
 
double p_overrideTimeScale
 Value set by caller to override computed time scale.
 
bool m_swapObserverTarget
 
double m_lt
 !< Swap traditional order
 

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

2003-12-01 Jeff Anderson Original Version derived from Spice class

2006-03-23 Jeff Anderson Added check to SetEphemeris to return if the time did not change. Should speed up line scan cameras

2007-07-10 Debbie A. Cook Added else to method SetAberrationCorrection to separate error section from the rest of the code

2007-08-24 Debbie A. Cook Added members p_coefficients, enums PartialType and Coefficient, and methods ReloadCache, SetPolynomial, GetPolynomial, ComputeBaseTime, DPolynomial, and CoordinatePartial to support solving for instrument position in BundleAdjust

2008-01-10 Debbie A. Cook The function was changed from Parabola to Poly1D. New methods GetBaseTime and SetOverrideBaseTime were added

2008-02-07 Jeannie Walldren Poly1D was changed to PolynomialUnivariate

2008-02-10 Debbie A. Cook Removed the (-1.) in case A of DPolynomial since it was not actually part of the position derivation but how it was being applied in BundleAdjust.

2008-06-18 Steven Lambright Fixed documentation

2008-06-25 Debbie A. Cook Added members p_velocity, p_cacheVelocity, and p_hasVelocity; also added methods Velocity() and HasVelocity() to support miniRF.

2008-06-26 Debbie A. Cook Replaced Naif error handling calls with Isis NaifStatus

2009-08-03 Jeannie Walldren - Added methods ReloadCache( table ), HermiteIndices(), Memcache2HermiteCache() to allow for downsized instrument position tables in labels of Isis cubes and added methods SetEphemerisTimeSpice(), SetEphemerisTimeHermiteCache(), and SetEphemerisTimeMemcache() to make software more readable.

2009-08-14 Debbie A. Cook - Corrected loop counter in HermiteIndices

2009-08-27 Jeannie Walldren - Added documentation.

2009-10-20 Debbie A. Cook - Corrected calculation of extremum in ReloadCache

2009-11-06 Debbie A. Cook - Added velocity partial derivative method

2009-12-15 Debbie A. Cook - Changed enumerated partial types and argument list for CoordinatePartial and VelocityPartial

2010-03-19 Debbie A. Cook - Added argument coeffIndex to method Velocity Partial

2010-12-07 Steven Lambright - Moved the cubic hermite splines to member scope for efficiency. It was a significant overhead to keep reconstructing these. Created ClearCache() to help increase code reusability.

2011-02-12 Debbie A. Cook - Added PolyFunction to source types and new method SetEphemerisTimePolyFunction to support the new type. Also added method SetPolynomialDegree to allow the degree of the polynomials fit to the position coordinates to be changed. The fit polynomial was changed from a fixed 2nd order polynomial to an nth degree polynomial with one independent variable. Added private class members p_fullCacheStartTime, p_fullCacheEndTime, and p_fullCacheSize. Added p_timeScale, GetBaseTime(), SetOverrideBaseTime, GetTimeScale(), Extrapolate(double timeEt), and ComputeVelocityInTime(PartialType var). The function was changed from a Parabola to Polynomial1Variable, now called PolynomialUnivariate.

2011-02-14 Debbie A. Cook - Corrected previous update for Hermite case. Created a special member, p_overrideScale, for Hermite case. Also removed obsolete enum Coefficient.

2011-02-17 Debbie A. Cook - Corrected missed problem with degree forced to be 2 and corrected calculation of velocity partial

2011-02-22 Debbie A. Cook - Corrected Extrapolation method

2011-02-28 Debbie A. Cook - Fixed typo in LoadCache and potential memory problem in SetPolynomial

2011-04-10 Debbie A. Cook - Added GetSource method to support spkwriter.

2011-05-27 Debbie A. Cook - Renamed old ReloadCache method for converting polynomial functions to Hermite cubic splines to LoadHermiteCache for use with spkwriter.

2012-03-31 Debbie A. Cook - Programmer notes: Added new interpolation algorithm, PolyFunctionOverHermiteConstant, and new supporting method, SetEphemerisTimePolyOverHermiteConstant. Also added argument to SetPolynomial methods for type. LineCache and ReloadCache were modified to allow other function types beyond the PolyFunction. Added new method to return Hermite coordinate.

2012-10-25 Jeannie Backer - Brought class closer to Isis3 standards: Ordered includes in cpp file, replaced quotation marks with angle braces in 3rd party includes, fixed history indentation. References #1181.

2012-10-31 Kris Becker - Added implementation for swapping of observer/target and light time correction to surface. Added a new, special protected constructor; virtualized the SetEphemerisTimeSpice method which allows for specialized determination of light time corrections; generalized the retrieval of body state vectors (computeStateVector); added many new, mostly protected, ssupport routines for this generalization. Fixes (mostly) #0909, #1136 and #1223.

2013-03-28 Debbie A. Cook - Fixed bug causing jigsaw to fail with a singular definite matrix on radar data. The bug was in the method VelocityPartial and occurred when the et = baseTime and the coeffIndex was 0. This caused the derivative equation to be 0 * 0 ** -1. This update fixes issue #1582.

Definition at line 174 of file SpicePosition.h.

Member Enumeration Documentation

This enum indicates the status of the object.

The class expects functions to be after MemCache in the list.

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.

PolyFunctionOverHermiteConstant 

Object is reading from splined.

Definition at line 181 of file SpicePosition.h.

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
targetCodeValid naif body name.
observerCodeValid naif body name.

Definition at line 29 of file SpicePosition.cpp.

References init().

Isis::SpicePosition::~SpicePosition ( )
virtual

Destructor.

Free the memory allocated by this SpicePosition instance.

Definition at line 126 of file SpicePosition.cpp.

References ClearCache().

Isis::SpicePosition::SpicePosition ( int  targetCode,
int  observerCode,
bool  swapObserverTarget 
)
protected

Constructor for swapping observer/target order.

This specialized constructor is provided to expressly support swap of observer/target order in the NAIF spkez_c/spkezp_c routines when determining the state of spacecraft and target body.

Traditionally, ISIS3 (and ISIS2!) has had the target/observer order swapped improperly to determine state vector between the two bodies. This constructor provides a transaitional path to correct this in a controlled way. See SpacecraftPosition for details.

Note that the targetCode and observerCode are always provided in the same order as the orignal constructor (i.e., targetCode=s/c, observerCode=planet) as the swapObserverTarget boolean parameter provides the means to properly implement the swap internally.

It is not as simple as simply swapping the codes in the constructor as internal representation of the state vector is in body fixed state of the target body relative to observer. If the swap is invoked, the state vector must be inverted.

Author
2012-10-29 Kris Becker
Parameters
targetCodeTraditional s/c code
observerCodeTraditional target/planet code
swapObserverTargetTrue indicates implement the observer/target swap, false will invoke preexisting behavior

Definition at line 62 of file SpicePosition.cpp.

References init().

Member Function Documentation

Table Isis::SpicePosition::Cache ( const QString &  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
tableNameName of the table to create and return
History:

2009-08-03 Jeannie Walldren - Added CacheType keyword to output table. This is based on p_source and will be read by LoadCache(Table).

2011-01-05 Debbie A. Cook - Added PolyFunction

Definition at line 479 of file SpicePosition.cpp.

References _FILEINFO_, CacheLabel(), Isis::TableField::Double, Isis::IException::Io, LineCache(), Memcache2HermiteCache(), p_baseTime, p_cache, p_cacheTime, p_cacheVelocity, p_coefficients, p_degree, p_fullCacheSize, p_hasVelocity, p_source, p_timeScale, PolyFunction, and PolyFunctionOverHermiteConstant.

Referenced by LineCache(), and LoadHermiteCache().

void Isis::SpicePosition::CacheLabel ( Table table)
private

Add labels to a SpicePosition table.

Return a table containing the labels defining the position table.

Parameters
TableTable to receive labels

Definition at line 589 of file SpicePosition.cpp.

References HermiteCache, Isis::Blob::Label(), Memcache, p_fullCacheEndTime, p_fullCacheSize, p_fullCacheStartTime, p_source, and Isis::toString().

Referenced by Cache().

void Isis::SpicePosition::ClearCache ( )
private

Removes the entire cache from memory.

Definition at line 1600 of file SpicePosition.cpp.

References p_cache, p_cacheTime, p_cacheVelocity, p_xhermite, p_yhermite, and p_zhermite.

Referenced by LoadHermiteCache(), ReloadCache(), and ~SpicePosition().

void Isis::SpicePosition::ComputeBaseTime ( )

Compute the base time using cached times.

Definition at line 1035 of file SpicePosition.cpp.

References p_baseTime, p_cacheTime, p_override, p_overrideBaseTime, p_overrideTimeScale, and p_timeScale.

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

void Isis::SpicePosition::computeStateVector ( double  et,
int  target,
int  observer,
const QString &  refFrame,
const QString &  abcorr,
double  state[6],
bool &  hasVelocity,
double &  lightTime 
) const
protected

Computes the state vector of the target w.r.t observer.

This method computes the state vector of the target that is relative of the observer. It first attempts to retrieve with velocity vectors. If that fails, it makes an additional attempt to get the state w/o velocity vectors. The final result is indicated by the hasVelocity parameter.

The parameters to this routine are the same as the NAIF spkez_c/spkezp_c routines. See http://naif.jpl.nasa.gov/pub/naif/toolkit_docs/C/cspice/spkez_c.html for complete description.

Note that this routine does not actually affect the internals of this object (hence the constness of the method). It is up to the caller to apply the results and potentially handle swapping of observer/target and light time correction. See SpacecraftPosition for additional details on this application of the results.

Author
2012-10-29 Kris Becker
Parameters
etTime to compute state vector for
targetNAIF target code
observerNAIF observer code
refFrameReference frame to express coordinates in (e.g., J2000)
abcorrStellar aberration correction option
stateReturns the 6-element state vector in target body fixed coordinates
hasVelocityReturns knowledge of whether the velocity vector is valid
lightTimeReturns the light time correct resulting from the request if applicable

Definition at line 1870 of file SpicePosition.cpp.

Referenced by Isis::SpacecraftPosition::SetEphemerisTimeSpice(), and SetEphemerisTimeSpice().

double Isis::SpicePosition::ComputeVelocityInTime ( PartialType  var)
private

Compute the velocity with respect to time instead of scaled time.

Parameters
coefA SpicePosition polynomial function partial type
History:
2011-02-12 Debbie A. Cook - Original version.

Definition at line 1731 of file SpicePosition.cpp.

References p_baseTime, p_coefficients, p_degree, p_et, and p_timeScale.

Referenced by SetEphemerisTimePolyFunction().

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
partialVarDesignated variable of the partial derivative
Returns
Derivative of j2000 vector calculated with polynomial with respect to partialVar

Definition at line 1086 of file SpicePosition.cpp.

References DPolynomial().

Referenced by Isis::CameraGroundMap::GetdXYdPosition(), and Isis::RadarGroundMap::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
coeffIndexIndex of coefficient to differentiate with respect to
Returns
The derivative evaluated at the current time

Definition at line 1163 of file SpicePosition.cpp.

References _FILEINFO_, p_baseTime, p_degree, p_et, p_timeScale, and Isis::IException::Programmer.

Referenced by CoordinatePartial().

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

Return the current ephemeris time.

Definition at line 205 of file SpicePosition.h.

References p_et.

Referenced by getAdjustedEphemerisTime().

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]timeEtThe time of the position to be extrapolated
[out]Anextrapolated position at the input time

Definition at line 1754 of file SpicePosition.cpp.

References p_coordinate, p_et, p_hasVelocity, and p_velocity.

QString Isis::SpicePosition::GetAberrationCorrection ( ) const
virtual

Returns current state of stellar aberration correction.

The aberration correction is the value of the parameter that will be provided to the spkez_c/spkezp_c routines when determining the targer/observer vector state. See SetAberrationCorrection() for valid values.

Author
2012-10-29 Kris Becker

Reimplemented in Isis::SpacecraftPosition.

Definition at line 206 of file SpicePosition.cpp.

References p_aberrationCorrection.

Referenced by SetEphemerisTimeSpice().

double Isis::SpicePosition::getAdjustedEphemerisTime ( ) const
protected

Returns adjusted ephemeris time.

This method returns the actual ephemeris time adjusted by a specifed time bias. The bias typically comes explicitly from the camera model but could be adjusted by the environment they are utlized in.

Author
2012-10-29 Kris Becker
Returns
double Adjusted ephemeris time with time bias applied

Definition at line 1833 of file SpicePosition.cpp.

References EphemerisTime(), and GetTimeBias().

Referenced by Isis::SpacecraftPosition::SetEphemerisTimeSpice(), and SetEphemerisTimeSpice().

double Isis::SpicePosition::GetBaseTime ( )
inline

Return the base time for the position.

Definition at line 263 of file SpicePosition.h.

References p_baseTime.

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

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

Return the light time coorection value.

This method returns the light time correction of the last call to SetEphemerisTimeSpice. This is the actual time after adjustments that has resulting int the adjustment of the target body. The observer position should be unchanged.

Author
2012-10-29 Kris Becker
Returns
double Returns the light time determined from the last call to SetEphemerisTime[Spice].

Definition at line 223 of file SpicePosition.cpp.

References m_lt.

int Isis::SpicePosition::getObserverCode ( ) const
protected

Returns observer code.

This methods returns the proper observer code as specified in constructor. Code has been subjected to swapping if requested.

Author
2012-10-29 Kris Becker
Returns
int NAIF code of observer

Definition at line 1803 of file SpicePosition.cpp.

References p_observerCode.

Referenced by Isis::SpacecraftPosition::SetEphemerisTimeSpice(), and SetEphemerisTimeSpice().

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]XCCoefficients of fit to first coordinate of position
[out]YCCoefficients of fit to second coordinate of position
[out]ZCCoefficients of fit to third coordinate of position

Definition at line 1022 of file SpicePosition.cpp.

References p_coefficients.

Referenced by Isis::BundleAdjust::applyParameterCorrections_CHOLMOD(), Isis::BundleAdjust::applyParameterCorrections_SPECIALK(), Isis::BundleAdjust::OutputImagesCSV(), Isis::BundleAdjust::OutputNoErrorPropagation(), Isis::BundleAdjust::OutputWithErrorPropagation(), Isis::BundleAdjust::Solve(), Isis::BundleAdjust::SolveCholesky(), and Isis::BundleAdjust::Update().

Source Isis::SpicePosition::GetSource ( )
inline

Return the source of the position.

Definition at line 256 of file SpicePosition.h.

References p_source.

int Isis::SpicePosition::getTargetCode ( ) const
protected

Returns target code.

This methods returns the proper target code as specified in constructor. Code has been subjected to swapping if requested.

Author
2012-10-29 Kris Becker
Returns
int NAIF code for target body

Definition at line 1817 of file SpicePosition.cpp.

References p_targetCode.

Referenced by Isis::SpacecraftPosition::SetEphemerisTimeSpice(), and SetEphemerisTimeSpice().

double Isis::SpicePosition::GetTimeBias ( ) const

Returns the value of the time bias added to ET.

Author
2012-10-29 Kris Becker
Returns
double Returns time bias for current object

Definition at line 154 of file SpicePosition.cpp.

References p_timeBias.

Referenced by getAdjustedEphemerisTime().

double Isis::SpicePosition::GetTimeScale ( )
inline

Return the time scale for the position.

Definition at line 270 of file SpicePosition.h.

References p_timeScale.

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

bool Isis::SpicePosition::HasVelocity ( )
inline

Return the flag indicating whether the velocity exists.

Definition at line 220 of file SpicePosition.h.

References p_hasVelocity.

std::vector< double > Isis::SpicePosition::HermiteCoordinate ( )

This method returns the Hermite coordinate for the current time for PolyFunctionOverHermiteConstant functions.

See Also
SetEphemerisTime()
History:
2012-02-05 Debbie A. Cook - Original version

Definition at line 1774 of file SpicePosition.cpp.

References _FILEINFO_, p_coordinate, p_source, PolyFunctionOverHermiteConstant, Isis::IException::Programmer, and SetEphemerisTimeHermiteCache().

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
toleranceMaximum error allowed between NAIF kernel coordinate values and values interpolated by the Hermite spline.
indexListVector 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.
History:

2009-08-03 Jeannie Walldren - Original version.

2009-08-14 Debbie A. Cook - Corrected indexing error in loop.

Definition at line 1531 of file SpicePosition.cpp.

References _FILEINFO_, Isis::NumericalApproximation::AddCubicHermiteDeriv(), Isis::NumericalApproximation::AddData(), Isis::NumericalApproximation::CubicHermite, Isis::NumericalApproximation::Evaluate(), Isis::IException::Io, p_baseTime, p_cache, p_cacheTime, p_cacheVelocity, p_hasVelocity, and p_timeScale.

Referenced by Memcache2HermiteCache().

void Isis::SpicePosition::init ( int  targetCode,
int  observerCode,
const bool &  swapObserverTarget = false 
)
private

Internal initialization of the object support observer/target swap.

This initailizer provides options to swap observer/target properly.

Author
2012-10-29 Kris Becker
Parameters
targetCodeTraditional s/c code
observerCodeTraditional target/planet code
swapObserverTargetTrue indicates implement the observer/target swap, false will invoke preexisting behavior

Definition at line 79 of file SpicePosition.cpp.

References m_lt, p_aberrationCorrection, p_baseTime, p_coefficients, p_coordinate, p_degree, p_degreeApplied, p_et, p_fullCacheEndTime, p_fullCacheSize, p_fullCacheStartTime, p_hasVelocity, p_observerCode, p_override, p_source, p_targetCode, p_timeBias, p_timeScale, p_velocity, p_xhermite, p_yhermite, p_zhermite, and Spice.

Referenced by SpicePosition().

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

Is this position cached.

Definition at line 237 of file SpicePosition.h.

References p_cache.

Table Isis::SpicePosition::LineCache ( const QString &  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
tableNameName of the table to create and return
History:
2012-01-25 Debbie A. Cook - Modified error checking for p_source to allow all function sources (>=HermiteCache)

Definition at line 635 of file SpicePosition.cpp.

References _FILEINFO_, Cache(), HermiteCache, Memcache, p_source, Isis::IException::Programmer, and ReloadCache().

Referenced by 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
startTimeStarting ephemeris time in seconds for the cache
endTimeEnding ephemeris time in seconds for the cache
sizeNumber of coordinates/positions to keep in the cache

Definition at line 287 of file SpicePosition.cpp.

References _FILEINFO_, HermiteCache, LoadTimeCache(), Memcache, p_cache, p_cacheTime, p_cacheVelocity, p_coordinate, p_fullCacheEndTime, p_fullCacheSize, p_fullCacheStartTime, p_hasVelocity, p_source, p_velocity, Isis::IException::Programmer, and SetEphemerisTime().

Referenced by LoadCache(), and ReloadCache().

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
timesingle ephemeris time in seconds to cache

Definition at line 335 of file SpicePosition.cpp.

References LoadCache().

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
tableAn ISIS table blob containing valid J2000 coordinate/time values
History:

2009-08-03 Jeannie Walldren - Reads CacheType keyword from table and sets p_source. If no CacheType keyword, we know this is an older image, so set p_source to Memcache.

2011-01-05 Debbie A. Cook - Added PolyFunction type

2011-04-08 Debbie A. Cook - Corrected loop counter in PolyFunction section to only go up to table.Records() - 1

Definition at line 362 of file SpicePosition.cpp.

References _FILEINFO_, Isis::TableRecord::Fields(), Isis::PvlObject::findKeyword(), Isis::PvlObject::hasKeyword(), HermiteCache, Isis::IException::Io, Isis::Blob::Label(), Memcache, Isis::Blob::Name(), p_cache, p_cacheTime, p_cacheVelocity, p_fullCacheEndTime, p_fullCacheSize, p_fullCacheStartTime, p_hasVelocity, p_override, p_overrideTimeScale, p_source, PolyFunction, Isis::IException::Programmer, Isis::Table::Records(), SetOverrideBaseTime(), SetPolynomial(), SetPolynomialDegree(), and Isis::toDouble().

Table Isis::SpicePosition::LoadHermiteCache ( const QString &  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.

Definition at line 720 of file SpicePosition.cpp.

References _FILEINFO_, Cache(), ClearCache(), Isis::BasisFunction::Coefficient(), Isis::BasisFunction::Evaluate(), HermiteCache, p_baseTime, p_cache, p_cacheTime, p_cacheVelocity, p_coefficients, p_coordinate, p_degree, p_et, p_fullCacheEndTime, p_fullCacheSize, p_fullCacheStartTime, p_hasVelocity, p_source, p_velocity, PolyFunction, Isis::IException::Programmer, Isis::BasisFunction::SetCoefficients(), and SetEphemerisTime().

void Isis::SpicePosition::LoadTimeCache ( )
private

Load the time cache.

This method should works with the LoadCache(startTime, endTime, size) method to load the time cache.

Definition at line 1697 of file SpicePosition.cpp.

References p_cacheTime, p_fullCacheEndTime, p_fullCacheSize, and p_fullCacheStartTime.

Referenced by LoadCache(), and ReloadCache().

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
toleranceMaximum error allowed between NAIF kernel coordinate values and values interpolated by the Hermite spline.
History:
2009-08-03 Jeannie Walldren - Original version.

Definition at line 1471 of file SpicePosition.cpp.

References _FILEINFO_, ComputeBaseTime(), HermiteCache, HermiteIndices(), Memcache, p_cache, p_cacheTime, p_cacheVelocity, p_override, p_overrideTimeScale, p_source, and Isis::IException::Programmer.

Referenced by Cache().

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.

History:
2012-01-25 Debbie A. Cook - Modified error checking for type to allow all function types (>=HermiteCache)

Definition at line 660 of file SpicePosition.cpp.

References _FILEINFO_, Isis::NaifStatus::CheckErrors(), HermiteCache, LoadTimeCache(), Memcache, p_cache, p_cacheTime, p_cacheVelocity, p_coordinate, p_degree, p_et, p_fullCacheSize, p_source, p_velocity, Isis::IException::Programmer, and SetEphemerisTime().

Referenced by LineCache().

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
tableAn ISIS table blob containing valid J2000 coordinate/time values.
History:
2009-08-03 Jeannie Walldren - Original version.

Definition at line 1687 of file SpicePosition.cpp.

References ClearCache(), LoadCache(), p_source, and Spice.

void Isis::SpicePosition::SetAberrationCorrection ( const QString &  correction)
virtual

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
correctionThis value must be one of the following: "NONE", "LT", "LT+S", "CN", "CN+S", "XLT", "XLT+S", "XCN" and "XCN+S" where LT is a correction for planetary aberration (light time), CN is converged Newtonian, XLT is transmission case using Newtonian formulation, XCN is transmission case using converged Newtonian formuation and S a correction for stellar aberration. See http://naif.jpl.nasa.gov/pub/naif/toolkit_docs/C/cspice/spkezr_c.html for additional information on the "abcorr" parameters. If never called the default is "LT+S".
History:
2012-10-10 Kris Becker - Added support for all current abcorr parameters options in the NAIF spkezX functions.

Reimplemented in Isis::SpacecraftPosition.

Definition at line 179 of file SpicePosition.cpp.

References _FILEINFO_, p_aberrationCorrection, and Isis::IException::Programmer.

Referenced by Isis::LoHighCamera::LoHighCamera(), Isis::LoMediumCamera::LoMediumCamera(), Isis::Mariner10Camera::Mariner10Camera(), and Isis::SpacecraftPosition::SetAberrationCorrection().

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
etephemeris time in seconds
History:
2009-08-03 Jeannie Walldren - Moved code to individual methods for each Source type to make software more readable.

Definition at line 242 of file SpicePosition.cpp.

References Isis::NaifStatus::CheckErrors(), HermiteCache, Memcache, p_coordinate, p_et, p_source, PolyFunction, PolyFunctionOverHermiteConstant, SetEphemerisTimeHermiteCache(), SetEphemerisTimeMemcache(), SetEphemerisTimePolyFunction(), SetEphemerisTimePolyFunctionOverHermiteConstant(), and SetEphemerisTimeSpice().

Referenced by GetCenterCoordinate(), LoadCache(), LoadHermiteCache(), ReloadCache(), and SetPolynomial().

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()
History:
2009-08-03 Jeannie Walldren - Original version

Definition at line 1274 of file SpicePosition.cpp.

References _FILEINFO_, Isis::NumericalApproximation::AddCubicHermiteDeriv(), Isis::NumericalApproximation::AddData(), ComputeBaseTime(), Isis::NumericalApproximation::CubicHermite, Isis::NumericalApproximation::Evaluate(), Isis::NumericalApproximation::EvaluateCubicHermiteFirstDeriv(), Isis::NumericalApproximation::Extrapolate, Isis::IException::Io, p_baseTime, p_cache, p_cacheTime, p_cacheVelocity, p_coordinate, p_et, p_hasVelocity, p_override, p_overrideTimeScale, p_timeScale, p_velocity, p_xhermite, p_yhermite, and p_zhermite.

Referenced by HermiteCoordinate(), SetEphemerisTime(), and SetEphemerisTimePolyFunctionOverHermiteConstant().

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()
History:
2009-08-03 Jeannie Walldren - Original version (moved code from SetEphemerisTime() to its own method)

Definition at line 1211 of file SpicePosition.cpp.

References p_cache, p_cacheTime, p_cacheVelocity, p_coordinate, p_et, p_hasVelocity, and p_velocity.

Referenced by 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()
History:
2011-01-05 Debbie A. Cook - Original version

Definition at line 1362 of file SpicePosition.cpp.

References ComputeVelocityInTime(), Isis::BasisFunction::Evaluate(), p_baseTime, p_cacheVelocity, p_coefficients, p_coordinate, p_degree, p_et, p_hasVelocity, p_timeScale, p_velocity, and Isis::BasisFunction::SetCoefficients().

Referenced by SetEphemerisTime(), and SetEphemerisTimePolyFunctionOverHermiteConstant().

void Isis::SpicePosition::SetEphemerisTimePolyFunctionOverHermiteConstant ( )
protected

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

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

See Also
SetEphemerisTime()
History:
2012-01-25 Debbie A. Cook - Original version

Definition at line 1409 of file SpicePosition.cpp.

References p_coordinate, p_velocity, SetEphemerisTimeHermiteCache(), and SetEphemerisTimePolyFunction().

Referenced by SetEphemerisTime().

void Isis::SpicePosition::SetEphemerisTimeSpice ( )
protectedvirtual

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()
History:
2009-08-03 Jeannie Walldren - Original version (moved code from SetEphemerisTime() to its own method)

Reimplemented in Isis::SpacecraftPosition.

Definition at line 1440 of file SpicePosition.cpp.

References computeStateVector(), GetAberrationCorrection(), getAdjustedEphemerisTime(), getObserverCode(), getTargetCode(), setLightTime(), and setStateVector().

Referenced by SetEphemerisTime(), and Isis::SpacecraftPosition::SetEphemerisTimeSpice().

void Isis::SpicePosition::setLightTime ( const double &  lightTime)
protected

Inheritors can set the light time if indicated.

Definition at line 1951 of file SpicePosition.cpp.

References m_lt.

Referenced by Isis::SpacecraftPosition::SetEphemerisTimeSpice(), and SetEphemerisTimeSpice().

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]baseTimeThe baseTime to use and override the computed base time
History:
2011-04-08 Debbie A. Cook - Corrected p_override set to BaseAndScale

Definition at line 1066 of file SpicePosition.cpp.

References p_override, p_overrideBaseTime, and p_overrideTimeScale.

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

void Isis::SpicePosition::SetPolynomial ( const Source  type = PolyFunction)

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.

< Basis function fit to X

< Basis function fit to Y

< Basis function fit to Z

Definition at line 850 of file SpicePosition.cpp.

References Isis::LeastSquares::AddKnown(), Isis::BasisFunction::Coefficient(), Isis::BasisFunction::Coefficients(), ComputeBaseTime(), Isis::LineEquation::Intercept(), p_baseTime, p_cache, p_cacheTime, p_coordinate, p_degree, p_source, p_timeScale, PolyFunction, PolyFunctionOverHermiteConstant, SetEphemerisTime(), Isis::LineEquation::Slope(), and Isis::LeastSquares::Solve().

Referenced by Isis::BundleAdjust::applyParameterCorrections_CHOLMOD(), Isis::BundleAdjust::applyParameterCorrections_SPECIALK(), LoadCache(), SetPolynomialDegree(), Isis::BundleAdjust::Solve(), Isis::BundleAdjust::SolveCholesky(), and Isis::BundleAdjust::Update().

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

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]XCCoefficients of fit to X coordinate
[in]YCCoefficients of fit to Y coordinate
[in]ZCCoefficients of fit to Z coordinate
History:
2012-02-05 Debbie A. Cook - Added type argument

Definition at line 973 of file SpicePosition.cpp.

References ComputeBaseTime(), p_coefficients, p_degree, p_degreeApplied, p_et, p_source, Isis::BasisFunction::SetCoefficients(), and SetEphemerisTime().

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]degreeDegree of the polynomial to be fit

Definition at line 1631 of file SpicePosition.cpp.

References p_coefficients, p_degree, p_degreeApplied, p_fullCacheSize, and SetPolynomial().

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

void Isis::SpicePosition::setStateVector ( const double  state[6],
const bool &  hasVelocity 
)
protected

Sets the state of target relative to observer.

This method sets the state of the target (vector) relative to the observer. Note that is the only place where the swap of observer/target adjustment to the state vector is handled. All contributors to this computation must compute the state vector representing the position and velocity of the target relative to the observer where the first three components of state[] are the x-, y- and z-component cartesian coordinates of the target's position; the last three are the corresponding velocity vector.

This routine maintains the directional integrity of the state vector should the observer/target be swapped. See the documentation for the spkez_c at http://naif.jpl.nasa.gov/pub/naif/toolkit_docs/C/cspice/spkez_c.html.

Author
2012-10-29 Kris Becker
Parameters
stateState vector. First three components of this 6 element array is the body fixed coordinates of the vector from the target to the observer. The last three components are the velocity state.
hasVelocityIf true, then the velocity components of the state vector are valid, otherwise they should be ignored.

Definition at line 1921 of file SpicePosition.cpp.

References p_coordinate, p_hasVelocity, and p_velocity.

Referenced by Isis::SpacecraftPosition::SetEphemerisTimeSpice(), and SetEphemerisTimeSpice().

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
timeBiastime bias in seconds

Definition at line 143 of file SpicePosition.cpp.

References p_timeBias.

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

Definition at line 1186 of file SpicePosition.cpp.

References _FILEINFO_, p_hasVelocity, p_velocity, and Isis::IException::Programmer.

Referenced by Isis::RadarGroundMap::GetXY(), Isis::MiniRF::MiniRF(), Isis::RadarGroundMap::SetFocalPlane(), and Isis::RadarGroundMap::SetGround().

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 = (1/p_timeScale) * (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
partialVarDesignated variable of the partial derivative
Returns
Derivative of j2000 velocity vector calculated with respect to partialVar

Definition at line 1127 of file SpicePosition.cpp.

References p_baseTime, p_et, and p_timeScale.

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

Member Data Documentation

double Isis::SpicePosition::m_lt
private

!< Swap traditional order

Definition at line 354 of file SpicePosition.h.

Referenced by GetLightTime(), init(), and setLightTime().

QString Isis::SpicePosition::p_aberrationCorrection
private

Light time correction to apply.

Definition at line 319 of file SpicePosition.h.

Referenced by GetAberrationCorrection(), init(), and SetAberrationCorrection().

double Isis::SpicePosition::p_baseTime
private
std::vector<std::vector<double> > Isis::SpicePosition::p_cache
private
std::vector<double> Isis::SpicePosition::p_cacheTime
private
std::vector<std::vector<double> > Isis::SpicePosition::p_cacheVelocity
private
std::vector<double> Isis::SpicePosition::p_coefficients[3]
private

Coefficients of polynomials fit to 3 coordinates.

Definition at line 336 of file SpicePosition.h.

Referenced by Cache(), ComputeVelocityInTime(), GetPolynomial(), init(), LoadHermiteCache(), SetEphemerisTimePolyFunction(), SetPolynomial(), and SetPolynomialDegree().

int Isis::SpicePosition::p_degree
private

Degree of polynomial function fit to the coordinates of the position.

Definition at line 343 of file SpicePosition.h.

Referenced by Cache(), ComputeVelocityInTime(), DPolynomial(), init(), LoadHermiteCache(), ReloadCache(), SetEphemerisTimePolyFunction(), SetPolynomial(), and SetPolynomialDegree().

bool Isis::SpicePosition::p_degreeApplied
private

Flag indicating whether or not a polynomial.

Definition at line 340 of file SpicePosition.h.

Referenced by init(), SetPolynomial(), and SetPolynomialDegree().

double Isis::SpicePosition::p_fullCacheEndTime
private

Original end time of the complete cache after spiceinit.

Definition at line 345 of file SpicePosition.h.

Referenced by CacheLabel(), GetCenterCoordinate(), init(), LoadCache(), LoadHermiteCache(), and LoadTimeCache().

double Isis::SpicePosition::p_fullCacheSize
private

Orignial size of the complete cache after spiceinit.

Definition at line 346 of file SpicePosition.h.

Referenced by Cache(), CacheLabel(), init(), LoadCache(), LoadHermiteCache(), LoadTimeCache(), ReloadCache(), and SetPolynomialDegree().

double Isis::SpicePosition::p_fullCacheStartTime
private

Original start time of the complete cache after spiceinit.

Definition at line 344 of file SpicePosition.h.

Referenced by CacheLabel(), GetCenterCoordinate(), init(), LoadCache(), LoadHermiteCache(), and LoadTimeCache().

bool Isis::SpicePosition::p_hasVelocity
private
int Isis::SpicePosition::p_observerCode
private

observer body code

Definition at line 316 of file SpicePosition.h.

Referenced by getObserverCode(), and init().

OverrideType Isis::SpicePosition::p_override
private

Time base and scale override options;.

Definition at line 348 of file SpicePosition.h.

Referenced by ComputeBaseTime(), init(), LoadCache(), Memcache2HermiteCache(), SetEphemerisTimeHermiteCache(), and SetOverrideBaseTime().

double Isis::SpicePosition::p_overrideBaseTime
private

Value set by caller to override computed base time.

Definition at line 349 of file SpicePosition.h.

Referenced by ComputeBaseTime(), and SetOverrideBaseTime().

double Isis::SpicePosition::p_overrideTimeScale
private

Value set by caller to override computed time scale.

Definition at line 350 of file SpicePosition.h.

Referenced by ComputeBaseTime(), LoadCache(), Memcache2HermiteCache(), SetEphemerisTimeHermiteCache(), and SetOverrideBaseTime().

Source Isis::SpicePosition::p_source
private

Enumerated value for the location of the SPK information used.

Definition at line 332 of file SpicePosition.h.

Referenced by Cache(), CacheLabel(), GetSource(), HermiteCoordinate(), init(), LineCache(), LoadCache(), LoadHermiteCache(), Memcache2HermiteCache(), ReloadCache(), SetEphemerisTime(), and SetPolynomial().

int Isis::SpicePosition::p_targetCode
private

target body code

Definition at line 315 of file SpicePosition.h.

Referenced by getTargetCode(), and init().

double Isis::SpicePosition::p_timeBias
private

iTime bias when reading kernels

Definition at line 318 of file SpicePosition.h.

Referenced by GetTimeBias(), init(), and SetTimeBias().

double Isis::SpicePosition::p_timeScale
private
NumericalApproximation* Isis::SpicePosition::p_xhermite
private

Hermite spline for x coordinate if Source == HermiteCache.

Definition at line 326 of file SpicePosition.h.

Referenced by ClearCache(), init(), and SetEphemerisTimeHermiteCache().

NumericalApproximation* Isis::SpicePosition::p_yhermite
private

Hermite spline for y coordinate if Source == HermiteCache.

Definition at line 328 of file SpicePosition.h.

Referenced by ClearCache(), init(), and SetEphemerisTimeHermiteCache().

NumericalApproximation* Isis::SpicePosition::p_zhermite
private

Hermite spline for z coordinate if Source == HermiteCache.

Definition at line 330 of file SpicePosition.h.

Referenced by ClearCache(), init(), and SetEphemerisTimeHermiteCache().


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