Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Eigen::JacobiRotation Class Reference

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

#include <ForwardDeclarations.h>

Public Types

typedef NumTraits< Scalar >::Real RealScalar
 

Public Member Functions

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

Protected Member Functions

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

Protected Attributes

Scalar m_c
 
Scalar m_s
 

Detailed Description

Rotation given by a cosine-sine pair.

\jacobi_module

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 283 of file ForwardDeclarations.h.

Member Typedef Documentation

◆ RealScalar

Definition at line 37 of file Jacobi.h.

Constructor & Destructor Documentation

◆ JacobiRotation() [1/2]

EIGEN_DEVICE_FUNC Eigen::JacobiRotation::JacobiRotation ( )
inline

Default constructor without any initialization.

Definition at line 41 of file Jacobi.h.

◆ JacobiRotation() [2/2]

EIGEN_DEVICE_FUNC Eigen::JacobiRotation::JacobiRotation ( const Scalar c,
const Scalar s 
)
inline

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

Definition at line 45 of file Jacobi.h.

Member Function Documentation

◆ adjoint()

EIGEN_DEVICE_FUNC JacobiRotation Eigen::JacobiRotation::adjoint ( ) const
inline

Returns the adjoint transformation

Definition at line 67 of file Jacobi.h.

◆ c() [1/2]

EIGEN_DEVICE_FUNC Scalar& Eigen::JacobiRotation::c ( )
inline

Definition at line 47 of file Jacobi.h.

◆ c() [2/2]

EIGEN_DEVICE_FUNC Scalar Eigen::JacobiRotation::c ( ) const
inline

Definition at line 48 of file Jacobi.h.

◆ makeGivens() [1/3]

EIGEN_DEVICE_FUNC void Eigen::JacobiRotation::makeGivens ( const Scalar p,
const Scalar q,
Scalar r,
internal::false_type   
)
protected

Definition at line 231 of file Jacobi.h.

◆ makeGivens() [2/3]

EIGEN_DEVICE_FUNC void Eigen::JacobiRotation::makeGivens ( const Scalar p,
const Scalar q,
Scalar r,
internal::true_type   
)
protected

Definition at line 171 of file Jacobi.h.

◆ makeGivens() [3/3]

EIGEN_DEVICE_FUNC void Eigen::JacobiRotation::makeGivens ( const Scalar p,
const Scalar q,
Scalar r = 0 
)

◆ makeJacobi() [1/3]

template<typename Derived >
EIGEN_DEVICE_FUNC bool Eigen::JacobiRotation::makeJacobi ( const MatrixBase< Derived > &  m,
Index  p,
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:

Matrix2f m = Matrix2f::Random();
m = (m + m.adjoint()).eval();
JacobiRotation<float> J;
J.makeJacobi(m, 0, 1);
cout << "Here is the matrix m:" << endl << m << endl;
m.applyOnTheLeft(0, 1, J.adjoint());
m.applyOnTheRight(0, 1, J);
cout << "Here is the matrix J' * m * J:" << endl << m << endl;

Output:

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

Definition at line 139 of file Jacobi.h.

◆ makeJacobi() [2/3]

EIGEN_DEVICE_FUNC bool Eigen::JacobiRotation::makeJacobi ( const RealScalar x,
const Scalar y,
const RealScalar z 
)

◆ makeJacobi() [3/3]

template<typename Scalar >
EIGEN_DEVICE_FUNC bool Eigen::JacobiRotation::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 94 of file Jacobi.h.

◆ operator*()

EIGEN_DEVICE_FUNC JacobiRotation Eigen::JacobiRotation::operator* ( const JacobiRotation other)
inline

Concatenates two planar rotation

Definition at line 54 of file Jacobi.h.

◆ s() [1/2]

EIGEN_DEVICE_FUNC Scalar& Eigen::JacobiRotation::s ( )
inline

Definition at line 49 of file Jacobi.h.

◆ s() [2/2]

EIGEN_DEVICE_FUNC Scalar Eigen::JacobiRotation::s ( ) const
inline

Definition at line 50 of file Jacobi.h.

◆ transpose()

EIGEN_DEVICE_FUNC JacobiRotation Eigen::JacobiRotation::transpose ( ) const
inline

Returns the transposed transformation

Definition at line 63 of file Jacobi.h.

Member Data Documentation

◆ m_c

Scalar Eigen::JacobiRotation::m_c
protected

Definition at line 84 of file Jacobi.h.

◆ m_s

Scalar Eigen::JacobiRotation::m_s
protected

Definition at line 84 of file Jacobi.h.


The documentation for this class was generated from the following files:
J
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
m
Matrix3f m
Definition: AngleAxis_mimic_euler.cpp:1
v
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
eval
internal::nested_eval< T, 1 >::type eval(const T &xpr)
Definition: sparse_permutations.cpp:38


gtsam
Author(s):
autogenerated on Sun Dec 22 2024 04:19:32