Shapeworks Studio  2.1
Shape analysis software suite
List of all members | Public Types | Public Member Functions
Eigen::GeneralizedSelfAdjointEigenSolver< _MatrixType > Class Template Reference

Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem. More...

#include <GeneralizedSelfAdjointEigenSolver.h>

+ Inheritance diagram for Eigen::GeneralizedSelfAdjointEigenSolver< _MatrixType >:
+ Collaboration diagram for Eigen::GeneralizedSelfAdjointEigenSolver< _MatrixType >:

Public Types

typedef Base::Index Index
 
typedef _MatrixType MatrixType
 
- Public Types inherited from Eigen::SelfAdjointEigenSolver< _MatrixType >
enum  { Size = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 Scalar type for matrices of type _MatrixType.
 
typedef MatrixType::Index Index
 
typedef NumTraits< Scalar >::Real RealScalar
 Real scalar type for _MatrixType. More...
 
typedef internal::plain_col_type< MatrixType, RealScalar >::type RealVectorType
 Type for vector of eigenvalues as returned by eigenvalues(). More...
 
typedef Tridiagonalization< MatrixType > TridiagonalizationType
 

Public Member Functions

 GeneralizedSelfAdjointEigenSolver ()
 Default constructor for fixed-size matrices. More...
 
 GeneralizedSelfAdjointEigenSolver (Index size)
 Constructor, pre-allocates memory for dynamic-size matrices. More...
 
 GeneralizedSelfAdjointEigenSolver (const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
 Constructor; computes generalized eigendecomposition of given matrix pencil. More...
 
GeneralizedSelfAdjointEigenSolvercompute (const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
 Computes generalized eigendecomposition of given matrix pencil. More...
 
- Public Member Functions inherited from Eigen::SelfAdjointEigenSolver< _MatrixType >
 SelfAdjointEigenSolver ()
 Default constructor for fixed-size matrices. More...
 
 SelfAdjointEigenSolver (Index size)
 Constructor, pre-allocates memory for dynamic-size matrices. More...
 
 SelfAdjointEigenSolver (const MatrixType &matrix, int options=ComputeEigenvectors)
 Constructor; computes eigendecomposition of given matrix. More...
 
SelfAdjointEigenSolvercompute (const MatrixType &matrix, int options=ComputeEigenvectors)
 Computes eigendecomposition of given matrix. More...
 
SelfAdjointEigenSolvercomputeDirect (const MatrixType &matrix, int options=ComputeEigenvectors)
 Computes eigendecomposition of given matrix using a direct algorithm. More...
 
const MatrixType & eigenvectors () const
 Returns the eigenvectors of given matrix. More...
 
const RealVectorTypeeigenvalues () const
 Returns the eigenvalues of given matrix. More...
 
MatrixType operatorSqrt () const
 Computes the positive-definite square root of the matrix. More...
 
MatrixType operatorInverseSqrt () const
 Computes the inverse square root of the matrix. More...
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 

Additional Inherited Members

- Static Public Attributes inherited from Eigen::SelfAdjointEigenSolver< _MatrixType >
static const int m_maxIterations = 30
 Maximum number of iterations. More...
 
- Protected Attributes inherited from Eigen::SelfAdjointEigenSolver< _MatrixType >
MatrixType m_eivec
 
RealVectorType m_eivalues
 
TridiagonalizationType::SubDiagonalType m_subdiag
 
ComputationInfo m_info
 
bool m_isInitialized
 
bool m_eigenvectorsOk
 

Detailed Description

template<typename _MatrixType>
class Eigen::GeneralizedSelfAdjointEigenSolver< _MatrixType >

Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.

Template Parameters
_MatrixTypethe type of the matrix of which we are computing the eigendecomposition; this is expected to be an instantiation of the Matrix class template.

This class solves the generalized eigenvalue problem $ Av = \lambda Bv $. In this case, the matrix $ A $ should be selfadjoint and the matrix $ B $ should be positive definite.

Only the lower triangular part of the input matrix is referenced.

Call the function compute() to compute the eigenvalues and eigenvectors of a given matrix. Alternatively, you can use the GeneralizedSelfAdjointEigenSolver(const MatrixType&, const MatrixType&, int) constructor which computes the eigenvalues and eigenvectors at construction time. Once the eigenvalue and eigenvectors are computed, they can be retrieved with the eigenvalues() and eigenvectors() functions.

The documentation for GeneralizedSelfAdjointEigenSolver(const MatrixType&, const MatrixType&, int) contains an example of the typical use of this class.

See also
class SelfAdjointEigenSolver, class EigenSolver, class ComplexEigenSolver

Definition at line 48 of file GeneralizedSelfAdjointEigenSolver.h.

Constructor & Destructor Documentation

template<typename _MatrixType >
Eigen::GeneralizedSelfAdjointEigenSolver< _MatrixType >::GeneralizedSelfAdjointEigenSolver ( )
inline

Default constructor for fixed-size matrices.

The default constructor is useful in cases in which the user intends to perform decompositions via compute(). This constructor can only be used if _MatrixType is a fixed-size matrix; use GeneralizedSelfAdjointEigenSolver(Index) for dynamic-size matrices.

Definition at line 63 of file GeneralizedSelfAdjointEigenSolver.h.

63 : Base() {}
template<typename _MatrixType >
Eigen::GeneralizedSelfAdjointEigenSolver< _MatrixType >::GeneralizedSelfAdjointEigenSolver ( Index  size)
inline

Constructor, pre-allocates memory for dynamic-size matrices.

Parameters
[in]sizePositive integer, size of the matrix whose eigenvalues and eigenvectors will be computed.

This constructor is useful for dynamic-size matrices, when the user intends to perform decompositions via compute(). The size parameter is only used as a hint. It is not an error to give a wrong size, but it may impair performance.

See also
compute() for an example

Definition at line 77 of file GeneralizedSelfAdjointEigenSolver.h.

78  : Base(size)
79  {}
template<typename _MatrixType >
Eigen::GeneralizedSelfAdjointEigenSolver< _MatrixType >::GeneralizedSelfAdjointEigenSolver ( const MatrixType &  matA,
const MatrixType &  matB,
int  options = ComputeEigenvectors|Ax_lBx 
)
inline

Constructor; computes generalized eigendecomposition of given matrix pencil.

Parameters
[in]matASelfadjoint matrix in matrix pencil. Only the lower triangular part of the matrix is referenced.
[in]matBPositive-definite matrix in matrix pencil. Only the lower triangular part of the matrix is referenced.
[in]optionsA or-ed set of flags {#ComputeEigenvectors,#EigenvaluesOnly} | {#Ax_lBx,#ABx_lx,#BAx_lx}. Default is #ComputeEigenvectors|#Ax_lBx.

This constructor calls compute(const MatrixType&, const MatrixType&, int) to compute the eigenvalues and (if requested) the eigenvectors of the generalized eigenproblem $ Ax = \lambda B x $ with matA the selfadjoint matrix $ A $ and matB the positive definite matrix $ B $. Each eigenvector $ x $ satisfies the property $ x^* B x = 1 $. The eigenvectors are computed if options contains ComputeEigenvectors.

In addition, the two following variants can be solved via options:

  • ABx_lx: $ ABx = \lambda x $
  • BAx_lx: $ BAx = \lambda x $

Example:

Output:

See also
compute(const MatrixType&, const MatrixType&, int)

Definition at line 107 of file GeneralizedSelfAdjointEigenSolver.h.

109  : Base(matA.cols())
110  {
111  compute(matA, matB, options);
112  }
GeneralizedSelfAdjointEigenSolver & compute(const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
Computes generalized eigendecomposition of given matrix pencil.

Member Function Documentation

template<typename MatrixType >
GeneralizedSelfAdjointEigenSolver< MatrixType > & Eigen::GeneralizedSelfAdjointEigenSolver< MatrixType >::compute ( const MatrixType &  matA,
const MatrixType &  matB,
int  options = ComputeEigenvectors|Ax_lBx 
)

Computes generalized eigendecomposition of given matrix pencil.

Parameters
[in]matASelfadjoint matrix in matrix pencil. Only the lower triangular part of the matrix is referenced.
[in]matBPositive-definite matrix in matrix pencil. Only the lower triangular part of the matrix is referenced.
[in]optionsA or-ed set of flags {#ComputeEigenvectors,#EigenvaluesOnly} | {#Ax_lBx,#ABx_lx,#BAx_lx}. Default is #ComputeEigenvectors|#Ax_lBx.
Returns
Reference to *this

Accoring to options, this function computes eigenvalues and (if requested) the eigenvectors of one of the following three generalized eigenproblems:

  • Ax_lBx: $ Ax = \lambda B x $
  • ABx_lx: $ ABx = \lambda x $
  • BAx_lx: $ BAx = \lambda x $ with matA the selfadjoint matrix $ A $ and matB the positive definite matrix $ B $. In addition, each eigenvector $ x $ satisfies the property $ x^* B x = 1 $.

The eigenvalues() function can be used to retrieve the eigenvalues. If options contains ComputeEigenvectors, then the eigenvectors are also computed and can be retrieved by calling eigenvectors().

The implementation uses LLT to compute the Cholesky decomposition $ B = LL^* $ and computes the classical eigendecomposition of the selfadjoint matrix $ L^{-1} A (L^*)^{-1} $ if options contains Ax_lBx and of $ L^{*} A L $ otherwise. This solves the generalized eigenproblem, because any solution of the generalized eigenproblem $ Ax = \lambda B x $ corresponds to a solution $ L^{-1} A (L^*)^{-1} (L^* x) = \lambda (L^* x) $ of the eigenproblem for $ L^{-1} A (L^*)^{-1} $. Similar statements can be made for the two other variants.

Example:

Output:

See also
GeneralizedSelfAdjointEigenSolver(const MatrixType&, const MatrixType&, int)

Definition at line 164 of file GeneralizedSelfAdjointEigenSolver.h.

165 {
166  eigen_assert(matA.cols()==matA.rows() && matB.rows()==matA.rows() && matB.cols()==matB.rows());
167  eigen_assert((options&~(EigVecMask|GenEigMask))==0
168  && (options&EigVecMask)!=EigVecMask
169  && ((options&GenEigMask)==0 || (options&GenEigMask)==Ax_lBx
170  || (options&GenEigMask)==ABx_lx || (options&GenEigMask)==BAx_lx)
171  && "invalid option parameter");
172 
173  bool computeEigVecs = ((options&EigVecMask)==0) || ((options&EigVecMask)==ComputeEigenvectors);
174 
175  // Compute the cholesky decomposition of matB = L L' = U'U
176  LLT<MatrixType> cholB(matB);
177 
178  int type = (options&GenEigMask);
179  if(type==0)
180  type = Ax_lBx;
181 
182  if(type==Ax_lBx)
183  {
184  // compute C = inv(L) A inv(L')
185  MatrixType matC = matA.template selfadjointView<Lower>();
186  cholB.matrixL().template solveInPlace<OnTheLeft>(matC);
187  cholB.matrixU().template solveInPlace<OnTheRight>(matC);
188 
189  Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly );
190 
191  // transform back the eigen vectors: evecs = inv(U) * evecs
192  if(computeEigVecs)
193  cholB.matrixU().solveInPlace(Base::m_eivec);
194  }
195  else if(type==ABx_lx)
196  {
197  // compute C = L' A L
198  MatrixType matC = matA.template selfadjointView<Lower>();
199  matC = matC * cholB.matrixL();
200  matC = cholB.matrixU() * matC;
201 
202  Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);
203 
204  // transform back the eigen vectors: evecs = inv(U) * evecs
205  if(computeEigVecs)
206  cholB.matrixU().solveInPlace(Base::m_eivec);
207  }
208  else if(type==BAx_lx)
209  {
210  // compute C = L' A L
211  MatrixType matC = matA.template selfadjointView<Lower>();
212  matC = matC * cholB.matrixL();
213  matC = cholB.matrixU() * matC;
214 
215  Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);
216 
217  // transform back the eigen vectors: evecs = L * evecs
218  if(computeEigVecs)
219  Base::m_eivec = cholB.matrixL() * Base::m_eivec;
220  }
221 
222  return *this;
223 }
SelfAdjointEigenSolver & compute(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.

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