Shapeworks Studio  2.1
Shape analysis software suite
List of all members | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes
Eigen::JacobiRotation< Scalar > Class Template Reference

Rotation given by a cosine-sine pair. More...

#include <Jacobi.h>

Public Types

typedef NumTraits< Scalar >::Real RealScalar
 

Public Member Functions

 JacobiRotation ()
 
 JacobiRotation (const Scalar &c, const Scalar &s)
 
Scalar & c ()
 
Scalar c () const
 
Scalar & s ()
 
Scalar s () const
 
JacobiRotation operator* (const JacobiRotation &other)
 
JacobiRotation transpose () const
 
JacobiRotation adjoint () const
 
template<typename Derived >
bool makeJacobi (const MatrixBase< Derived > &, typename Derived::Index p, typename Derived::Index q)
 
bool makeJacobi (const RealScalar &x, const Scalar &y, const RealScalar &z)
 
void makeGivens (const Scalar &p, const Scalar &q, Scalar *z=0)
 

Protected Member Functions

void makeGivens (const Scalar &p, const Scalar &q, Scalar *z, internal::true_type)
 
void makeGivens (const Scalar &p, const Scalar &q, Scalar *z, internal::false_type)
 

Protected Attributes

Scalar m_c
 
Scalar m_s
 

Detailed Description

template<typename Scalar>
class Eigen::JacobiRotation< Scalar >

Rotation given by a cosine-sine pair.

This class represents a Jacobi or Givens rotation. This is a 2D rotation in the plane J of angle $ \theta $ defined by its cosine c and sine s as follow: $ J = \left ( \begin{array}{cc} c & \overline s \\ -s & \overline c \end{array} \right ) $

You can apply the respective counter-clockwise rotation to a column vector v by applying its adjoint on the left: $ v = J^* v $ that translates to the following Eigen code:

v.applyOnTheLeft(J.adjoint());
See also
MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()

Definition at line 228 of file ForwardDeclarations.h.

Constructor & Destructor Documentation

template<typename Scalar>
Eigen::JacobiRotation< Scalar >::JacobiRotation ( )
inline

Default constructor without any initialization.

Definition at line 40 of file Jacobi.h.

40 {}
template<typename Scalar>
Eigen::JacobiRotation< Scalar >::JacobiRotation ( const Scalar &  c,
const Scalar &  s 
)
inline

Construct a planar rotation from a cosine-sine pair (c, s).

Definition at line 43 of file Jacobi.h.

43 : m_c(c), m_s(s) {}

Member Function Documentation

template<typename Scalar>
JacobiRotation Eigen::JacobiRotation< Scalar >::adjoint ( ) const
inline

Returns the adjoint transformation

Definition at line 62 of file Jacobi.h.

62 { using numext::conj; return JacobiRotation(conj(m_c), -m_s); }
template<typename Scalar >
void Eigen::JacobiRotation< Scalar >::makeGivens ( const Scalar &  p,
const Scalar &  q,
Scalar *  z = 0 
)

Makes *this as a Givens rotation G such that applying $ G^* $ to the left of the vector $ V = \left ( \begin{array}{c} p \\ q \end{array} \right )$ yields: $ G^* V = \left ( \begin{array}{c} r \\ 0 \end{array} \right )$.

The value of z is returned if z is not null (the default is null). Also note that G is built such that the cosine is always real.

Example:

Output:

This function implements the continuous Givens rotation generation algorithm found in Anderson (2000), Discontinuous Plane Rotations and the Symmetric Eigenvalue Problem. LAPACK Working Note 150, University of Tennessee, UT-CS-00-454, December 4, 2000.

See also
MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()

Definition at line 148 of file Jacobi.h.

149 {
150  makeGivens(p, q, z, typename internal::conditional<NumTraits<Scalar>::IsComplex, internal::true_type, internal::false_type>::type());
151 }
void makeGivens(const Scalar &p, const Scalar &q, Scalar *z=0)
Definition: Jacobi.h:148
template<typename Scalar >
template<typename Derived >
bool Eigen::JacobiRotation< Scalar >::makeJacobi ( const MatrixBase< Derived > &  m,
typename Derived::Index  p,
typename Derived::Index  q 
)
inline

Makes *this as a Jacobi rotation J such that applying J on both the right and left sides of the 2x2 selfadjoint matrix $ B = \left ( \begin{array}{cc} \text{this}_{pp} & \text{this}_{pq} \\ (\text{this}_{pq})^* & \text{this}_{qq} \end{array} \right )$ yields a diagonal matrix $ A = J^* B J $

Example:

Output:

See also
JacobiRotation::makeJacobi(RealScalar, Scalar, RealScalar), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()

Definition at line 126 of file Jacobi.h.

127 {
128  return makeJacobi(numext::real(m.coeff(p,p)), m.coeff(p,q), numext::real(m.coeff(q,q)));
129 }
bool makeJacobi(const MatrixBase< Derived > &, typename Derived::Index p, typename Derived::Index q)
Definition: Jacobi.h:126
template<typename Scalar >
bool Eigen::JacobiRotation< Scalar >::makeJacobi ( const RealScalar &  x,
const Scalar &  y,
const RealScalar &  z 
)

Makes *this as a Jacobi rotation J such that applying J on both the right and left sides of the selfadjoint 2x2 matrix $ B = \left ( \begin{array}{cc} x & y \\ \overline y & z \end{array} \right )$ yields a diagonal matrix $ A = J^* B J $

See also
MatrixBase::makeJacobi(const MatrixBase<Derived>&, Index, Index), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()

Definition at line 83 of file Jacobi.h.

84 {
85  using std::sqrt;
86  using std::abs;
87  typedef typename NumTraits<Scalar>::Real RealScalar;
88  if(y == Scalar(0))
89  {
90  m_c = Scalar(1);
91  m_s = Scalar(0);
92  return false;
93  }
94  else
95  {
96  RealScalar tau = (x-z)/(RealScalar(2)*abs(y));
97  RealScalar w = sqrt(numext::abs2(tau) + RealScalar(1));
98  RealScalar t;
99  if(tau>RealScalar(0))
100  {
101  t = RealScalar(1) / (tau + w);
102  }
103  else
104  {
105  t = RealScalar(1) / (tau - w);
106  }
107  RealScalar sign_t = t > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
108  RealScalar n = RealScalar(1) / sqrt(numext::abs2(t)+RealScalar(1));
109  m_s = - sign_t * (numext::conj(y) / abs(y)) * abs(t) * n;
110  m_c = n;
111  return true;
112  }
113 }
template<typename Scalar>
JacobiRotation Eigen::JacobiRotation< Scalar >::operator* ( const JacobiRotation< Scalar > &  other)
inline

Concatenates two planar rotation

Definition at line 51 of file Jacobi.h.

52  {
53  using numext::conj;
54  return JacobiRotation(m_c * other.m_c - conj(m_s) * other.m_s,
55  conj(m_c * conj(other.m_s) + conj(m_s) * conj(other.m_c)));
56  }
template<typename Scalar>
JacobiRotation Eigen::JacobiRotation< Scalar >::transpose ( ) const
inline

Returns the transposed transformation

Definition at line 59 of file Jacobi.h.

59 { using numext::conj; return JacobiRotation(m_c, -conj(m_s)); }

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