JSBSim Flight Dynamics Model  1.0 (02 March 2017)
An Open Source Flight Dynamics and Control Software Library in C++
FGQuaternion Class Reference

Models the Quaternion representation of rotations. More...

#include <FGQuaternion.h>

+ Inheritance diagram for FGQuaternion:
+ Collaboration diagram for FGQuaternion:

Public Member Functions

 FGQuaternion ()
 Default initializer. More...
 
 FGQuaternion (const FGQuaternion &q)
 Copy constructor. More...
 
 FGQuaternion (double phi, double tht, double psi)
 Initializer by euler angles. More...
 
 FGQuaternion (FGColumnVector3 vOrient)
 Initializer by euler angle vector. More...
 
 FGQuaternion (int idx, double angle)
 Initializer by one euler angle. More...
 
 FGQuaternion (double angle, const FGColumnVector3 &axis)
 Initializer by a rotation axis and an angle. More...
 
 FGQuaternion (const FGMatrix33 &m)
 Initializer by matrix. More...
 
 ~FGQuaternion ()
 Destructor.
 
FGQuaternion Conjugate (void) const
 Conjugate of the quaternion. More...
 
std::string Dump (const std::string &delimiter) const
 
double Entry (unsigned int idx) const
 Read access the entries of the vector. More...
 
double & Entry (unsigned int idx)
 Write access the entries of the vector. More...
 
double GetCosEuler (int i) const
 Retrieves cosine of the given euler angle. More...
 
const FGColumnVector3GetEuler (void) const
 Retrieves the Euler angles. More...
 
double GetEuler (int i) const
 Retrieves the Euler angles. More...
 
double GetEulerDeg (int i) const
 Retrieves the Euler angles. More...
 
FGColumnVector3 const GetEulerDeg (void) const
 Retrieves the Euler angle vector. More...
 
FGQuaternion GetQDot (const FGColumnVector3 &PQR) const
 Quaternion derivative for given angular rates. More...
 
double GetSinEuler (int i) const
 Retrieves sine of the given euler angle. More...
 
const FGMatrix33GetT (void) const
 Transformation matrix. More...
 
const FGMatrix33GetTInv (void) const
 Backward transformation matrix. More...
 
FGQuaternion Inverse (void) const
 Inverse of the quaternion. More...
 
double Magnitude (void) const
 Length of the vector. More...
 
void Normalize (void)
 Normalize. More...
 
 operator FGMatrix33 () const
 Conversion from Quat to Matrix.
 
bool operator!= (const FGQuaternion &q) const
 Comparison operator "!=". More...
 
double operator() (unsigned int idx) const
 Read access the entries of the vector. More...
 
double & operator() (unsigned int idx)
 Write access the entries of the vector. More...
 
FGQuaternion operator* (const FGQuaternion &q) const
 Arithmetic operator "*". More...
 
const FGQuaternionoperator*= (double scalar)
 Arithmetic operator "*=". More...
 
const FGQuaternionoperator*= (const FGQuaternion &q)
 Arithmetic operator "*=". More...
 
FGQuaternion operator+ (const FGQuaternion &q) const
 Arithmetic operator "+". More...
 
const FGQuaternionoperator+= (const FGQuaternion &q)
 
FGQuaternion operator- (const FGQuaternion &q) const
 Arithmetic operator "-". More...
 
const FGQuaternionoperator-= (const FGQuaternion &q)
 Arithmetic operator "-=". More...
 
const FGQuaternionoperator/= (double scalar)
 Arithmetic operator "/=". More...
 
const FGQuaternionoperator= (const FGQuaternion &q)
 Assignment operator "=". More...
 
bool operator== (const FGQuaternion &q) const
 Comparison operator "==". More...
 
double SqrMagnitude (void) const
 Square of the length of the vector. More...
 
- Public Member Functions inherited from FGJSBBase
 FGJSBBase ()
 Constructor for FGJSBBase.
 
virtual ~FGJSBBase ()
 Destructor for FGJSBBase.
 
void disableHighLighting (void)
 Disables highlighting in the console output.
 
std::string GetVersion (void)
 Returns the version number of JSBSim. More...
 
void PutMessage (const Message &msg)
 Places a Message structure on the Message queue. More...
 
void PutMessage (const std::string &text)
 Creates a message with the given text and places it on the queue. More...
 
void PutMessage (const std::string &text, bool bVal)
 Creates a message with the given text and boolean value and places it on the queue. More...
 
void PutMessage (const std::string &text, int iVal)
 Creates a message with the given text and integer value and places it on the queue. More...
 
void PutMessage (const std::string &text, double dVal)
 Creates a message with the given text and double value and places it on the queue. More...
 
int SomeMessages (void)
 Reads the message on the queue (but does not delete it). More...
 
void ProcessMessage (void)
 Reads the message on the queue and removes it from the queue. More...
 
MessageProcessNextMessage (void)
 Reads the next message on the queue and removes it from the queue. More...
 

Static Public Member Functions

static FGQuaternion zero (void)
 Zero quaternion vector. More...
 
- Static Public Member Functions inherited from FGJSBBase
static double CelsiusToFahrenheit (double celsius)
 Converts from degrees Celsius to degrees Fahrenheit. More...
 
static double CelsiusToKelvin (double celsius)
 Converts from degrees Celsius to degrees Kelvin. More...
 
static double CelsiusToRankine (double celsius)
 Converts from degrees Celsius to degrees Rankine. More...
 
static double Constrain (double min, double value, double max)
 Constrain a value between a minimum and a maximum value.
 
static bool EqualToRoundoff (double a, double b)
 Finite precision comparison. More...
 
static bool EqualToRoundoff (float a, float b)
 Finite precision comparison. More...
 
static bool EqualToRoundoff (float a, double b)
 Finite precision comparison. More...
 
static bool EqualToRoundoff (double a, float b)
 Finite precision comparison. More...
 
static double FahrenheitToCelsius (double fahrenheit)
 Converts from degrees Fahrenheit to degrees Celsius. More...
 
static double FeetToMeters (double measure)
 Converts from feet to meters. More...
 
static double GaussianRandomNumber (void)
 
static double KelvinToCelsius (double kelvin)
 Converts from degrees Kelvin to degrees Celsius. More...
 
static double KelvinToFahrenheit (double kelvin)
 Converts from degrees Kelvin to degrees Fahrenheit. More...
 
static double KelvinToRankine (double kelvin)
 Converts from degrees Kelvin to degrees Rankine. More...
 
static double MachFromVcalibrated (double vcas, double p, double psl, double rhosl)
 Calculate the Mach number from the calibrated airspeed. More...
 
static double PitotTotalPressure (double mach, double p)
 Compute the total pressure in front of the Pitot tube. More...
 
static double RankineToCelsius (double rankine)
 Converts from degrees Rankine to degrees Celsius. More...
 
static double RankineToKelvin (double rankine)
 Converts from degrees Rankine to degrees Kelvin. More...
 
static double sign (double num)
 
static double VcalibratedFromMach (double mach, double p, double psl, double rhosl)
 Calculate the calibrated airspeed from the Mach number. More...
 

Friends

FGQuaternion operator* (double, const FGQuaternion &)
 Scalar multiplication. More...
 
FGQuaternion QExp (const FGColumnVector3 &omega)
 Quaternion exponential. More...
 

Additional Inherited Members

- Public Types inherited from FGJSBBase
enum  { eL = 1, eM, eN }
 Moments L, M, N.
 
enum  { eP = 1, eQ, eR }
 Rates P, Q, R.
 
enum  { eU = 1, eV, eW }
 Velocities U, V, W.
 
enum  { eX = 1, eY, eZ }
 Positions X, Y, Z.
 
enum  { ePhi = 1, eTht, ePsi }
 Euler angles Phi, Theta, Psi.
 
enum  { eDrag = 1, eSide, eLift }
 Stability axis forces, Drag, Side force, Lift.
 
enum  { eRoll = 1, ePitch, eYaw }
 Local frame orientation Roll, Pitch, Yaw.
 
enum  { eNorth = 1, eEast, eDown }
 Local frame position North, East, Down.
 
enum  { eLat = 1, eLong, eRad }
 Locations Radius, Latitude, Longitude.
 
enum  {
  inNone = 0, inDegrees, inRadians, inMeters,
  inFeet
}
 Conversion specifiers.
 
- Static Public Attributes inherited from FGJSBBase
static short debug_lvl = 1
 
static char highint [5] = {27, '[', '1', 'm', '\0' }
 highlights text
 
static char halfint [5] = {27, '[', '2', 'm', '\0' }
 low intensity text
 
static char normint [6] = {27, '[', '2', '2', 'm', '\0' }
 normal intensity text
 
static char reset [5] = {27, '[', '0', 'm', '\0' }
 resets text properties
 
static char underon [5] = {27, '[', '4', 'm', '\0' }
 underlines text
 
static char underoff [6] = {27, '[', '2', '4', 'm', '\0' }
 underline off
 
static char fgblue [6] = {27, '[', '3', '4', 'm', '\0' }
 blue text
 
static char fgcyan [6] = {27, '[', '3', '6', 'm', '\0' }
 cyan text
 
static char fgred [6] = {27, '[', '3', '1', 'm', '\0' }
 red text
 
static char fggreen [6] = {27, '[', '3', '2', 'm', '\0' }
 green text
 
static char fgdef [6] = {27, '[', '3', '9', 'm', '\0' }
 default text
 
- Protected Member Functions inherited from FGJSBBase
void Debug (int)
 
- Static Protected Member Functions inherited from FGJSBBase
static std::string CreateIndexedPropertyName (const std::string &Property, int index)
 
- Static Protected Attributes inherited from FGJSBBase
static const double degtorad = 0.017453292519943295769236907684886
 
static const double fpstokts = 1.0/ktstofps
 
static const double fttom = 0.3048
 
static int gaussian_random_number_phase = 0
 
static const double hptoftlbssec = 550.0
 
static const double in3tom3 = 1.638706E-5
 
static const double inchtoft = 0.08333333
 
static const double inhgtopa = 3386.38
 
static const std::string JSBSim_version = "1.0 " __DATE__ " " __TIME__
 
static const double kgtolb = 2.20462
 
static const double kgtoslug = 0.06852168
 
static const double ktstofps = 1.68781
 
static const double lbtoslug = 1.0/slugtolb
 
static Message localMsg
 
static const double m3toft3 = 1.0/(fttom*fttom*fttom)
 
static double Mair = 28.9645
 
static unsigned int messageId = 0
 
static std::queue< MessageMessages
 
static const std::string needed_cfg_version = "2.0"
 
static const double psftoinhg = 0.014138
 
static const double psftopa = 47.88
 
static const double radtodeg = 57.295779513082320876798154814105
 
static double Reng = 1716.56
 
static double Rstar = 1545.348
 
static const double SHRatio = 1.40
 
static const double slugtolb = 32.174049
 

Detailed Description

Models the Quaternion representation of rotations.

FGQuaternion is a representation of an arbitrary rotation through a quaternion. It has vector properties. This class also contains access functions to the euler angle representation of rotations and access to transformation matrices for 3D vectors. Transformations and euler angles are therefore computed once they are requested for the first time. Then they are cached for later usage as long as the class is not accessed trough a nonconst member function.

Note: The order of rotations used in this class corresponds to a 3-2-1 sequence, or Y-P-R, or Z-Y-X, if you prefer.

See also
Cooke, Zyda, Pratt, and McGhee, "NPSNET: Flight Simulation Dynamic Modeling Using Quaternions", Presence, Vol. 1, No. 4, pp. 404-420 Naval Postgraduate School, January 1994
D. M. Henderson, "Euler Angles, Quaternions, and Transformation Matrices", JSC 12960, July 1977
Richard E. McFarland, "A Standard Kinematic Model for Flight Simulation at NASA-Ames", NASA CR-2497, January 1975
Barnes W. McCormick, "Aerodynamics, Aeronautics, and Flight Mechanics", Wiley & Sons, 1979 ISBN 0-471-03032-5
Bernard Etkin, "Dynamics of Flight, Stability and Control", Wiley & Sons, 1982 ISBN 0-471-08936-2
Author
Mathias Froehlich, extended FGColumnVector4 originally by Tony Peden and Jon Berndt

Definition at line 92 of file FGQuaternion.h.

Constructor & Destructor Documentation

◆ FGQuaternion() [1/7]

FGQuaternion ( )
inline

Default initializer.

Default initializer, initializes the class with the identity rotation.

Definition at line 96 of file FGQuaternion.h.

96  : mCacheValid(false) {
97  data[0] = 1.0;
98  data[1] = data[2] = data[3] = 0.0;
99  }
+ Here is the caller graph for this function:

◆ FGQuaternion() [2/7]

FGQuaternion ( const FGQuaternion q)

Copy constructor.

Copy constructor, initializes the quaternion.

Parameters
qa constant reference to another FGQuaternion instance

Definition at line 67 of file FGQuaternion.cpp.

67  : mCacheValid(q.mCacheValid)
68 {
69  data[0] = q(1);
70  data[1] = q(2);
71  data[2] = q(3);
72  data[3] = q(4);
73  if (mCacheValid) {
74  mT = q.mT;
75  mTInv = q.mTInv;
76  mEulerAngles = q.mEulerAngles;
77  mEulerSines = q.mEulerSines;
78  mEulerCosines = q.mEulerCosines;
79  }
80 }

◆ FGQuaternion() [3/7]

FGQuaternion ( double  phi,
double  tht,
double  psi 
)

Initializer by euler angles.

Initialize the quaternion with the euler angles.

Parameters
phiThe euler X axis (roll) angle in radians
thtThe euler Y axis (attitude) angle in radians
psiThe euler Z axis (heading) angle in radians

Definition at line 85 of file FGQuaternion.cpp.

85  : mCacheValid(false)
86 {
87  InitializeFromEulerAngles(phi, tht, psi);
88 }

◆ FGQuaternion() [4/7]

Initializer by euler angle vector.

Initialize the quaternion with the euler angle vector.

Parameters
vOrientThe euler axis angle vector in radians (phi, tht, psi)

Definition at line 92 of file FGQuaternion.cpp.

92  : mCacheValid(false)
93 {
94  double phi = vOrient(ePhi);
95  double tht = vOrient(eTht);
96  double psi = vOrient(ePsi);
97 
98  InitializeFromEulerAngles(phi, tht, psi);
99 }
+ Here is the call graph for this function:

◆ FGQuaternion() [5/7]

FGQuaternion ( int  idx,
double  angle 
)
inline

Initializer by one euler angle.

Initialize the quaternion with the single euler angle where its index is given in the first argument.

Parameters
idxIndex of the euler angle to initialize
angleThe euler angle in radians

Definition at line 123 of file FGQuaternion.h.

124  : mCacheValid(false) {
125 
126  double angle2 = 0.5*angle;
127 
128  double Sangle2 = sin(angle2);
129  double Cangle2 = cos(angle2);
130 
131  if (idx == ePhi) {
132  data[0] = Cangle2;
133  data[1] = Sangle2;
134  data[2] = 0.0;
135  data[3] = 0.0;
136 
137  } else if (idx == eTht) {
138  data[0] = Cangle2;
139  data[1] = 0.0;
140  data[2] = Sangle2;
141  data[3] = 0.0;
142 
143  } else {
144  data[0] = Cangle2;
145  data[1] = 0.0;
146  data[2] = 0.0;
147  data[3] = Sangle2;
148 
149  }
150  }

◆ FGQuaternion() [6/7]

FGQuaternion ( double  angle,
const FGColumnVector3 axis 
)
inline

Initializer by a rotation axis and an angle.

Initialize the quaternion to represent the rotation around a given angle and an arbitrary axis.

Parameters
angleThe angle in radians
axisThe rotation axis

Definition at line 158 of file FGQuaternion.h.

159  : mCacheValid(false) {
160 
161  double angle2 = 0.5 * angle;
162 
163  double length = axis.Magnitude();
164  double Sangle2 = sin(angle2) / length;
165  double Cangle2 = cos(angle2);
166 
167  data[0] = Cangle2;
168  data[1] = Sangle2 * axis(1);
169  data[2] = Sangle2 * axis(2);
170  data[3] = Sangle2 * axis(3);
171  }
+ Here is the call graph for this function:

◆ FGQuaternion() [7/7]

FGQuaternion ( const FGMatrix33 m)

Initializer by matrix.

Initialize the quaternion with the matrix representing a transform from one frame to another using the standard aerospace sequence, Yaw-Pitch-Roll (3-2-1).

Parameters
mthe rotation matrix

Definition at line 142 of file FGQuaternion.cpp.

142  : mCacheValid(false)
143 {
144  data[0] = 0.50*sqrt(1.0 + m(1,1) + m(2,2) + m(3,3));
145  double t = 0.25/data[0];
146  data[1] = t*(m(2,3) - m(3,2));
147  data[2] = t*(m(3,1) - m(1,3));
148  data[3] = t*(m(1,2) - m(2,1));
149 
150  Normalize();
151 }
void Normalize(void)
Normalize.
+ Here is the call graph for this function:

Member Function Documentation

◆ Conjugate()

FGQuaternion Conjugate ( void  ) const
inline

Conjugate of the quaternion.

Compute and return the conjugate of the quaternion. This one is equal to the inverse iff the quaternion is normalized.

Definition at line 456 of file FGQuaternion.h.

456  {
457  return FGQuaternion( data[0], -data[1], -data[2], -data[3] );
458  }
FGQuaternion()
Default initializer.
Definition: FGQuaternion.h:96
+ Here is the call graph for this function:

◆ Entry() [1/2]

double Entry ( unsigned int  idx) const
inline

Read access the entries of the vector.

Parameters
idxthe component index.

Return the value of the matrix entry at the given index. Indices are counted starting with 1.

This function is just a shortcut for the double operator()(unsigned int idx) const function. It is used internally to access the elements in a more convenient way.

Note that the index given in the argument is unchecked.

Definition at line 291 of file FGQuaternion.h.

291 { return data[idx-1]; }

◆ Entry() [2/2]

double& Entry ( unsigned int  idx)
inline

Write access the entries of the vector.

Parameters
idxthe component index.

Return a reference to the vector entry at the given index. Indices are counted starting with 1.

This function is just a shortcut for the double& operator()(unsigned int idx) function. It is used internally to access the elements in a more convenient way.

Note that the index given in the argument is unchecked.

Definition at line 306 of file FGQuaternion.h.

306  {
307  mCacheValid = false;
308  return data[idx-1];
309  }

◆ GetCosEuler()

double GetCosEuler ( int  i) const
inline

Retrieves cosine of the given euler angle.

Returns
the sine of the Euler angle theta (pitch attitude) corresponding to this quaternion rotation.

Definition at line 251 of file FGQuaternion.h.

251  {
252  ComputeDerived();
253  return mEulerCosines(i);
254  }
+ Here is the caller graph for this function:

◆ GetEuler() [1/2]

const FGColumnVector3& GetEuler ( void  ) const
inline

Retrieves the Euler angles.

Returns
a reference to the triad of Euler angles corresponding to this quaternion rotation. units radians

Definition at line 205 of file FGQuaternion.h.

205  {
206  ComputeDerived();
207  return mEulerAngles;
208  }
+ Here is the caller graph for this function:

◆ GetEuler() [2/2]

double GetEuler ( int  i) const
inline

Retrieves the Euler angles.

Parameters
ithe Euler angle index. units radians.
Returns
a reference to the i-th euler angles corresponding to this quaternion rotation.

Definition at line 216 of file FGQuaternion.h.

216  {
217  ComputeDerived();
218  return mEulerAngles(i);
219  }

◆ GetEulerDeg() [1/2]

double GetEulerDeg ( int  i) const
inline

Retrieves the Euler angles.

Parameters
ithe Euler angle index.
Returns
a reference to the i-th euler angles corresponding to this quaternion rotation. units degrees

Definition at line 226 of file FGQuaternion.h.

226  {
227  ComputeDerived();
228  return radtodeg*mEulerAngles(i);
229  }
+ Here is the caller graph for this function:

◆ GetEulerDeg() [2/2]

FGColumnVector3 const GetEulerDeg ( void  ) const
inline

Retrieves the Euler angle vector.

Returns
an Euler angle column vector corresponding to this quaternion rotation. units degrees

Definition at line 235 of file FGQuaternion.h.

235  {
236  ComputeDerived();
237  return radtodeg*mEulerAngles;
238  }

◆ GetQDot()

FGQuaternion GetQDot ( const FGColumnVector3 PQR) const

Quaternion derivative for given angular rates.

Returns the derivative of the quaternion corresponding to the angular velocities PQR.

Computes the quaternion derivative which results from the given angular velocities

Parameters
PQRa constant reference to a rotation rate vector
Returns
the quaternion derivative
See also
Stevens and Lewis, "Aircraft Control and Simulation", Second Edition, Equation 1.3-36.

See Stevens and Lewis, "Aircraft Control and Simulation", Second Edition, Equation 1.3-36. Also see Jack Kuipers, "Quaternions and Rotation Sequences", Equation 11.12.

Definition at line 161 of file FGQuaternion.cpp.

162 {
163  return FGQuaternion(
164  -0.5*( data[1]*PQR(eP) + data[2]*PQR(eQ) + data[3]*PQR(eR)),
165  0.5*( data[0]*PQR(eP) - data[3]*PQR(eQ) + data[2]*PQR(eR)),
166  0.5*( data[3]*PQR(eP) + data[0]*PQR(eQ) - data[1]*PQR(eR)),
167  0.5*(-data[2]*PQR(eP) + data[1]*PQR(eQ) + data[0]*PQR(eR))
168  );
169 }
FGQuaternion()
Default initializer.
Definition: FGQuaternion.h:96
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetSinEuler()

double GetSinEuler ( int  i) const
inline

Retrieves sine of the given euler angle.

Returns
the sine of the Euler angle theta (pitch attitude) corresponding to this quaternion rotation.

Definition at line 243 of file FGQuaternion.h.

243  {
244  ComputeDerived();
245  return mEulerSines(i);
246  }
+ Here is the caller graph for this function:

◆ GetT()

const FGMatrix33& GetT ( void  ) const
inline

Transformation matrix.

Returns
a reference to the transformation/rotation matrix corresponding to this quaternion rotation.

Definition at line 194 of file FGQuaternion.h.

194 { ComputeDerived(); return mT; }
+ Here is the caller graph for this function:

◆ GetTInv()

const FGMatrix33& GetTInv ( void  ) const
inline

Backward transformation matrix.

Returns
a reference to the inverse transformation/rotation matrix corresponding to this quaternion rotation.

Definition at line 199 of file FGQuaternion.h.

199 { ComputeDerived(); return mTInv; }
+ Here is the caller graph for this function:

◆ Inverse()

FGQuaternion Inverse ( void  ) const
inline

Inverse of the quaternion.

Compute and return the inverse of the quaternion so that the orientation represented with *this multiplied with the returned value is equal to the identity orientation.

Definition at line 442 of file FGQuaternion.h.

442  {
443  double norm = SqrMagnitude();
444  if (norm == 0.0)
445  return *this;
446  double rNorm = 1.0/norm;
447  return FGQuaternion( data[0]*rNorm, -data[1]*rNorm,
448  -data[2]*rNorm, -data[3]*rNorm );
449  }
FGQuaternion()
Default initializer.
Definition: FGQuaternion.h:96
double SqrMagnitude(void) const
Square of the length of the vector.
Definition: FGQuaternion.h:472
+ Here is the call graph for this function:

◆ Magnitude()

double Magnitude ( void  ) const
inline

Length of the vector.

Compute and return the euclidean norm of this vector.

Definition at line 466 of file FGQuaternion.h.

466 { return sqrt(SqrMagnitude()); }
double SqrMagnitude(void) const
Square of the length of the vector.
Definition: FGQuaternion.h:472
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ Normalize()

void Normalize ( void  )

Normalize.

Normalize the vector to have the Magnitude() == 1.0. If the vector is equal to zero it is left untouched.

Definition at line 173 of file FGQuaternion.cpp.

174 {
175  // Note: this does not touch the cache since it does not change the orientation
176  double norm = Magnitude();
177  if (norm == 0.0 || fabs(norm - 1.000) < 1e-10) return;
178 
179  double rnorm = 1.0/norm;
180 
181  data[0] *= rnorm;
182  data[1] *= rnorm;
183  data[2] *= rnorm;
184  data[3] *= rnorm;
185 }
double Magnitude(void) const
Length of the vector.
Definition: FGQuaternion.h:466
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ operator!=()

bool operator!= ( const FGQuaternion q) const
inline

Comparison operator "!=".

Parameters
qa quaternion reference
Returns
true if both quaternions do not represent the same rotation.

Definition at line 349 of file FGQuaternion.h.

349 { return ! operator==(q); }
bool operator==(const FGQuaternion &q) const
Comparison operator "==".
Definition: FGQuaternion.h:341
+ Here is the call graph for this function:

◆ operator()() [1/2]

double operator() ( unsigned int  idx) const
inline

Read access the entries of the vector.

Parameters
idxthe component index.

Return the value of the matrix entry at the given index. Indices are counted starting with 1.

Note that the index given in the argument is unchecked.

Definition at line 265 of file FGQuaternion.h.

265 { return data[idx-1]; }

◆ operator()() [2/2]

double& operator() ( unsigned int  idx)
inline

Write access the entries of the vector.

Parameters
idxthe component index.

Return a reference to the vector entry at the given index. Indices are counted starting with 1.

Note that the index given in the argument is unchecked.

Definition at line 276 of file FGQuaternion.h.

276 { mCacheValid = false; return data[idx-1]; }

◆ operator*()

FGQuaternion operator* ( const FGQuaternion q) const
inline

Arithmetic operator "*".

Multiplication of two quaternions is like performing successive rotations.

Parameters
qa quaternion to be multiplied.
Returns
a quaternion representing Q, where Q = Q * q.

Definition at line 412 of file FGQuaternion.h.

412  {
413  return FGQuaternion(data[0]*q.data[0]-data[1]*q.data[1]-data[2]*q.data[2]-data[3]*q.data[3],
414  data[0]*q.data[1]+data[1]*q.data[0]+data[2]*q.data[3]-data[3]*q.data[2],
415  data[0]*q.data[2]-data[1]*q.data[3]+data[2]*q.data[0]+data[3]*q.data[1],
416  data[0]*q.data[3]+data[1]*q.data[2]-data[2]*q.data[1]+data[3]*q.data[0]);
417  }
FGQuaternion()
Default initializer.
Definition: FGQuaternion.h:96
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ operator*=() [1/2]

const FGQuaternion& operator*= ( double  scalar)
inline

Arithmetic operator "*=".

Parameters
scalara multiplicative value.
Returns
a quaternion reference representing Q, where Q = Q * scalar.

Definition at line 376 of file FGQuaternion.h.

376  {
377  data[0] *= scalar;
378  data[1] *= scalar;
379  data[2] *= scalar;
380  data[3] *= scalar;
381  mCacheValid = false;
382  return *this;
383  }
+ Here is the caller graph for this function:

◆ operator*=() [2/2]

const FGQuaternion& operator*= ( const FGQuaternion q)
inline

Arithmetic operator "*=".

Multiplication of two quaternions is like performing successive rotations.

Parameters
qa quaternion to be multiplied.
Returns
a quaternion reference representing Q, where Q = Q * q.

Definition at line 423 of file FGQuaternion.h.

423  {
424  double q0 = data[0]*q.data[0]-data[1]*q.data[1]-data[2]*q.data[2]-data[3]*q.data[3];
425  double q1 = data[0]*q.data[1]+data[1]*q.data[0]+data[2]*q.data[3]-data[3]*q.data[2];
426  double q2 = data[0]*q.data[2]-data[1]*q.data[3]+data[2]*q.data[0]+data[3]*q.data[1];
427  double q3 = data[0]*q.data[3]+data[1]*q.data[2]-data[2]*q.data[1]+data[3]*q.data[0];
428  data[0] = q0;
429  data[1] = q1;
430  data[2] = q2;
431  data[3] = q3;
432  mCacheValid = false;
433  return *this;
434  }

◆ operator+()

FGQuaternion operator+ ( const FGQuaternion q) const
inline

Arithmetic operator "+".

Parameters
qa quaternion to be summed.
Returns
a quaternion representing Q, where Q = Q + q.

Definition at line 395 of file FGQuaternion.h.

395  {
396  return FGQuaternion(data[0]+q.data[0], data[1]+q.data[1],
397  data[2]+q.data[2], data[3]+q.data[3]);
398  }
FGQuaternion()
Default initializer.
Definition: FGQuaternion.h:96
+ Here is the call graph for this function:

◆ operator-()

FGQuaternion operator- ( const FGQuaternion q) const
inline

Arithmetic operator "-".

Parameters
qa quaternion to be subtracted.
Returns
a quaternion representing Q, where Q = Q - q.

Definition at line 403 of file FGQuaternion.h.

403  {
404  return FGQuaternion(data[0]-q.data[0], data[1]-q.data[1],
405  data[2]-q.data[2], data[3]-q.data[3]);
406  }
FGQuaternion()
Default initializer.
Definition: FGQuaternion.h:96
+ Here is the call graph for this function:

◆ operator-=()

const FGQuaternion& operator-= ( const FGQuaternion q)
inline

Arithmetic operator "-=".

Parameters
qa quaternion reference.
Returns
a quaternion reference representing Q, where Q = Q - q.

Definition at line 363 of file FGQuaternion.h.

363  {
364  // Copy the master values ...
365  data[0] -= q.data[0];
366  data[1] -= q.data[1];
367  data[2] -= q.data[2];
368  data[3] -= q.data[3];
369  mCacheValid = false;
370  return *this;
371  }

◆ operator/=()

const FGQuaternion& operator/= ( double  scalar)
inline

Arithmetic operator "/=".

Parameters
scalara divisor value.
Returns
a quaternion reference representing Q, where Q = Q / scalar.

Definition at line 388 of file FGQuaternion.h.

388  {
389  return operator*=(1.0/scalar);
390  }
const FGQuaternion & operator*=(double scalar)
Arithmetic operator "*=".
Definition: FGQuaternion.h:376
+ Here is the call graph for this function:

◆ operator=()

const FGQuaternion& operator= ( const FGQuaternion q)
inline

Assignment operator "=".

Assign the value of q to the current object. Cached values are conserved.

Parameters
qreference to an FGQuaternion instance
Returns
reference to a quaternion object

Definition at line 316 of file FGQuaternion.h.

316  {
317  // Copy the master values ...
318  data[0] = q.data[0];
319  data[1] = q.data[1];
320  data[2] = q.data[2];
321  data[3] = q.data[3];
322  ComputeDerived();
323  // .. and copy the derived values if they are valid
324  mCacheValid = q.mCacheValid;
325  if (mCacheValid) {
326  mT = q.mT;
327  mTInv = q.mTInv;
328  mEulerAngles = q.mEulerAngles;
329  mEulerSines = q.mEulerSines;
330  mEulerCosines = q.mEulerCosines;
331  }
332  return *this;
333  }

◆ operator==()

bool operator== ( const FGQuaternion q) const
inline

Comparison operator "==".

Parameters
qa quaternion reference
Returns
true if both quaternions represent the same rotation.

Definition at line 341 of file FGQuaternion.h.

341  {
342  return data[0] == q.data[0] && data[1] == q.data[1]
343  && data[2] == q.data[2] && data[3] == q.data[3];
344  }
+ Here is the caller graph for this function:

◆ SqrMagnitude()

double SqrMagnitude ( void  ) const
inline

Square of the length of the vector.

Compute and return the square of the euclidean norm of this vector.

Definition at line 472 of file FGQuaternion.h.

472  {
473  return data[0]*data[0] + data[1]*data[1]
474  + data[2]*data[2] + data[3]*data[3];
475  }
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ zero()

static FGQuaternion zero ( void  )
inlinestatic

Zero quaternion vector.

Does not represent any orientation. Useful for initialization of increments

Definition at line 486 of file FGQuaternion.h.

486 { return FGQuaternion( 0.0, 0.0, 0.0, 0.0 ); }
FGQuaternion()
Default initializer.
Definition: FGQuaternion.h:96
+ Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator*

FGQuaternion operator* ( double  scalar,
const FGQuaternion q 
)
friend

Scalar multiplication.

Parameters
scalarscalar value to multiply with.
qVector to multiply.

Multiply the Vector with a scalar value.

Definition at line 545 of file FGQuaternion.h.

545  {
546  return FGQuaternion(scalar*q.data[0], scalar*q.data[1], scalar*q.data[2], scalar*q.data[3]);
547 }
FGQuaternion()
Default initializer.
Definition: FGQuaternion.h:96

◆ QExp

FGQuaternion QExp ( const FGColumnVector3 omega)
friend

Quaternion exponential.

Parameters
omegarotation velocity Calculate the unit quaternion which is the result of the exponentiation of the vector 'omega'.

Definition at line 554 of file FGQuaternion.h.

554  {
555  FGQuaternion qexp;
556  double angle = omega.Magnitude();
557  double sina_a = angle > 0.0 ? sin(angle)/angle : 1.0;
558 
559  qexp.data[0] = cos(angle);
560  qexp.data[1] = omega(1) * sina_a;
561  qexp.data[2] = omega(2) * sina_a;
562  qexp.data[3] = omega(3) * sina_a;
563 
564  return qexp;
565 }
FGQuaternion()
Default initializer.
Definition: FGQuaternion.h:96

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