Isis 3 Programmer Reference
SpiceRotation.h
Go to the documentation of this file.
1 #ifndef SpiceRotation_h
2 #define SpiceRotation_h
3 
25 #include <string>
26 #include <vector>
27 
28 //#include <SpiceUsr.h>
29 //#include <SpiceZfc.h>
30 //#include <SpiceZmc.h>
31 
32 #include "Angle.h"
33 #include "Table.h"
34 #include "PolynomialUnivariate.h"
35 #include "Quaternion.h"
36 
37 #define J2000Code 1
38 
39 namespace Isis {
224  public:
225  // Constructors
226  SpiceRotation(int frameCode);
227  /* SpiceRotation( int NaifCode );
228  We would like to call refchg instead to avoid the strings. Currently Naif does
229  not have refchg_c, but only the f2c'd refchg.c.*/
230  SpiceRotation(int frameCode, int targetCode);
231  SpiceRotation(const SpiceRotation &rotToCopy);
232 
233  // Destructor
234  virtual ~SpiceRotation();
235 
236  // Change the frame (has no effect if cached)
237  void SetFrame(int frameCode);
238  int Frame();
239 
240  void SetTimeBias(double timeBias);
241 
258  enum Source {
265  };
266 
271  enum PartialType {
275  };
276 
281  Yes,
283  No
284  };
285 
289  enum FrameType {
290  UNKNOWN = 0,
291  INERTL = 1,
292  PCK = 2,
293  CK = 3,
294  TK = 4,
295  DYN = 5,
296  BPC = 6,
298  };
299 
300  void SetEphemerisTime(double et);
301  double EphemerisTime() const;
302 
303  std::vector<double> GetCenterAngles();
304 
305  std::vector<double> Matrix();
306  std::vector<double> AngularVelocity();
307 
308  // TC
309  std::vector<double> ConstantRotation();
310  std::vector<double> &ConstantMatrix();
311  void SetConstantMatrix(std::vector<double> constantMatrix);
312 
313  // CJ
314  std::vector<double> TimeBasedRotation();
315  std::vector<double> &TimeBasedMatrix();
316  void SetTimeBasedMatrix(std::vector<double> timeBasedMatrix);
317 
318  std::vector<double> J2000Vector(const std::vector<double> &rVec);
319 
320  std::vector<Angle> poleRaCoefs();
321 
322  std::vector<Angle> poleDecCoefs();
323 
324  std::vector<Angle> pmCoefs();
325 
326  std::vector<double> poleRaNutPrecCoefs();
327 
328  std::vector<double> poleDecNutPrecCoefs();
329 
330  std::vector<double> pmNutPrecCoefs();
331 
332  std::vector<Angle> sysNutPrecConstants();
333 
334  std::vector<Angle> sysNutPrecCoefs();
335 
336  std::vector<double> ReferenceVector(const std::vector<double> &jVec);
337 
338  std::vector<double> EvaluatePolyFunction();
339 
340  void loadPCFromSpice(int CenterBodyCode);
341  void loadPCFromTable(const PvlObject &Label);
342 
343  void MinimizeCache(DownsizeStatus status);
344 
345  void LoadCache(double startTime, double endTime, int size);
346 
347  void LoadCache(double time);
348 
349  void LoadCache(Table &table);
350 
351  Table LineCache(const QString &tableName);
352 
353  void ReloadCache();
354 
355  Table Cache(const QString &tableName);
356  void CacheLabel(Table &table);
357 
358  void LoadTimeCache();
359 
360  std::vector<double> Angles(int axis3, int axis2, int axis1);
361  void SetAngles(std::vector<double> angles, int axis3, int axis2, int axis1);
362 
363  bool IsCached() const;
364 
365  void SetPolynomial(const Source type=PolyFunction);
366 
367  void SetPolynomial(const std::vector<double> &abcAng1,
368  const std::vector<double> &abcAng2,
369  const std::vector<double> &abcAng3,
370  const Source type = PolyFunction);
371 
372  void usePckPolynomial();
373  void setPckPolynomial(const std::vector<Angle> &raCoeff,
374  const std::vector<Angle> &decCoeff,
375  const std::vector<Angle> &pmCoeff);
376 
377  void GetPolynomial(std::vector<double> &abcAng1,
378  std::vector<double> &abcAng2,
379  std::vector<double> &abcAng3);
380 
381  void getPckPolynomial(std::vector<Angle> &raCoeff,
382  std::vector<Angle> &decCoeff,
383  std::vector<Angle> &pmCoeff);
384 
385  // Set the polynomial degree
386  void SetPolynomialDegree(int degree);
387  Source GetSource();
388  void SetSource(Source source);
389  void ComputeBaseTime();
391  double GetBaseTime();
392  double GetTimeScale();
393 
394  void SetOverrideBaseTime(double baseTime, double timeScale);
395  void SetCacheTime(std::vector<double> cacheTime);
396 
397  // Derivative methods
398  double DPolynomial(const int coeffIndex);
399  double DPckPolynomial(PartialType partialVar, const int coeffIndex);
400 
401  std::vector<double> toJ2000Partial(const std::vector<double> &lookT,
402  PartialType partialVar, int coeffIndex);
403  std::vector<double> ToReferencePartial(std::vector<double> &lookJ,
404  PartialType partialVar, int coeffIndex);
405  void DCJdt(std::vector<double> &dRJ);
406 
407  double WrapAngle(double compareAngle, double angle);
408  void SetAxes(int axis1, int axis2, int axis3);
409  std::vector<double> GetFullCacheTime();
410  void FrameTrace(double et);
411 
412  // Return the frame chain for the constant part of the rotation (ends in target)
413  std::vector<int> ConstantFrameChain();
414  std::vector<int> TimeFrameChain();
415  void InitConstantRotation(double et);
416  bool HasAngularVelocity();
417 
418  void ComputeAv();
419  std::vector<double> Extrapolate(double timeEt);
420 
421  void checkForBinaryPck();
422 
423  protected:
424  void SetFullCacheParameters(double startTime, double endTime, int cacheSize);
426  void setEphemerisTimeNadir();
427  void setEphemerisTimeSpice();
431  std::vector<double> p_cacheTime;
432  std::vector<std::vector<double> > p_cache;
437  int p_degree;
438  int p_axis1;
439  int p_axis2;
440  int p_axis3;
441 
442  private:
443  // method
444  void setFrameType();
445  std::vector<int> p_constantFrames;
448  std::vector<int> p_timeFrames;
451  double p_timeBias;
452 
453  double p_et;
457  bool p_matrixSet;
459 
460 
463  int p_axisP;
465  int p_axisV;
468 
469  double p_baseTime;
470  double p_timeScale;
474  std::vector<double> p_coefficients[3];
483  std::vector<double> p_TC;
486  std::vector<double> p_CJ;
488  std::vector<std::vector<double> > p_cacheAv;
490  std::vector<double> p_av;
493  std::vector<double> StateTJ();
495  // The remaining items are only used for PCK frame types. In this case the
496  // rotation is stored as a cache, but the coefficients are available for display
497  // or comparison, and the first three coefficient sets can be solved for and
498  // updated in jigsaw. The initial coefficient values are read from a Naif PCK.
499  //
500  // The general equation for the right ascension of the pole is
501  //
502  // raPole = raPole[0] + raPole[1]*Time + raPole[2]*Time**2 + raNutPrec,
503  // where
504  // raNutPrec = raNutPrec1[0]*sin(sysNutPrec[0][0] + sysNutPrec[0][1]*Time) +
505  // raNutPrec1[1]*sin(sysNutPrec[1][0] + sysNutPrec[1][1]*Time) + ...
506  // raNutPrec1[N-1]*sin(sysNutPrec[N-1][0] + sysNutPrec[N-1][1]*Time) +
507  // (optional for multiples of nutation precession angles)
508  // raNutPrec2[0]*sin(2*(sysNutPrec[0][0] + sysNutPrec[0][1]*Time)) +
509  // raNutPrec2[1]*sin(2*(sysNutPrec[1][0] + sysNutPrec[1][1]*Time)) + ...
510  // raNutPrec2[N-1]*sin(2*(sysNutPrec[N-1][0] + sysNutPrec[N-1][1]*Time)) +
511  // raNutPrecM[0]*sin(M*(sysNutPrec[0][0] + sysNutPrec[0][1]*Time)) +
512  // raNutPrecM[1]*sin(M*(sysNutPrec[1][0] + sysNutPrec[1][1]*Time)) + ...
513  // raNutPrecM[N-1]*sin(M*(sysNutPrec[N-1][0] + sysNutPrec[N-1][1]*Time)) +
514  //
515  // The general equation for the declination of the pole is
516  //
517  // decPole = p_decPole[0] + p_decPole[1]*Time + p_decPole[2]*Time**2 + decNutPrec,
518  // where
519  // decNutPrec = decNutPrec1[0]*cos(sysNutPrec[0][0] + sysNutPrec[0][1]*Time) +
520  // decNutPrec1[1]*cos(sysNutPrec[1][0] + sysNutPrec[1][1]*Time) + ...
521  // decNutPrec1[N-1]*cos(sysNutPrec[N-1][0] + sysNutPrec[N-1][1]*Time) +
522  // decNutPrec2[0]*cos(2*(sysNutPrec[0][0] + sysNutPrec[0][1]*Time)) +
523  // decNutPrec2[1]*cos(2*(sysNutPrec[1][0] + sysNutPrec[1][1]*Time)) + ...
524  // decNutPrec2[N-1]*cos(2*(sysNutPrec[N-1][0] + sysNutPrec[N-1][1]*Time)) +
525  // (optional for multiples of nutation precession angles)
526  // decNutPrecM[0]*sin(M*(sysNutPrec[0][0] + sysNutPrec[0][1]*Time)) +
527  // decNutPrecM[1]*sin(M*(sysNutPrec[1][0] + sysNutPrec[1][1]*Time)) + ...
528  // decNutPrecM[N-1]*sin(M*(sysNutPrec[N-1][0] + sysNutPrec[N-1][1]*Time))
529  //
530  // and Time is julian centuries since J2000.
531  //
532  // The general equation for the prime meridian rotation is
533  //
534  // pm = p_pm[0] + p_pm[1]*Dtime + p_pm[2]*Dtime**2 + pmNutPrec,
535  // where
536  // pmNutPrec = pmNutPrec1[0]*sin(sysNutPrec[0][0] + sysNutPrec[0][1]*Time) +
537  // pmNutPrec1[1]*sin(sysNutPrec[1][0] + sysNutPrec[1][1]*Time) + ...
538  // pmNutPrec1[N-1]*sin(sysNutPrec[N-1][0] + sysNutPrec[N-1][1]*Time) +
539  // (optional for multiples of nutation precession angles)
540  // pmNutPrec2[0]*sin(2*(sysNutPrec[0][0] + sysNutPrec[0][1]*Time)) +
541  // pmNutPrec2[1]*sin(2*(sysNutPrec[1][0] + sysNutPrec[1][1]*Time)) + ...
542  // pmNutPrec2[N-1]*sin(2*(sysNutPrec[N-1][0] + sysNutPrec[N-1][1]*Time)) +
543  // pmNutPrecM[0]*sin(M*(sysNutPrec[0][0] + sysNutPrec[0][1]*Time)) +
544  // pmNutPrecM[1]*sin(M*(sysNutPrec[1][0] + sysNutPrec[1][1]*Time)) + ...
545  // pmNutPrecM[N-1]*sin(M*(sysNutPrec[N-1][0] + sysNutPrec[N-1][1]*Time))
546  //
547  // Time is interval in Julian centuries since the standard epoch,
548  // dTime is interval in days from the standard epoch (J2000),
549  //
550  // N is the number of nutation/precession terms for the planetary system of the target
551  // body, (possibly including multiple angles as unique terms,
552  // ie. 2*sysNutPrec[0][0] + sysNutPrec[][1]*Time).
553  //
554  // Many of the constants in this equation are 0. for a given body.
555  //
556  // M is included as an option for future improvements. M = highest multiple (period)
557  // of any of the nutation/precession angles included in the equations.
558  //
559  // ***NOTE*** Currently Naif stores multiples (amplitudes) as if they were additional
560  // nutation/precession terms (periods) in the equation. This method works as
561  // long as jigsaw does not solve for those values. In order to solve for
562  // those values, the multiples will need to be known so that the partial
563  // derivatives can be correctly calculated. Some possible ways of doing this
564  // are 1) Convince Naif to change their data format indicating the relation
565  // 2) Make an Isis version of the PCK data and have Isis software to
566  // calculate the rotation and partials.
567  // 3) Have an Isis addendum file that identifies the repeated periods
568  // and software to apply them when calculating the rotation and partials.
569  //
570  // For now this software will handle any terms with the same period and different
571  // amplitudes as unique terms in the equation (raNutPrec, decNutPrec,
572  // and pmNutPrec).
573  //
574  // The next three vectors will have length 3 (for a quadratic polynomial) if used.
575  std::vector<Angle>m_raPole;
576  std::vector<Angle>m_decPole;
577  std::vector<Angle>m_pm ;
578  //
579  // Currently multiples (terms with periods matching other terms but varying amplitudes)
580  // are handled as additional terms added to the end of the vector as Naif does (see
581  // comments in any of the standard Naif PCK.
582  std::vector<double>m_raNutPrec;
583  std::vector<double>m_decNutPrec;
584  std::vector<double>m_pmNutPrec;
585 
586  // The periods of bodies in the same system are modeled with a linear equation
587  std::vector<Angle>m_sysNutPrec0;
588  std::vector<Angle>m_sysNutPrec1;
589 
590  // The following scalars are used in the IAU equations to convert p_et to the appropriate time
591  // units for calculating target body ra, dec, and w. These need to be initialized in every
592  // constructor.
594  static const double m_centScale;
596  static const double m_dayScale;
597  };
598 };
599 
600 #endif
Table Cache(const QString &tableName)
Return a table with J2000 to reference rotations.
Isis specific code for binary pck.
Provide operations for quaternion arithmetic.
Definition: Quaternion.h:52
void SetOverrideBaseTime(double baseTime, double timeScale)
Set an override base time to be used with observations on scanners to allow all images in an observat...
void MinimizeCache(DownsizeStatus status)
Set the downsize status to minimize cache.
std::vector< double > EvaluatePolyFunction()
Evaluate the polynomial fit function for the three pointing angles for the current ephemeris time...
double DPckPolynomial(PartialType partialVar, const int coeffIndex)
Evaluate the derivative of the fit polynomial defined by the given coefficients with respect to the c...
std::vector< double > p_cacheTime
iTime for corresponding rotation
std::vector< double > p_CJ
Rotation matrix from J2000 to first constant rotation.
std::vector< double > Matrix()
Return the full rotation TJ as a matrix.
std::vector< Angle > m_pm
Coefficients of a quadratic polynomial fitting pole pm.
void SetPolynomialDegree(int degree)
Set the degree of the polynomials to be fit to the three camera angles for the time period covered by...
std::vector< double > m_decNutPrec
Coefficients of pole decliniation nut/prec terms.
int p_axis3
Axis of rotation for angle 3 of rotation.
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 co...
void loadPCFromTable(const PvlObject &Label)
Initialize planetary orientation constants from an cube body rotation label.
int p_targetCode
For computing Nadir rotation only.
void SetAxes(int axis1, int axis2, int axis3)
Set the axes of rotation for decomposition of a rotation matrix into 3 angles.
bool m_tOrientationAvailable
Target orientation constants are available.
std::vector< double > m_raNutPrec
Coefficients of pole right ascension nut/prec terms.
int p_degree
Degree of fit polynomial for angles.
void LoadCache(double startTime, double endTime, int size)
Cache J2000 rotation quaternion over a time range.
void setEphemerisTimePolyFunction()
When setting the ephemeris time, updates the rotation according to a polynomial that defines the thre...
int p_axis1
Axis of rotation for angle 1 of rotation.
void SetConstantMatrix(std::vector< double > constantMatrix)
Set the constant 3x3 rotation TC matrix from a vector of length 9.
void usePckPolynomial()
Set the coefficients of a polynomial fit to each of the three planet angles for the time period cover...
bool p_hasAngularVelocity
Flag indicating whether the rotation includes angular velocity.
void SetEphemerisTime(double et)
Return the J2000 to reference frame quaternion at given time.
std::vector< double > p_TC
Rotation matrix from first constant rotation (after all time-based rotations in frame chain from J200...
std::vector< Angle > m_sysNutPrec0
Constants of planetary system nut/prec periods.
std::vector< double > AngularVelocity()
Accessor method to get the angular velocity.
DownsizeStatus p_minimizeCache
Status of downsizing the cache (set to No to ignore)
void LoadTimeCache()
Load the time cache.
int Frame()
Accessor method that returns the frame code.
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.
std::vector< double > GetCenterAngles()
Return the camera angles at the center time of the observation.
std::vector< double > TimeBasedRotation()
Return time-based 3x3 rotation CJ matrix as a quaternion.
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 cover...
double WrapAngle(double compareAngle, double angle)
Wrap the input angle to keep it within 2pi radians of the angle to compare.
std::vector< double > & TimeBasedMatrix()
Return time-based 3x3 rotation CJ matrix as a vector of length 9.
std::vector< Angle > poleRaCoefs()
Return the coefficients used to calculate the target body pole ra.
SpiceRotation(int frameCode)
Construct an empty SpiceRotation class using a valid Naif frame code to set up for getting rotation f...
Isis specific code for unknown frame type.
int p_axis2
Axis of rotation for angle 2 of rotation.
std::vector< double > Extrapolate(double timeEt)
Extrapolate pointing for a given time assuming a constant angular velocity.
double p_fullCacheStartTime
Initial requested starting time of cache.
PCK frame not referenced to J2000.
std::vector< double > p_av
Angular velocity for rotation at time p_et.
FrameType
Enumeration for the frame type of the rotation.
void SetFrame(int frameCode)
Change the frame to the given frame code.
void loadPCFromSpice(int CenterBodyCode)
Initialize planetary orientation constants from Spice PCK.
double p_overrideTimeScale
Value set by caller to override computed time scale.
DownsizeStatus
Status of downsizing the cache.
bool p_noOverride
Flag to compute base time;.
double GetBaseTime()
Accessor method to get the rotation base time.
std::vector< double > poleDecNutPrecCoefs()
Return the coefficients used to calculate the target body pole dec nut/prec coefficients.
bool p_degreeApplied
Flag indicating whether or not a polynomial of degree p_degree has been created and used to fill the ...
void checkForBinaryPck()
Check loaded pck to see if any are binary and set frame type to indicate binary pck.
bool IsCached() const
Checks if the cache is empty.
void CacheLabel(Table &table)
Add labels to a SpiceRotation table.
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...
std::vector< double > ReferenceVector(const std::vector< double > &jVec)
Given a direction vector in J2000, return a reference frame direction.
double p_baseTime
Base time used in fit equations.
int p_axisP
The axis defined by the spacecraft vector for defining a nadir rotation.
bool p_matrixSet
Flag indicating p_TJ has been set.
void SetTimeBias(double timeBias)
Apply a time bias when invoking SetEphemerisTime method.
Source p_source
The source of the rotation data.
std::vector< double > GetFullCacheTime()
Return full listing (cache) of original time coverage requested.
double p_et
Current ephemeris time.
std::vector< int > TimeFrameChain()
Accessor method to get the frame chain for the rotation (begins in J2000).
static const double m_centScale
Seconds per Julian century for scaling time in seconds.
With respect to Twist or Prime Meridian Rotation.
double p_overrideBaseTime
Value set by caller to override computed base time.
See Naif Frames.req document for.
void SetTimeBasedMatrix(std::vector< double > timeBasedMatrix)
Set the time-based 3x3 rotation CJ matrix from a vector of length 9.
std::vector< std::vector< double > > p_cache
Cached rotations, stored as rotation matrix from J2000 to 1st constant frame (CJ) or coefficients of ...
void SetFullCacheParameters(double startTime, double endTime, int cacheSize)
Set the full cache time parameters.
Do not downsize the cache.
void setEphemerisTimePolyFunctionOverSpice()
When setting the ephemeris time, updates the rotation state based on a polynomial fit over spice kern...
Kernels plus nth degree polynomial.
std::vector< Angle > sysNutPrecConstants()
Return the constants used to calculate the target body system nut/prec angles.
Source
The rotation can come from one of 3 places for an Isis cube.
std::vector< std::vector< double > > p_cacheAv
Cached angular velocities for corresponding rotactions in p_cache.
FrameType getFrameType()
Accessor method to get the rotation frame type.
std::vector< Angle > pmCoefs()
Return the coefficients used to calculate the target body prime meridian.
double EphemerisTime() const
Accessor method to get current ephemeris time.
void FrameTrace(double et)
Compute frame trace chain from target frame to J2000.
Obtain SPICE rotation information for a body.
Quaternion p_quaternion
Quaternion for J2000 to reference rotation at et.
With respect to Declination.
PartialType
This enumeration indicates whether the partial derivative is taken with respect to Right Ascension...
double p_timeBias
iTime bias when reading kernels
std::vector< int > ConstantFrameChain()
Accessor method to get the frame chain for the constant part of the rotation (ends in target) ...
std::vector< double > pmNutPrecCoefs()
Return the coefficients used to calculate the target body pm nut/prec coefficients.
void ComputeBaseTime()
Compute the base time using cached times.
void setEphemerisTimeNadir()
When setting the ephemeris time, uses spacecraft nadir source to update the rotation state...
int p_fullCacheSize
Initial requested cache size.
void InitConstantRotation(double et)
Initialize the constant rotation.
double DPolynomial(const int coeffIndex)
Evaluate the derivative of the fit polynomial defined by the given coefficients with respect to the c...
void ReloadCache()
Cache J2000 rotation over existing cached time range using polynomials.
void ComputeAv()
Compute the angular velocity from the time-based functions fit to the pointing angles This method com...
Directly from the kernels.
void DCJdt(std::vector< double > &dRJ)
Compute the derivative of the 3x3 rotation matrix CJ with respect to time.
void SetSource(Source source)
Resets the source of the rotation to the given value.
std::vector< double > m_pmNutPrec
Coefficients of prime meridian nut/prec terms.
Downsize the cache.
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 co...
std::vector< double > StateTJ()
State matrix (6x6) for rotating state vectors from J2000 to target frame.
void setFrameType()
Set the frame type (m_frameType).
Class for storing Table blobs information.
Definition: Table.h:77
std::vector< double > J2000Vector(const std::vector< double > &rVec)
Given a direction vector in the reference frame, return a J2000 direction.
std::vector< Angle > m_sysNutPrec1
Linear terms of planetary system nut/prec periods.
Cache is downsized.
std::vector< Angle > sysNutPrecCoefs()
Return the coefficients used to calculate the target body system nut/prec angles. ...
FrameType m_frameType
The type of rotation frame.
Namespace for ISIS/Bullet specific routines.
Definition: Apollo.h:31
std::vector< double > p_coefficients[3]
Coefficients defining functions fit to 3 pointing angles.
From nth degree polynomial.
std::vector< double > & ConstantMatrix()
Return the constant 3x3 rotation TC matrix as a vector of length 9.
static const double m_dayScale
Seconds per day for scaling time in seconds to get target body w.
void setEphemerisTimeSpice()
When setting the ephemeris time, updates the rotation state based on data read directly from NAIF ker...
int p_axisV
The axis defined by the velocity vector for defining a nadir rotation.
std::vector< Angle > poleDecCoefs()
Return the coefficients used to calculate the target body pole dec.
double GetTimeScale()
Accessor method to get the rotation time scale.
std::vector< double > poleRaNutPrecCoefs()
Return the coefficients used to calculate the target body pole ra nut/prec coefficients.
void setEphemerisTimeMemcache()
Updates rotation state based on the rotation cache.
double p_fullCacheEndTime
Initial requested ending time of cache.
void setEphemerisTimePckPolyFunction()
When setting the ephemeris time, updates the rotation state based on the PcK polynomial.
void SetPolynomial(const Source type=PolyFunction)
Set the coefficients of a polynomial fit to each of the three camera angles for the time period cover...
Contains Pvl Groups and Pvl Objects.
Definition: PvlObject.h:74
std::vector< int > p_constantFrames
Chain of Naif frame codes in constant rotation TC.
double p_timeScale
Time scale used in fit equations.
bool HasAngularVelocity()
Checks whether the rotation has angular velocities.
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...
Table LineCache(const QString &tableName)
Return a table with J2000 to reference rotations.
std::vector< int > p_timeFrames
Chain of Naif frame codes in time-based rotation CJ.
std::vector< double > ConstantRotation()
Return the constant 3x3 rotation TC matrix as a quaternion.
Quadratic polynomial function with linear trignometric terms.
virtual ~SpiceRotation()
Destructor for SpiceRotation object.
Source GetSource()
Accessor method to get the rotation source.
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 o...
std::vector< Angle > m_raPole
Coefficients of a quadratic polynomial fitting pole ra.
With respect to Right Ascension.
std::vector< Angle > m_decPole
Coefficients of a quadratic polynomial fitting pole dec.