Isis Developer Reference
Isis::SpiceRotation Class Reference

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

#include <SpiceRotation.h>

Collaboration diagram for Isis::SpiceRotation:
Collaboration graph

Public Types

enum  Source {
  Spice , Nadir , Memcache , PolyFunction ,
  PolyFunctionOverSpice , PckPolyFunction
}
 The rotation can come from one of 3 places for an Isis cube. More...
 
enum  PartialType { WRT_RightAscension , WRT_Declination , WRT_Twist }
 This enumeration indicates whether the partial derivative is taken with respect to Right Ascension, Declination, or Twist (or Rotation). More...
 
enum  DownsizeStatus { Yes , Done , No }
 Status of downsizing the cache. More...
 
enum  FrameType {
  UNKNOWN = 0 , INERTL = 1 , PCK = 2 , CK = 3 ,
  TK = 4 , DYN = 5 , BPC = 6 , NOTJ2000PCK = 7
}
 Enumeration for the frame type of the rotation. More...
 

Public Member Functions

 SpiceRotation (int frameCode)
 Construct an empty SpiceRotation class using a valid Naif frame code to set up for getting rotation from Spice kernels.
 
 SpiceRotation (int frameCode, int targetCode)
 Construct an empty SpiceRotation object using valid Naif frame code and.
 
 SpiceRotation (const SpiceRotation &rotToCopy)
 Construct a SpiceRotation object by copying from an existing one.
 
virtual ~SpiceRotation ()
 Destructor for SpiceRotation object.
 
void SetFrame (int frameCode)
 Change the frame to the given frame code.
 
int Frame ()
 Accessor method that returns the frame code.
 
void SetTimeBias (double timeBias)
 Apply a time bias when invoking SetEphemerisTime method.
 
void SetEphemerisTime (double et)
 Return the J2000 to reference frame quaternion at given time.
 
double EphemerisTime () const
 Accessor method to get current ephemeris time.
 
std::vector< double > GetCenterAngles ()
 Return the camera angles at the center time of the observation.
 
std::vector< double > Matrix ()
 Return the full rotation TJ as a matrix.
 
std::vector< double > AngularVelocity ()
 Accessor method to get the angular velocity.
 
std::vector< double > ConstantRotation ()
 Return the constant 3x3 rotation TC matrix as a quaternion.
 
std::vector< double > & ConstantMatrix ()
 Return the constant 3x3 rotation TC matrix as a vector of length 9.
 
void SetConstantMatrix (std::vector< double > constantMatrix)
 Set the constant 3x3 rotation TC matrix from a vector of length 9.
 
std::vector< double > TimeBasedRotation ()
 Return time-based 3x3 rotation CJ matrix as a quaternion.
 
std::vector< double > & TimeBasedMatrix ()
 Return time-based 3x3 rotation CJ matrix as a vector of length 9.
 
void SetTimeBasedMatrix (std::vector< double > timeBasedMatrix)
 Set the time-based 3x3 rotation CJ matrix from a vector of length 9.
 
std::vector< double > J2000Vector (const std::vector< double > &rVec)
 Given a direction vector in the reference frame, return a J2000 direction.
 
std::vector< AnglepoleRaCoefs ()
 Return the coefficients used to calculate the target body pole ra.
 
std::vector< AnglepoleDecCoefs ()
 Return the coefficients used to calculate the target body pole dec.
 
std::vector< AnglepmCoefs ()
 Return the coefficients used to calculate the target body prime meridian.
 
std::vector< double > poleRaNutPrecCoefs ()
 Return the coefficients used to calculate the target body pole ra nut/prec coefficients.
 
std::vector< double > poleDecNutPrecCoefs ()
 Return the coefficients used to calculate the target body pole dec nut/prec coefficients.
 
std::vector< double > pmNutPrecCoefs ()
 Return the coefficients used to calculate the target body pm nut/prec coefficients.
 
std::vector< AnglesysNutPrecConstants ()
 Return the constants used to calculate the target body system nut/prec angles.
 
std::vector< AnglesysNutPrecCoefs ()
 Return the coefficients used to calculate the target body system nut/prec angles.
 
std::vector< double > ReferenceVector (const std::vector< double > &jVec)
 Given a direction vector in J2000, return a reference frame direction.
 
std::vector< double > EvaluatePolyFunction ()
 Evaluate the polynomial fit function for the three pointing angles for the current ephemeris time.
 
void loadPCFromSpice (int CenterBodyCode)
 Initialize planetary orientation constants from Spice PCK.
 
void loadPCFromTable (const PvlObject &Label)
 Initialize planetary orientation constants from an cube body rotation label.
 
void MinimizeCache (DownsizeStatus status)
 Set the downsize status to minimize cache.
 
void LoadCache (double startTime, double endTime, int size)
 Cache J2000 rotation quaternion over a time range.
 
void LoadCache (double time)
 Cache J2000 to frame rotation for a time.
 
void LoadCache (Table &table)
 Cache J2000 rotations using a table file.
 
void LoadCache (nlohmann::json &isd)
 Load the cached data from an ALE ISD.
 
Table LineCache (const QString &tableName)
 Return a table with J2000 to reference rotations.
 
void ReloadCache ()
 Cache J2000 rotation over existing cached time range using polynomials.
 
Table Cache (const QString &tableName)
 Return a table with J2000 to reference rotations.
 
void CacheLabel (Table &table)
 Add labels to a SpiceRotation table.
 
void LoadTimeCache ()
 Load the time cache.
 
std::vector< double > Angles (int axis3, int axis2, int axis1)
 Return the camera angles (right ascension, declination, and twist) for the time-based matrix CJ.
 
void SetAngles (std::vector< double > angles, int axis3, int axis2, int axis1)
 Set the rotation angles (phi, delta, and w) for the current time to define the time-based matrix CJ.
 
bool IsCached () const
 Checks if the cache is empty.
 
void SetPolynomial (const Source type=PolyFunction)
 Set the coefficients of a polynomial fit to each of the three camera angles for the time period covered by the cache, angle = a + bt + ct**2, where t = (time - p_baseTime)/ p_timeScale.
 
void SetPolynomial (const std::vector< double > &abcAng1, const std::vector< double > &abcAng2, const std::vector< double > &abcAng3, const Source type=PolyFunction)
 Set the coefficients of a polynomial fit to each of the three camera angles 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, and n = p_degree.
 
void usePckPolynomial ()
 Set the coefficients of a polynomial fit to each of the three planet angles for the time period covered by the cache, ra = ra0 + ra1*t + ra2*t**2 + raTrig, dec = dec0 + dec1*t + dec2*t**2 + decTrig, pm = pm0 + pm1*d + pm2*t**2 + pmTrig, where t = time / (seconds per day).
 
void setPckPolynomial (const std::vector< Angle > &raCoeff, const std::vector< Angle > &decCoeff, const std::vector< Angle > &pmCoeff)
 Set the coefficients of a polynomial fit to each of the three planet angles for the time period covered by the cache, ra = ra0 + ra1*t + ra2*t**2 + raTrig, dec = dec0 + dec1*t + dec2*t**2 + decTrig, pm = pm0 + pm1*d + pm2*t**2 + pmTrig, where t = time / (seconds per day).
 
void GetPolynomial (std::vector< double > &abcAng1, std::vector< double > &abcAng2, std::vector< double > &abcAng3)
 Return the coefficients of a polynomial fit to each of the three camera angles 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 and n = p_degree.
 
void getPckPolynomial (std::vector< Angle > &raCoeff, std::vector< Angle > &decCoeff, std::vector< Angle > &pmCoeff)
 Return the coefficients of a polynomial fit to each of the three planet angles.
 
void SetPolynomialDegree (int degree)
 Set the degree of the polynomials to be fit to the three camera angles 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, and n = p_degree.
 
Source GetSource ()
 Accessor method to get the rotation source.
 
void SetSource (Source source)
 Resets the source of the rotation to the given value.
 
void ComputeBaseTime ()
 Compute the base time using cached times.
 
FrameType getFrameType ()
 Accessor method to get the rotation frame type.
 
double GetBaseTime ()
 Accessor method to get the rotation base time.
 
double GetTimeScale ()
 Accessor method to get the rotation time scale.
 
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 save base time and polynomials for the angles.
 
void SetCacheTime (std::vector< double > cacheTime)
 
double DPolynomial (const int coeffIndex)
 Evaluate the derivative of the fit polynomial defined by the given coefficients with respect to the coefficient at the given index, at the current time.
 
double DPckPolynomial (PartialType partialVar, const int coeffIndex)
 Evaluate the derivative of the fit polynomial defined by the given coefficients with respect to the coefficient at the given index, at the current time.
 
std::vector< double > toJ2000Partial (const std::vector< double > &lookT, PartialType partialVar, int coeffIndex)
 Given a direction vector in the reference frame, compute the derivative with respect to one of the coefficients in the angle polynomial fit equation of a vector rotated from the reference frame to J2000.
 
std::vector< double > ToReferencePartial (std::vector< double > &lookJ, PartialType partialVar, int coeffIndex)
 Compute the derivative with respect to one of the coefficients in the angle polynomial fit equation of a vector rotated from J2000 to a reference frame.
 
void DCJdt (std::vector< double > &dRJ)
 Compute the derivative of the 3x3 rotation matrix CJ with respect to time.
 
double WrapAngle (double compareAngle, double angle)
 Wrap the input angle to keep it within 2pi radians of the angle to compare.
 
void SetAxes (int axis1, int axis2, int axis3)
 Set the axes of rotation for decomposition of a rotation matrix into 3 angles.
 
std::vector< double > GetFullCacheTime ()
 Return full listing (cache) of original time coverage requested.
 
void FrameTrace (double et)
 Compute frame trace chain from target frame to J2000.
 
std::vector< int > ConstantFrameChain ()
 Accessor method to get the frame chain for the constant part of the rotation (ends in target)
 
std::vector< int > TimeFrameChain ()
 Accessor method to get the frame chain for the rotation (begins in J2000).
 
void InitConstantRotation (double et)
 Initialize the constant rotation.
 
bool HasAngularVelocity ()
 Checks whether the rotation has angular velocities.
 
void ComputeAv ()
 Compute the angular velocity from the time-based functions fit to the pointing angles This method computes omega = angular velocity matrix, and extracts the angular velocity.
 
std::vector< double > Extrapolate (double timeEt)
 Extrapolate pointing for a given time assuming a constant angular velocity.
 
void checkForBinaryPck ()
 Check loaded pck to see if any are binary and set frame type to indicate binary pck.
 
int cacheSize ()
 

Protected Member Functions

void SetFullCacheParameters (double startTime, double endTime, int cacheSize)
 Set the full cache time parameters.
 
void setEphemerisTimeMemcache ()
 Updates rotation state based on the rotation cache.
 
void setEphemerisTimeNadir ()
 When setting the ephemeris time, uses spacecraft nadir source to update the rotation state.
 
void setEphemerisTimeSpice ()
 When setting the ephemeris time, updates the rotation state based on data read directly from NAIF kernels using NAIF Spice routines.
 
void setEphemerisTimePolyFunction ()
 When setting the ephemeris time, updates the rotation according to a polynomial that defines the three camera angles and angular velocity, if available.
 
void setEphemerisTimePolyFunctionOverSpice ()
 When setting the ephemeris time, updates the rotation state based on a polynomial fit over spice kernel data.
 
void setEphemerisTimePckPolyFunction ()
 When setting the ephemeris time, updates the rotation state based on the PcK polynomial.
 

Protected Attributes

std::vector< double > p_cacheTime
 iTime for corresponding rotation
 
int p_degree
 Degree of fit polynomial for angles.
 
int p_axis1
 Axis of rotation for angle 1 of rotation.
 
int p_axis2
 Axis of rotation for angle 2 of rotation.
 
int p_axis3
 Axis of rotation for angle 3 of rotation.
 
ale::Orientations * m_orientation
 

Detailed Description

Obtain SPICE rotation information for a body.

This class will obtain the rotation from J2000 to a particular reference frame, for example the rotation from J2000 to MOC NA.

It is essentially used to convert position vectors from one frame to another, making it is a C++ wrapper to the NAIF routines pxform_c and mxv or mtxv. Therefore, appropriate NAIF kernels are expected to be loaded prior to using this class. A position can be returned in either the J2000 frame or the selected reference frame. See NAIF required reading for more information regarding this subject at ftp://naif.jpl.nasa.gov/pub/naif/toolkit_docs/C/req/spk.html

An important functionality of this class is the ability to cache the rotations so they do not have to be constantly read from the NAIF kernels and they can be more conveniently updated. Once the data is cached, the NAIF kernels can be unloaded. If the rotation has a fixed part and a time- based part, the rotation is computed and stored in those two parts.

Author
2005-12-01 Debbie A. Cook

Member Enumeration Documentation

◆ DownsizeStatus

Status of downsizing the cache.

Enumerator
Yes 

Downsize the cache.

Done 

Cache is downsized.

No 

Do not downsize the cache.

◆ FrameType

Enumeration for the frame type of the rotation.

Enumerator
UNKNOWN 

Isis specific code for unknown frame type.

INERTL 

See Naif Frames.req document for.

PCK 

definitions

CK 
TK 
DYN 
BPC 

Isis specific code for binary pck.

NOTJ2000PCK 

PCK frame not referenced to J2000.

◆ PartialType

This enumeration indicates whether the partial derivative is taken with respect to Right Ascension, Declination, or Twist (or Rotation).

Enumerator
WRT_RightAscension 

With respect to Right Ascension.

WRT_Declination 

With respect to Declination.

WRT_Twist 

With respect to Twist or Prime Meridian Rotation.

◆ Source

The rotation can come from one of 3 places for an Isis cube.

The class expects function to be after Memcache. Spice - the rotation is calculated by Naif Spice routines with data read directly from Naif kernels. Nadir - the rotation is calculated using the Naif routine twovec with the position and velocity vectors of the spacecraft. Memcache - the rotation is linearly interpolated from time-based values in a table. PolyFunction - the rotation is calculated from an nth degree polynomial in one variable (time in scaled seconds) PolyFunctionOverSpice - the rotation is calculated from an nth degree polynomial fit over the Naif Spice results. PckPolyFunction - The rotation is calculated using the IAU fit polynomials in one variable (time in Julian centuries and days).

Enumerator
Spice 

Directly from the kernels.

Nadir 

Nadir pointing.

Memcache 

From cached table.

PolyFunction 

From nth degree polynomial.

PolyFunctionOverSpice 

Kernels plus nth degree polynomial.

PckPolyFunction 

Quadratic polynomial function with linear trignometric terms.

Constructor & Destructor Documentation

◆ SpiceRotation() [1/3]

Isis::SpiceRotation::SpiceRotation ( int frameCode)

Construct an empty SpiceRotation class using a valid Naif frame code to set up for getting rotation from Spice kernels.

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

Parameters
frameCodeValid naif frame code.

References m_orientation, No, p_axis1, p_axis2, p_axis3, p_degree, Spice, and UNKNOWN.

◆ SpiceRotation() [2/3]

Isis::SpiceRotation::SpiceRotation ( int frameCode,
int targetCode )

Construct an empty SpiceRotation object using valid Naif frame code and.

body code to set up for computing nadir rotation. See required reading ftp://naif.jpl.nasa.gov/pub/naif/toolkit_docs/C/req/naif_ids.html

Parameters
frameCodeValid naif frame code.
targetCodeValid naif body code.
Exceptions
IException::Io"Cannot find [key] in text kernels"

References _FILEINFO_, Isis::NaifStatus::CheckErrors(), DYN, Isis::IException::Io, m_orientation, Nadir, No, p_axis1, p_axis2, p_axis3, p_degree, and Isis::toString().

◆ SpiceRotation() [3/3]

Isis::SpiceRotation::SpiceRotation ( const SpiceRotation & rotToCopy)

Construct a SpiceRotation object by copying from an existing one.

Parameters
rotToCopyconst reference to other SpiceRotation to copy

References m_orientation, p_axis1, p_axis2, p_axis3, p_cacheTime, and p_degree.

◆ ~SpiceRotation()

Isis::SpiceRotation::~SpiceRotation ( )
virtual

Destructor for SpiceRotation object.

References m_orientation.

Member Function Documentation

◆ Angles()

std::vector< double > Isis::SpiceRotation::Angles ( int axis3,
int axis2,
int axis1 )

Return the camera angles (right ascension, declination, and twist) for the time-based matrix CJ.

Parameters
axis3The rotation axis for the third angle
axis2The rotation axis for the second angle
axis1The rotation axis for the first angle
Returns
vector<double> Camera angles (ra, dec, twist)

References Isis::NaifStatus::CheckErrors().

Referenced by DCJdt(), GetCenterAngles(), setEphemerisTimePolyFunctionOverSpice(), SetPolynomial(), toJ2000Partial(), and ToReferencePartial().

◆ AngularVelocity()

std::vector< double > Isis::SpiceRotation::AngularVelocity ( )

Accessor method to get the angular velocity.

Returns
vector<double> Angular velocity

◆ Cache()

Table Isis::SpiceRotation::Cache ( const QString & tableName)

Return a table with J2000 to reference rotations.

Return a table containing the cached pointing with the given name. The table will have either five columns (for a list cache) of J2000 to reference quaternions and times, eight columns (if angular velocity is available), or three columns (for a coefficient cache), of J2000 to reference frame rotation angles defined by coefficients of a polynomial function (see SetPolynommial). In the coefficient cache the last row of the table is the base time, time scale, and polynomial degree. Note: In the case of the coefficient cache, the angular velocity is not written to the table since it can be calculated from the polynomials.

Parameters
tableNameName of the table to create and return
Exceptions
IException::Programmer"To create table source of data must be either Memcache or PolyFunction"
Returns
Table Table with given name that contains the cached pointing

References _FILEINFO_, CacheLabel(), Isis::TableField::Double, LineCache(), LoadTimeCache(), m_orientation, Memcache, p_cacheTime, p_degree, PolyFunction, PolyFunctionOverSpice, Isis::IException::Programmer, and Yes.

Referenced by LineCache().

◆ CacheLabel()

void Isis::SpiceRotation::CacheLabel ( Table & table)

Add labels to a SpiceRotation table.

Return a table containing the labels defining the rotation.

Parameters
TableTable to receive labels

References Isis::NaifStatus::CheckErrors(), Isis::Table::Label(), PCK, and Isis::toString().

Referenced by Cache().

◆ cacheSize()

int Isis::SpiceRotation::cacheSize ( )
inline

References m_orientation.

Referenced by LoadCache(), and SetFullCacheParameters().

◆ checkForBinaryPck()

void Isis::SpiceRotation::checkForBinaryPck ( )

Check loaded pck to see if any are binary and set frame type to indicate binary pck.

This is strictly a local method to be called only when the source is Spice. Its purpose is to check loaded PCK to see if any are binary. If any loaded PCK is binary, the body rotation label keyword FrameTypeCode is set to BPC (6).

References BPC.

Referenced by loadPCFromSpice().

◆ ComputeAv()

void Isis::SpiceRotation::ComputeAv ( )

Compute the angular velocity from the time-based functions fit to the pointing angles This method computes omega = angular velocity matrix, and extracts the angular velocity.

See comments in the Naif Spicelib routine xf2rav_c.c.

       _                     _
      |                       |
      |   0    -av[2]  av[1]  |
      |                       |

omega = | av[2] 0 -av[0] | | | | -av[1] av[0] 0 | |_ _|

Exceptions
IException::Programmer"The SpiceRotation pointing angles must be fit to polynomials in order to compute angular velocity."
IException::Programmer"Planetary angular velocity must be fit computed with PCK polynomials"

References _FILEINFO_, BPC, Isis::NaifStatus::CheckErrors(), CK, DCJdt(), DYN, INERTL, NOTJ2000PCK, PCK, PckPolyFunction, PolyFunction, Isis::IException::Programmer, TK, and UNKNOWN.

Referenced by setEphemerisTimePolyFunction().

◆ ComputeBaseTime()

void Isis::SpiceRotation::ComputeBaseTime ( )

Compute the base time using cached times.

References p_cacheTime.

Referenced by SetPolynomial(), and SetPolynomial().

◆ ConstantFrameChain()

std::vector< int > Isis::SpiceRotation::ConstantFrameChain ( )

Accessor method to get the frame chain for the constant part of the rotation (ends in target)

Returns
vector<int> The frame chain for the constant part of the rotation.

◆ ConstantMatrix()

std::vector< double > & Isis::SpiceRotation::ConstantMatrix ( )

Return the constant 3x3 rotation TC matrix as a vector of length 9.

Returns
vector<double> Constant rotation matrix, TC.

◆ ConstantRotation()

std::vector< double > Isis::SpiceRotation::ConstantRotation ( )

Return the constant 3x3 rotation TC matrix as a quaternion.

Returns
vector<double> Constant rotation quaternion, TC.

References Isis::NaifStatus::CheckErrors().

◆ DCJdt()

void Isis::SpiceRotation::DCJdt ( std::vector< double > & dCJ)

Compute the derivative of the 3x3 rotation matrix CJ with respect to time.

The derivative is computed based on p_CJ (J2000 to first constant frame).

Parameters
[out]dCJDerivative of p_CJ

References Angles(), Isis::NaifStatus::CheckErrors(), Isis::PolynomialUnivariate::DerivativeVar(), p_axis1, p_axis2, p_axis3, p_degree, and Isis::BasisFunction::SetCoefficients().

Referenced by ComputeAv().

◆ DPckPolynomial()

double Isis::SpiceRotation::DPckPolynomial ( SpiceRotation::PartialType partialVar,
const int coeffIndex )

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

Parameters
partialVarVariable derivative is to be with respect to
coeffIndexThe index of the coefficient to differentiate
Exceptions
IException::Programmer"Unable to evaluate the derivative of the SPICE rotation fit polynomial for the given coefficient index. Index is negative or exceeds degree of polynomial"
Returns
double The derivative evaluated at the current time.

References _FILEINFO_, p_degree, Isis::IException::Programmer, Isis::toString(), WRT_Declination, WRT_RightAscension, and WRT_Twist.

Referenced by toJ2000Partial(), and ToReferencePartial().

◆ DPolynomial()

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

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

Parameters
coeffIndexThe index of the coefficient to differentiate
Exceptions
IException::Programmer"Unable to evaluate the derivative of the SPICE rotation fit polynomial for the given coefficient index. Index is negative or exceeds degree of polynomial"
Returns
double The derivative evaluated at the current time.

References _FILEINFO_, p_degree, Isis::IException::Programmer, and Isis::toString().

Referenced by toJ2000Partial(), and ToReferencePartial().

◆ EphemerisTime()

double Isis::SpiceRotation::EphemerisTime ( ) const

Accessor method to get current ephemeris time.

Returns
double The current ephemeris time.

Referenced by Isis::Spice::computeSolarLongitude(), and Isis::IsisBody::rotation().

◆ EvaluatePolyFunction()

std::vector< double > Isis::SpiceRotation::EvaluatePolyFunction ( )

Evaluate the polynomial fit function for the three pointing angles for the current ephemeris time.

Returns
vector<double> Vector containing the three rotation angles.

References p_degree.

Referenced by setEphemerisTimePolyFunctionOverSpice().

◆ Extrapolate()

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

Extrapolate pointing for a given time assuming a constant angular velocity.

The pointing and angular velocity at the current time will be used to extrapolate pointing at the input time. If angular velocity does not exist, the value at the current time will be output.

Parameters
[in]timeEtThe time of the pointing to be extrapolated
Returns
vector<double> A quaternion defining the rotation at the input time.

References Isis::NaifStatus::CheckErrors().

◆ Frame()

int Isis::SpiceRotation::Frame ( )

Accessor method that returns the frame code.

This is the first value of the constant frames member variable.

Returns
int An integer value indicating the frame code.

◆ FrameTrace()

void Isis::SpiceRotation::FrameTrace ( double et)

Compute frame trace chain from target frame to J2000.

Parameters
etEphemeris time
Exceptions
IException::Programmer"The frame is not supported by Naif"
IException::Programmer"The ck rotation from frame can not be found due to no pointing available at request time or a problem with the frame"
IException::Programmer"The frame has a type not supported by your version of Naif Spicelib. You need to update."

References _FILEINFO_, Isis::NaifStatus::CheckErrors(), CK, DYN, INERTL, J2000Code, Nadir, PCK, Isis::IException::Programmer, TK, and Isis::toString().

Referenced by InitConstantRotation().

◆ GetBaseTime()

double Isis::SpiceRotation::GetBaseTime ( )

Accessor method to get the rotation base time.

Returns
double The base time for the rotation.

◆ GetCenterAngles()

std::vector< double > Isis::SpiceRotation::GetCenterAngles ( )

Return the camera angles at the center time of the observation.

Returns
vector<double> Camera angles at center time

References Angles(), p_axis1, p_axis2, p_axis3, and SetEphemerisTime().

Referenced by Isis::IsisBundleObservation::bundleOutputFetchData().

◆ getFrameType()

SpiceRotation::FrameType Isis::SpiceRotation::getFrameType ( )

Accessor method to get the rotation frame type.

Returns
SpiceRotation::FrameType The frame type of the rotation.

◆ GetFullCacheTime()

std::vector< double > Isis::SpiceRotation::GetFullCacheTime ( )

Return full listing (cache) of original time coverage requested.

Exceptions
IException::User"Time cache not availabe -- rerun spiceinit"
Returns
vector<double> Cache of original time coverage.

References _FILEINFO_, and Isis::IException::User.

◆ getPckPolynomial()

void Isis::SpiceRotation::getPckPolynomial ( std::vector< Angle > & raCoeff,
std::vector< Angle > & decCoeff,
std::vector< Angle > & pmCoeff )

Return the coefficients of a polynomial fit to each of the three planet angles.

See SetPckPolynomial() for more information on the polynomial.

TODO??? Should the other coefficients be returned as well???

Parameters
[out]raCoeffQuadratic coefficients of fit to ra
[out]decCoeffQuadratic coefficients of fit to dec
[out]pmcoeffQuadratic coefficients of fit to pm

◆ GetPolynomial()

void Isis::SpiceRotation::GetPolynomial ( std::vector< double > & coeffAng1,
std::vector< double > & coeffAng2,
std::vector< double > & coeffAng3 )

Return the coefficients of a polynomial fit to each of the three camera angles 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 and n = p_degree.

Parameters
[out]coeffAng1Coefficients of fit to Angle 1
[out]coeffAng2Coefficients of fit to Angle 2
[out]coeffAng3Coefficients of fit to Angle 3

Referenced by Isis::IsisBundleObservation::applyParameterCorrections(), and Isis::IsisBundleObservation::bundleOutputFetchData().

◆ GetSource()

SpiceRotation::Source Isis::SpiceRotation::GetSource ( )

Accessor method to get the rotation source.

Returns
SpiceRotation::Source The source of the rotation.

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

◆ GetTimeScale()

double Isis::SpiceRotation::GetTimeScale ( )

Accessor method to get the rotation time scale.

Returns
double The time scale for the rotation.

◆ HasAngularVelocity()

bool Isis::SpiceRotation::HasAngularVelocity ( )

Checks whether the rotation has angular velocities.

Returns
bool Indicates whether the rotation has angular velocities.

◆ InitConstantRotation()

void Isis::SpiceRotation::InitConstantRotation ( double et)

Initialize the constant rotation.

Parameters
etEphemeris time.

References FrameTrace(), and refchg_().

Referenced by LoadCache(), setEphemerisTimeNadir(), and setEphemerisTimeSpice().

◆ IsCached()

bool Isis::SpiceRotation::IsCached ( ) const

Checks if the cache is empty.

Returns
bool Indicates whether this rotation is cached.

References m_orientation.

Referenced by Isis::Spice::computeSolarLongitude(), and Isis::Spice::createCache().

◆ J2000Vector()

std::vector< double > Isis::SpiceRotation::J2000Vector ( const std::vector< double > & rVec)

Given a direction vector in the reference frame, return a J2000 direction.

Parameters
[in]rVecA direction vector in the reference frame
Returns
vector<double> A direction vector in J2000 frame.

References Isis::NaifStatus::CheckErrors(), and invstm_().

◆ LineCache()

Table Isis::SpiceRotation::LineCache ( const QString & tableName)

Return a table with J2000 to reference rotations.

Return a table containing the cached pointing with the given name. The table will have eight columns, quaternion, angular velocity, and time of J2000 to reference frame rotation.

Parameters
tableNameName of the table to create and return
Exceptions
IException::Programmer"Only cached rotations can be returned as a line cache of quaternions and time"
Returns
Table Table with given name that contains the cached pointing

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

Referenced by Cache().

◆ LoadCache() [1/4]

void Isis::SpiceRotation::LoadCache ( double startTime,
double endTime,
int size )

Cache J2000 rotation quaternion over a time range.

This method will load an internal cache with frames over a time range. This prevents the NAIF kernels from being read over-and-over again and slowing an 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 frames to keep in the cache
Exceptions
IException::Programmer"Argument cacheSize must not be less than or equal to zero"
IException::Programmer"Argument startTime must be less than or equal to endTime"
IException::Programmer"Cache size must be more than 1 if startTime and endTime differ"
IException::Programmer"A SpiceRotation cache has already men

References _FILEINFO_, cacheSize(), InitConstantRotation(), LoadTimeCache(), m_orientation, Memcache, p_cacheTime, Isis::IException::Programmer, SetEphemerisTime(), and Yes.

Referenced by Isis::Spice::createCache(), and LoadCache().

◆ LoadCache() [2/4]

void Isis::SpiceRotation::LoadCache ( double time)

Cache J2000 to frame rotation for a time.

This method will load an internal cache with a rotation 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

References LoadCache().

◆ LoadCache() [3/4]

void Isis::SpiceRotation::LoadCache ( nlohmann::json & isd)

Load the cached data from an ALE ISD.

The SpiceRotation object must be set to a SPICE source before loading the cache.

Parameters
isdRotThe ALE ISD as a JSON object.

References _FILEINFO_, CK, m_orientation, Memcache, p_cacheTime, Isis::IException::Programmer, and SetEphemerisTime().

◆ LoadCache() [4/4]

void Isis::SpiceRotation::LoadCache ( Table & table)

Cache J2000 rotations using a table file.

This method will load either an internal cache with rotations (quaternions) or coefficients (for 3 polynomials defining the camera angles) from an ISIS table file. In the first case, the table must have 5 columns and at least one row. The 5 columns contain the following information, J2000 to reference quaternion (q0, q1, q2, q3) and the ephemeris time of that position. If there are multiple rows, it is assumed the quaternions between the rows can be interpolated. In the second case, the table must have three columns and at least two rows. The three columns contain the coefficients for each of the three camera angles. Row one of the table contains coefficient 0 (constant term) for angles 1, 2, and 3. If the degree of the fit equation is greater than 1, row 2 contains coefficient 1 (linear) for each of the three angles. Row n contains coefficient n-1 and the last row contains the time parameters, base time, and time scale, and the degree of the polynomial.

Parameters
tableAn ISIS table blob containing valid J2000 to reference quaternion/time values
Exceptions
IException::Programmer"Expecting either three, five, or eight fields in the SpiceRotation table"

References _FILEINFO_, Isis::PvlObject::findKeyword(), Isis::PvlObject::hasKeyword(), J2000Code, Isis::Table::Label(), loadPCFromTable(), m_orientation, Memcache, p_cacheTime, PCK, PolyFunction, Isis::IException::Programmer, Isis::Table::Records(), SetOverrideBaseTime(), SetPolynomial(), SetPolynomialDegree(), Isis::toDouble(), Isis::toInt(), and UNKNOWN.

◆ loadPCFromSpice()

void Isis::SpiceRotation::loadPCFromSpice ( int centerBody)

Initialize planetary orientation constants from Spice PCK.

Retrieve planetary orientation constants from a Spice PCK and store them in the class.

Parameters
centerBodyNAIF id for the planetary body to retrieve the PCK for

References Isis::NaifStatus::CheckErrors(), checkForBinaryPck(), Isis::Angle::Degrees, NOTJ2000PCK, PCK, and Isis::toString().

◆ loadPCFromTable()

void Isis::SpiceRotation::loadPCFromTable ( const PvlObject & label)

Initialize planetary orientation constants from an cube body rotation label.

Retrieve planetary orientation constants from a cube body rotation label if they are present.

Parameters
labelconst reference to the cube body rotation pvl label

References Isis::NaifStatus::CheckErrors(), Isis::Angle::Degrees, Isis::PvlObject::hasKeyword(), and Isis::toDouble().

Referenced by LoadCache().

◆ LoadTimeCache()

void Isis::SpiceRotation::LoadTimeCache ( )

Load the time cache.

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

Exceptions
IException::Programmer"Full cache size does NOT match cache size in LoadTimeCache -- should never happen"
IException::Programmer"Observation crosses segment boundary--unable to interpolate pointing"
IException::User"No camera kernels loaded...Unable to determine time cache to downsize"

References _FILEINFO_, Isis::NaifStatus::CheckErrors(), ck3sdn(), Done, m_orientation, Memcache, Nadir, No, p_cacheTime, Isis::IException::Programmer, Isis::IException::User, and Yes.

Referenced by Cache(), LoadCache(), and ReloadCache().

◆ Matrix()

std::vector< double > Isis::SpiceRotation::Matrix ( )

Return the full rotation TJ as a matrix.

Returns
vector<double> Returned matrix.

References Isis::NaifStatus::CheckErrors().

Referenced by Isis::Spice::computeSolarLongitude(), Isis::IsisBody::rotation(), and Isis::Spice::sunToBodyDist().

◆ MinimizeCache()

void Isis::SpiceRotation::MinimizeCache ( DownsizeStatus status)

Set the downsize status to minimize cache.

Parameters
statusThe DownsizeStatus enumeration value.

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

◆ pmCoefs()

std::vector< Angle > Isis::SpiceRotation::pmCoefs ( )

Return the coefficients used to calculate the target body prime meridian.

Return the coefficients used to calculate the target body prime meridian without nutation/precession. The model is a standard quadratic polynomial in time in days from the standard epoch (usually J2000). To match the Naif PCK prime meridian, (usually the same as the most recent IAU Report) the trignometric terms to account for nutation/precession need to be added.

pm = pm0 + pm1*d + pm2*d**2 + sum(pmcoef[i]i*sin(angle[i]))

Returns
vector<double> A vector of length 3 containing the prime meridian coefficients.

◆ pmNutPrecCoefs()

std::vector< double > Isis::SpiceRotation::pmNutPrecCoefs ( )

Return the coefficients used to calculate the target body pm nut/prec coefficients.

Return the coefficients used to calculate the target body prime meridian nutation/precession contribution. The model is the sum of the products of the coefficients returned by this method with the sin of the corresponding angles associated with the barycenter related to the target body.

prime meridian = pm0 + pm1*T + pm2*d**2 + sum(pmCoef[i]i*sin(angle[i]))

Returns
vector<double> A vector containing the pm nut/prec coeffcients.

◆ poleDecCoefs()

std::vector< Angle > Isis::SpiceRotation::poleDecCoefs ( )

Return the coefficients used to calculate the target body pole dec.

Return the coefficients used to calculate the target body declination without nutation/precession. The model is a standard quadratic polynomial in time in Julian centuries (36525 days) from the standard epoch (usually J2000). To match the Naif PCK pole declination (usually the same as the most recent IAU Report) the trignometric terms to account for nutation/precession need to be added.

pole dec = dec0 + dec1*T + dec2*T**2 + sum(racoef[i]i*sin(angle[i]))

Returns
vector<double> A vector of length 3 containing the pole declination coefficients.

◆ poleDecNutPrecCoefs()

std::vector< double > Isis::SpiceRotation::poleDecNutPrecCoefs ( )

Return the coefficients used to calculate the target body pole dec nut/prec coefficients.

Return the coefficients used to calculate the target body declination nutation/precession contribution. The model is the sum of the products of the coefficients returned by this method with the sin of the corresponding angles associated with the barycenter related to the target body.

pole dec = dec0 + dec1*T + dec2*T**2 + sum(decCoef[i]i*sin(angle[i]))

Returns
vector<double> A vector containing the pole dec nut/prec coeffcients.

◆ poleRaCoefs()

std::vector< Angle > Isis::SpiceRotation::poleRaCoefs ( )

Return the coefficients used to calculate the target body pole ra.

Return the coefficients used to calculate the target body right ascension without nutation/precession. The model is a standard quadratic polynomial in time in Julian centuries (36525 days) from the standard epoch (usually J2000). To match the Naif PCK pole right ascension (usually the same as the most recent IAU Report) the trignometric terms to account for nutation/precession need to be added.

pole ra = ra0 + ra1*T + ra2*T**2 + sum(racoef[i]i*sin(angle[i]))

Returns
vector<double> A vector of length 3 containing the pole ra coefficients

◆ poleRaNutPrecCoefs()

std::vector< double > Isis::SpiceRotation::poleRaNutPrecCoefs ( )

Return the coefficients used to calculate the target body pole ra nut/prec coefficients.

Return the coefficients used to calculate the target body right ascension nutation/precession contribution. The model is the sum of the products of the coefficients returned by this method with the sin of the corresponding angles associated with the barycenter related to the target body.

pole ra = ra0 + ra1*T + ra2*T**2 + sum(raCoef[i]i*sin(angle[i]))

Returns
vector<double> A vector containing the pole ra nut/prec coefficients.

◆ ReferenceVector()

std::vector< double > Isis::SpiceRotation::ReferenceVector ( const std::vector< double > & jVec)

Given a direction vector in J2000, return a reference frame direction.

Parameters
[in]jVecA direction vector in J2000
Returns
vector<double> A direction vector in reference frame.

References Isis::NaifStatus::CheckErrors().

Referenced by Isis::IsisBody::fixedVector(), Isis::Spice::instrumentBodyFixedPosition(), Isis::Spice::instrumentBodyFixedVelocity(), Isis::Spice::setTime(), Isis::Spice::subSpacecraftPoint(), and Isis::Spice::targetCenterDistance().

◆ ReloadCache()

void Isis::SpiceRotation::ReloadCache ( )

Cache J2000 rotation over existing cached time range using polynomials.

This method will reload an internal cache with matrices formed from rotation angles fit to functions over a time range.

Exceptions
IException::Programmer"The SpiceRotation has not yet been fit to a function"

References _FILEINFO_, LoadTimeCache(), m_orientation, Memcache, No, p_cacheTime, PolyFunction, PolyFunctionOverSpice, Isis::IException::Programmer, and SetEphemerisTime().

Referenced by LineCache().

◆ SetAngles()

void Isis::SpiceRotation::SetAngles ( std::vector< double > angles,
int axis3,
int axis2,
int axis1 )

Set the rotation angles (phi, delta, and w) for the current time to define the time-based matrix CJ.

This method was created for unitTests and should not be used otherwise. It only works for cached data with a cache size of 1.

Parameters
[in]anglesThe angles defining the rotation (phi, delta, and w) in radians
[in]axis3The rotation axis for the third angle
[in]axis2The rotation axis for the second angle
[in]axis1The rotation axis for the first angle

References m_orientation, p_cacheTime, and SetEphemerisTime().

◆ SetAxes()

void Isis::SpiceRotation::SetAxes ( int axis1,
int axis2,
int axis3 )

Set the axes of rotation for decomposition of a rotation matrix into 3 angles.

Parameters
[in]axis1Axes of rotation of first angle applied (right rotation)
[in]axis2Axes of rotation of second angle applied (center rotation)
[in]axis3Axes of rotation of third angle applied (left rotation)
Exceptions
IException::Programmer"A rotation axis is outside the valid range of 1 to 3"
Returns
double Wrapped angle.

References _FILEINFO_, p_axis1, p_axis2, p_axis3, and Isis::IException::Programmer.

◆ SetCacheTime()

void Isis::SpiceRotation::SetCacheTime ( std::vector< double > cacheTime)

References p_cacheTime.

◆ SetConstantMatrix()

void Isis::SpiceRotation::SetConstantMatrix ( std::vector< double > constantMatrix)

Set the constant 3x3 rotation TC matrix from a vector of length 9.

Parameters
constantMatrixConstant rotation matrix, TC.

◆ SetEphemerisTime()

void Isis::SpiceRotation::SetEphemerisTime ( double et)

Return the J2000 to reference frame quaternion at given time.

This method returns the J2000 to reference frame rotational matrix at a given et in seconds. The quaternion is obtained from either valid NAIF ck and/or fk, or alternatively from an internal cache loaded from an ISIS Table object. In the first case, the kernels must contain the rotation for the frame specified in the constructor at the given time (as well as all the intermediate frames going from the reference frame to J2000) and they must be loaded using the SpiceKernel class.

Parameters
etephemeris time in seconds

References Memcache, PckPolyFunction, PolyFunction, PolyFunctionOverSpice, setEphemerisTimeMemcache(), setEphemerisTimeNadir(), setEphemerisTimePckPolyFunction(), setEphemerisTimePolyFunction(), setEphemerisTimePolyFunctionOverSpice(), and setEphemerisTimeSpice().

Referenced by Isis::Spice::computeSolarLongitude(), GetCenterAngles(), LoadCache(), LoadCache(), ReloadCache(), Isis::IsisBody::rotation(), SetAngles(), SetPolynomial(), SetPolynomial(), and Isis::Spice::setTime().

◆ setEphemerisTimeMemcache()

void Isis::SpiceRotation::setEphemerisTimeMemcache ( )
protected

Updates rotation state based on the rotation cache.

When setting the ephemeris time, this method is used to update the rotation state by reading from the rotation cache

See also
SpiceRotation::SetEphemerisTime

References Isis::NaifStatus::CheckErrors(), m_orientation, and p_cacheTime.

Referenced by SetEphemerisTime(), and setEphemerisTimePolyFunctionOverSpice().

◆ setEphemerisTimeNadir()

void Isis::SpiceRotation::setEphemerisTimeNadir ( )
protected

When setting the ephemeris time, uses spacecraft nadir source to update the rotation state.

See also
SpiceRotation::SetEphemerisTime

References Isis::NaifStatus::CheckErrors(), and InitConstantRotation().

Referenced by SetEphemerisTime().

◆ setEphemerisTimePckPolyFunction()

void Isis::SpiceRotation::setEphemerisTimePckPolyFunction ( )
protected

When setting the ephemeris time, updates the rotation state based on the PcK polynomial.

See also
SpiceRotation::SetEphemerisTime

References Isis::DEG2RAD, Isis::Angle::Degrees, Isis::HALFPI, p_axis1, p_axis2, p_axis3, and Isis::Angle::Radians.

Referenced by SetEphemerisTime(), and setPckPolynomial().

◆ setEphemerisTimePolyFunction()

void Isis::SpiceRotation::setEphemerisTimePolyFunction ( )
protected

When setting the ephemeris time, updates the rotation according to a polynomial that defines the three camera angles and angular velocity, if available.

See also
SpiceRotatation::SetEphemerisTime

References Isis::NaifStatus::CheckErrors(), ComputeAv(), m_orientation, p_axis1, p_axis2, p_axis3, and p_degree.

Referenced by SetEphemerisTime(), and setEphemerisTimePolyFunctionOverSpice().

◆ setEphemerisTimePolyFunctionOverSpice()

void Isis::SpiceRotation::setEphemerisTimePolyFunctionOverSpice ( )
protected

When setting the ephemeris time, updates the rotation state based on a polynomial fit over spice kernel data.

See also
SpiceRotation::SetEphemerisTime

References Angles(), Isis::NaifStatus::CheckErrors(), EvaluatePolyFunction(), p_axis1, p_axis2, p_axis3, setEphemerisTimeMemcache(), and setEphemerisTimePolyFunction().

Referenced by SetEphemerisTime().

◆ setEphemerisTimeSpice()

void Isis::SpiceRotation::setEphemerisTimeSpice ( )
protected

When setting the ephemeris time, updates the rotation state based on data read directly from NAIF kernels using NAIF Spice routines.

Exceptions
IException::Io"[framecode] is an unrecognized reference frame code. Has the mission frames kernel been loaded?"
IException::Io"No pointing is availabe at request time for frame code"
See also
SpiceRotation::SetEphemerisTime

References _FILEINFO_, Isis::NaifStatus::CheckErrors(), frmchg_(), InitConstantRotation(), Isis::IException::Io, J2000Code, and refchg_().

Referenced by SetEphemerisTime().

◆ SetFrame()

void Isis::SpiceRotation::SetFrame ( int frameCode)

Change the frame to the given frame code.

This method has no effect if spice is cached.

Parameters
frameCodeThe integer-valued frame code

◆ SetFullCacheParameters()

void Isis::SpiceRotation::SetFullCacheParameters ( double startTime,
double endTime,
int cacheSize )
protected

Set the full cache time parameters.

Parameters
[in]startTimeThe earliest time of the full cache coverage
[in]endTimeThe latest time of the full cache coverage
[in]cacheSizeThe number of epochs in the full (line) cache

References cacheSize().

◆ SetOverrideBaseTime()

void Isis::SpiceRotation::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 save base time and polynomials for the angles.

Parameters
[in]baseTimeThe base time to use and override the computed base time
[in]timeScaleThe time scale to use and override the computed time scale

Referenced by LoadCache().

◆ setPckPolynomial()

void Isis::SpiceRotation::setPckPolynomial ( const std::vector< Angle > & raCoeff,
const std::vector< Angle > & decCoeff,
const std::vector< Angle > & pmCoeff )

Set the coefficients of a polynomial fit to each of the three planet angles for the time period covered by the cache, ra = ra0 + ra1*t + ra2*t**2 + raTrig, dec = dec0 + dec1*t + dec2*t**2 + decTrig, pm = pm0 + pm1*d + pm2*t**2 + pmTrig, where t = time / (seconds per day).

for time = et d = t / (days per Julian century), and the trig terms are of the form sum(i = 0, numTerms) tcoef[i] * sin(constant[i] + lcoef[i]*t) for ra and pm and sum(i = 0, numTerms) tcoef[i] * cos(constant[i] + lcoef[i]*t) for dec.

Parameters
[in]raCoeffpole ra coefficients to set
[in]decCoeffpole dec coefficients to set
[in]pmCoeffpole pm coefficients to set

References setEphemerisTimePckPolyFunction(), and usePckPolynomial().

◆ SetPolynomial() [1/2]

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

Set the coefficients of a polynomial fit to each of the three camera angles for the time period covered by the cache, angle = a + bt + ct**2, where t = (time - p_baseTime)/ p_timeScale.

Parameters
typeRotation source type

References Angles(), Isis::NaifStatus::CheckErrors(), ComputeBaseTime(), m_orientation, p_axis1, p_axis2, p_axis3, p_cacheTime, p_degree, PolyFunction, PolyFunctionOverSpice, SetEphemerisTime(), SetPolynomial(), and WrapAngle().

Referenced by Isis::IsisBundleObservation::applyParameterCorrections(), LoadCache(), SetPolynomial(), and SetPolynomialDegree().

◆ SetPolynomial() [2/2]

void Isis::SpiceRotation::SetPolynomial ( const std::vector< double > & coeffAng1,
const std::vector< double > & coeffAng2,
const std::vector< double > & coeffAng3,
const Source type = PolyFunction )

Set the coefficients of a polynomial fit to each of the three camera angles 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, and n = p_degree.

Parameters
[in]coeffAng1Coefficients of fit to Angle 1
[in]coeffAng2Coefficients of fit to Angle 2
[in]coeffAng3Coefficients of fit to Angle 3
[in]typeRotation source type

References Isis::NaifStatus::CheckErrors(), ComputeBaseTime(), p_degree, and SetEphemerisTime().

◆ SetPolynomialDegree()

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

Set the degree of the polynomials to be fit to the three camera angles 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, and n = p_degree.

Parameters
[in]degreeDegree of the polynomial to be fit

References p_degree, and SetPolynomial().

Referenced by LoadCache().

◆ SetSource()

void Isis::SpiceRotation::SetSource ( Source source)

Resets the source of the rotation to the given value.

Parameters
sourceThe rotation source to be set.

◆ SetTimeBasedMatrix()

void Isis::SpiceRotation::SetTimeBasedMatrix ( std::vector< double > timeBasedMatrix)

Set the time-based 3x3 rotation CJ matrix from a vector of length 9.

Parameters
timeBasedMatrixTime-based rotation matrix, TC.

◆ SetTimeBias()

void Isis::SpiceRotation::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

◆ sysNutPrecCoefs()

std::vector< Angle > Isis::SpiceRotation::sysNutPrecCoefs ( )

Return the coefficients used to calculate the target body system nut/prec angles.

Return the linear coefficients used to calculate the target body system (barycenter) nutation/precession angles (periods). The model for the angles is linear in time in Julian centuries since the standard epoch (usually J2000). angles associated with the barycenter related to the target body.

angle[i] = constant[i] + coef[i]*T

Returns
vector<Angle> A vector containing the system nut/prec linear coefficients.

◆ sysNutPrecConstants()

std::vector< Angle > Isis::SpiceRotation::sysNutPrecConstants ( )

Return the constants used to calculate the target body system nut/prec angles.

Return the constant terms used to calculate the target body system (barycenter) nutation/precession angles (periods). The model for the angles is linear in time in Julian centuries since the standard epoch (usually J2000). angles associated with the barycenter related to the target body.

angle[i] = constant[i] + coef[i]*T

Returns
vector<Angle> A vector containing the system nut/prec constant terms.

◆ TimeBasedMatrix()

std::vector< double > & Isis::SpiceRotation::TimeBasedMatrix ( )

Return time-based 3x3 rotation CJ matrix as a vector of length 9.

Returns
vector<double> Time-based rotation matrix, CJ.

◆ TimeBasedRotation()

std::vector< double > Isis::SpiceRotation::TimeBasedRotation ( )

Return time-based 3x3 rotation CJ matrix as a quaternion.

Returns
vector<double> Time-based rotation quaternion, CJ.

References Isis::NaifStatus::CheckErrors().

◆ TimeFrameChain()

std::vector< int > Isis::SpiceRotation::TimeFrameChain ( )

Accessor method to get the frame chain for the rotation (begins in J2000).

Returns
vector<int> The frame chain for the rotation.

◆ toJ2000Partial()

std::vector< double > Isis::SpiceRotation::toJ2000Partial ( const std::vector< double > & lookT,
SpiceRotation::PartialType partialVar,
int coeffIndex )

Given a direction vector in the reference frame, compute the derivative with respect to one of the coefficients in the angle polynomial fit equation of a vector rotated from the reference frame to J2000.

TODO - merge this method with ToReferencePartial

Parameters
[in]lookTA direction vector in the targeted reference frame
[in]partialVarVariable derivative is to be with respect to
[in]coeffIndexCoefficient index in the polynomial fit to the variable (angle)
Exceptions
IException::User"Body rotation uses a binary PCK. Solutions for this model are not supported"
IException::User"Body rotation uses a PCK not referenced to J2000. Solutions for this model are not supported"
IException::User"Solutions are not supported for this frame type."
Returns
vector<double> A direction vector rotated by derivative of reference to J2000 rotation.

References _FILEINFO_, Angles(), BPC, Isis::NaifStatus::CheckErrors(), CK, DPckPolynomial(), DPolynomial(), DYN, NOTJ2000PCK, p_axis1, p_axis2, p_axis3, PCK, and Isis::IException::User.

◆ ToReferencePartial()

std::vector< double > Isis::SpiceRotation::ToReferencePartial ( std::vector< double > & lookJ,
SpiceRotation::PartialType partialVar,
int coeffIndex )

Compute the derivative with respect to one of the coefficients in the angle polynomial fit equation of a vector rotated from J2000 to a reference frame.

The polynomial equation is of the form angle = c0 + c1*t + c2*t**2 + ... cn*t**n, where t = (time - p_basetime) / p_timeScale and n = p_degree (the degree of the equation)

Parameters
[in]lookJLook vector in J2000 frame
[in]partialVarVariable derivative is to be with respect to
[in]coeffIndexCoefficient index in the polynomial fit to the variable (angle)
Exceptions
IException::Programmer"Only CK, DYN, and PCK partials can be calculated"
Returns
vector<double> Vector rotated by derivative of J2000 to reference rotation.

References _FILEINFO_, Angles(), Isis::NaifStatus::CheckErrors(), CK, DPckPolynomial(), DPolynomial(), DYN, p_axis1, p_axis2, p_axis3, PCK, Isis::IException::Programmer, and UNKNOWN.

◆ usePckPolynomial()

void Isis::SpiceRotation::usePckPolynomial ( )

Set the coefficients of a polynomial fit to each of the three planet angles for the time period covered by the cache, ra = ra0 + ra1*t + ra2*t**2 + raTrig, dec = dec0 + dec1*t + dec2*t**2 + decTrig, pm = pm0 + pm1*d + pm2*t**2 + pmTrig, where t = time / (seconds per day).

for time = et d = t / (days per Julian century), and the trig terms are of the form sum(i = 0, numTerms) tcoef[i] * sin(constant[i] + lcoef[i]*t) for ra and pm and sum(i = 0, numTerms) tcoef[i] * cos(constant[i] + lcoef[i]*t) for dec.

Exceptions
IException::User"Target body orientation information not available. Rerun spiceinit."

References _FILEINFO_, PckPolyFunction, and Isis::IException::User.

Referenced by setPckPolynomial().

◆ WrapAngle()

double Isis::SpiceRotation::WrapAngle ( double compareAngle,
double angle )

Wrap the input angle to keep it within 2pi radians of the angle to compare.

Parameters
[in]compareAngleLook vector in J2000 frame
[in]angleAngle to be wrapped if needed
Returns
double Wrapped angle.

References Isis::NaifStatus::CheckErrors().

Referenced by SetPolynomial().

Member Data Documentation

◆ m_orientation

◆ p_axis1

◆ p_axis2

◆ p_axis3

◆ p_cacheTime

std::vector<double> Isis::SpiceRotation::p_cacheTime
protected

◆ p_degree


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