OPALX (Object Oriented Parallel Accelerator Library for Exascal) master (dc2a29eed580)
OPALX
Loading...
Searching...
No Matches
Quaternion Class Reference

Quaternion storage and rotation algebra used by OPALX geometry code. More...

#include <Quaternion.hpp>

Inherits ippl::Vector< double, 4 >.

Collaboration diagram for Quaternion:
Collaboration graph

Public Member Functions

 Quaternion ()
 Construct the identity rotation quaternion \((1, 0, 0, 0)\).
 
 Quaternion (const Quaternion &)
 
 Quaternion (const double &, const double &, const double &, const double &)
 Construct a quaternion from its four components \((w, x, y, z)\).
 
 Quaternion (const ippl::Vector< double, 3 > &)
 Construct a pure quaternion \((0, \mathbf{v})\) from a 3-vector.
 
 Quaternion (const double &, const ippl::Vector< double, 3 > &)
 Construct a quaternion from scalar part \(w\) and vector part \(\mathbf{v}\).
 
 Quaternion (const matrix3x3_t &)
 Construct a quaternion from a 3x3 rotation matrix.
 
Quaternion inverse () const
 Return the multiplicative inverse.
 
Quaternion conjugate () const
 Return the quaternion conjugate \(q^\ast = (w, -\mathbf{v})\).
 
Quaternion operator* (const double &) const
 
Quaternion operator* (const Quaternion &) const
 Hamilton product of two quaternions.
 
Quaternionoperator= (const Quaternion &)=default
 
Quaternionoperator*= (const Quaternion &)
 
Quaternion operator/ (const double &) const
 
double Norm () const
 Return \(\| q \|^2 = w^2 + x^2 + y^2 + z^2\).
 
double length () const
 Return the Euclidean norm \(\| q \|\).
 
Quaternionnormalize ()
 Normalize the quaternion in place.
 
bool isUnit () const
 Return true if \(\| q \|^2 \approx 1\) within a fixed tolerance.
 
bool isPure () const
 Return true if the scalar part is approximately zero.
 
bool isPureUnit () const
 Return true if the quaternion is both pure and unit length.
 
double real () const
 Return the scalar part \(w\).
 
ippl::Vector< double, 3 > imag () const
 Return the vector part \(\mathbf{v} = (x, y, z)\).
 
ippl::Vector< double, 3 > rotate (const ippl::Vector< double, 3 > &) const
 Rotate a spatial vector by a unit quaternion.
 
matrix3x3_t getRotationMatrix () const
 Convert the quaternion to a 3x3 rotation matrix.
 

Detailed Description

Quaternion storage and rotation algebra used by OPALX geometry code.

The quaternion is stored as

\[ q = (w, x, y, z) = \left(w, \mathbf{v}\right) \]

with scalar part \(w\) and vector part \(\mathbf{v} = (x, y, z)\).

For a unit quaternion, the active rotation of a spatial vector \(\mathbf{r}\) is computed as

\[ \mathbf{r}' = \left(q \, (0, \mathbf{r}) \, q^\ast\right)_{\mathrm{imag}}, \]

where \(q^\ast\) is the quaternion conjugate and the imaginary part is mapped back to \(\mathbb{R}^3\).

OPALX uses quaternions as the rotational kernel of CoordinateSystemTrafo. The conversion to a rotation matrix is defined for unit quaternions; non-unit inputs are normalized before the matrix is built.

Definition at line 29 of file Quaternion.hpp.

Constructor & Destructor Documentation

◆ Quaternion() [1/6]

Quaternion::Quaternion ( )
inline

Construct the identity rotation quaternion \((1, 0, 0, 0)\).

Definition at line 140 of file Quaternion.hpp.

Referenced by operator*=().

◆ Quaternion() [2/6]

Quaternion::Quaternion ( const Quaternion quat)
inline

Definition at line 142 of file Quaternion.hpp.

◆ Quaternion() [3/6]

Quaternion::Quaternion ( const double &  x0,
const double &  x1,
const double &  x2,
const double &  x3 
)
inline

Construct a quaternion from its four components \((w, x, y, z)\).

Definition at line 144 of file Quaternion.hpp.

◆ Quaternion() [4/6]

Quaternion::Quaternion ( const ippl::Vector< double, 3 > &  vec)
inline

Construct a pure quaternion \((0, \mathbf{v})\) from a 3-vector.

Definition at line 150 of file Quaternion.hpp.

◆ Quaternion() [5/6]

Quaternion::Quaternion ( const double &  realPart,
const ippl::Vector< double, 3 > &  vec 
)
inline

Construct a quaternion from scalar part \(w\) and vector part \(\mathbf{v}\).

Definition at line 155 of file Quaternion.hpp.

◆ Quaternion() [6/6]

Quaternion::Quaternion ( const matrix3x3_t M)

Construct a quaternion from a 3x3 rotation matrix.

Definition at line 22 of file Quaternion.cpp.

Member Function Documentation

◆ conjugate()

Quaternion Quaternion::conjugate ( ) const
inline

◆ getRotationMatrix()

matrix3x3_t Quaternion::getRotationMatrix ( ) const

Convert the quaternion to a 3x3 rotation matrix.

The quaternion is normalized first, so the returned matrix represents the same rotation for any non-zero scalar multiple of \(q\).

Definition at line 138 of file Quaternion.cpp.

References normalize().

Referenced by CoordinateSystemTrafo::operator*=(), TEST_F(), TEST_F(), and TEST_F().

Here is the call graph for this function:

◆ imag()

ippl::Vector< double, 3 > Quaternion::imag ( ) const
inline

Return the vector part \(\mathbf{v} = (x, y, z)\).

Definition at line 178 of file Quaternion.hpp.

Referenced by conjugate(), operator*(), operator*=(), operator/(), and TEST_F().

◆ inverse()

Quaternion Quaternion::inverse ( ) const

Return the multiplicative inverse.

The inverse is

\[ q^{-1} = \frac{q^\ast}{\| q \|^2}. \]

Definition at line 112 of file Quaternion.cpp.

References conjugate(), and Norm().

Referenced by TEST_F(), and TEST_F().

Here is the call graph for this function:

◆ isPure()

bool Quaternion::isPure ( ) const
inline

Return true if the scalar part is approximately zero.

Definition at line 166 of file Quaternion.hpp.

Referenced by isPureUnit(), TEST_F(), and TEST_F().

◆ isPureUnit()

bool Quaternion::isPureUnit ( ) const
inline

Return true if the quaternion is both pure and unit length.

Definition at line 168 of file Quaternion.hpp.

References isPure(), and isUnit().

Referenced by TEST_F().

Here is the call graph for this function:

◆ isUnit()

bool Quaternion::isUnit ( ) const
inline

Return true if \(\| q \|^2 \approx 1\) within a fixed tolerance.

Definition at line 164 of file Quaternion.hpp.

References Norm().

Referenced by isPureUnit(), rotate(), TEST_F(), TEST_F(), TEST_F(), TEST_F(), and TEST_F().

Here is the call graph for this function:

◆ length()

double Quaternion::length ( ) const
inline

Return the Euclidean norm \(\| q \|\).

Definition at line 162 of file Quaternion.hpp.

References Norm().

Referenced by normalize(), and TEST_F().

Here is the call graph for this function:

◆ Norm()

double Quaternion::Norm ( ) const
inline

Return \(\| q \|^2 = w^2 + x^2 + y^2 + z^2\).

Definition at line 160 of file Quaternion.hpp.

References dot().

Referenced by inverse(), isUnit(), length(), normalize(), rotate(), TEST_F(), TEST_F(), TEST_F(), and TEST_F().

Here is the call graph for this function:

◆ normalize()

Quaternion & Quaternion::normalize ( )

Normalize the quaternion in place.

After normalization,

\[ \| q \| = 1. \]

Definition at line 100 of file Quaternion.cpp.

References length(), and Norm().

Referenced by getRotationMatrix(), CoordinateSystemTrafo::operator*=(), TEST_F(), and TEST_F().

Here is the call graph for this function:

◆ operator*() [1/2]

Quaternion Quaternion::operator* ( const double &  d) const

Definition at line 71 of file Quaternion.cpp.

References imag(), and real().

Here is the call graph for this function:

◆ operator*() [2/2]

Quaternion Quaternion::operator* ( const Quaternion other) const

Hamilton product of two quaternions.

For \(q_1 = (w_1, \mathbf{v}_1)\) and \(q_2 = (w_2, \mathbf{v}_2)\),

\[ q_1 q_2 = \left(w_1 w_2 - \mathbf{v}_1 \cdot \mathbf{v}_2, w_1 \mathbf{v}_2 + w_2 \mathbf{v}_1 + \mathbf{v}_1 \times \mathbf{v}_2\right). \]

Definition at line 77 of file Quaternion.cpp.

◆ operator*=()

Quaternion & Quaternion::operator*= ( const Quaternion other)

Definition at line 82 of file Quaternion.cpp.

References cross(), imag(), and Quaternion().

Here is the call graph for this function:

◆ operator/()

Quaternion Quaternion::operator/ ( const double &  d) const

Definition at line 94 of file Quaternion.cpp.

References imag(), and real().

Here is the call graph for this function:

◆ operator=()

Quaternion & Quaternion::operator= ( const Quaternion )
default

◆ real()

double Quaternion::real ( ) const
inline

◆ rotate()

ippl::Vector< double, 3 > Quaternion::rotate ( const ippl::Vector< double, 3 > &  vec) const

Rotate a spatial vector by a unit quaternion.

This method requires a unit quaternion and evaluates

\[ \mathbf{r}' = \left(q \, (0, \mathbf{r}) \, q^\ast\right)_{\mathrm{imag}}. \]

Definition at line 125 of file Quaternion.cpp.

References conjugate(), isUnit(), and Norm().

Referenced by OpalBeamline::compileCompatibilityPlacement(), Util::getTaitBryantAngles(), CoordinateSystemTrafo::invert(), CoordinateSystemTrafo::operator*=(), CoordinateSystemTrafo::print(), TEST_F(), TEST_F(), TEST_F(), TEST_F(), TEST_F(), and TEST_F().

Here is the call graph for this function:

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