Eigen::GeneralizedSelfAdjointEigenSolver< MatrixType_ > Class Template Reference

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

+ Inheritance diagram for Eigen::GeneralizedSelfAdjointEigenSolver< MatrixType_ >:

Public Types

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

Public Member Functions

GeneralizedSelfAdjointEigenSolvercompute (const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
 Computes generalized eigendecomposition of given matrix pencil. More...
 
 GeneralizedSelfAdjointEigenSolver ()
 Default constructor for fixed-size matrices. More...
 
 GeneralizedSelfAdjointEigenSolver (const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
 Constructor; computes generalized eigendecomposition of given matrix pencil. More...
 
 GeneralizedSelfAdjointEigenSolver (Index size)
 Constructor, pre-allocates memory for dynamic-size matrices. More...
 
- Public Member Functions inherited from Eigen::SelfAdjointEigenSolver< MatrixType_ >
template<typename InputType >
SelfAdjointEigenSolver< MatrixType > & compute (const EigenBase< InputType > &a_matrix, int options)
 
template<typename InputType >
SelfAdjointEigenSolvercompute (const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
 Computes eigendecomposition of given matrix. More...
 
SelfAdjointEigenSolvercomputeDirect (const MatrixType &matrix, int options=ComputeEigenvectors)
 Computes eigendecomposition of given matrix using a closed-form algorithm. More...
 
SelfAdjointEigenSolvercomputeFromTridiagonal (const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
 Computes the eigen decomposition from a tridiagonal symmetric matrix. More...
 
const RealVectorTypeeigenvalues () const
 Returns the eigenvalues of given matrix. More...
 
const EigenvectorsTypeeigenvectors () const
 Returns the eigenvectors of given matrix. More...
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
MatrixType operatorInverseSqrt () const
 Computes the inverse square root of the matrix. More...
 
MatrixType operatorSqrt () const
 Computes the positive-definite square root of the matrix. More...
 
 SelfAdjointEigenSolver ()
 Default constructor for fixed-size matrices. More...
 
template<typename InputType >
 SelfAdjointEigenSolver (const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
 Constructor; computes eigendecomposition of given matrix. More...
 
 SelfAdjointEigenSolver (Index size)
 Constructor, pre-allocates memory for dynamic-size matrices. More...
 

Private Types

typedef SelfAdjointEigenSolver< MatrixType_ > Base
 

Additional Inherited Members

- Static Public Attributes inherited from Eigen::SelfAdjointEigenSolver< MatrixType_ >
static const int m_maxIterations
 Maximum number of iterations. More...
 
- Protected Attributes inherited from Eigen::SelfAdjointEigenSolver< MatrixType_ >
bool m_eigenvectorsOk
 
RealVectorType m_eivalues
 
EigenvectorsType m_eivec
 
TridiagonalizationType::CoeffVectorType m_hcoeffs
 
ComputationInfo m_info
 
bool m_isInitialized
 
TridiagonalizationType::SubDiagonalType m_subdiag
 
VectorType m_workspace
 

Detailed Description

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

Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.

This is defined in the Eigenvalues module.

Template Parameters
MatrixType_the 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 50 of file GeneralizedSelfAdjointEigenSolver.h.

Member Typedef Documentation

◆ Base

template<typename MatrixType_ >
typedef SelfAdjointEigenSolver<MatrixType_> Eigen::GeneralizedSelfAdjointEigenSolver< MatrixType_ >::Base
private

Definition at line 52 of file GeneralizedSelfAdjointEigenSolver.h.

◆ MatrixType

template<typename MatrixType_ >
typedef MatrixType_ Eigen::GeneralizedSelfAdjointEigenSolver< MatrixType_ >::MatrixType

Definition at line 55 of file GeneralizedSelfAdjointEigenSolver.h.

Constructor & Destructor Documentation

◆ GeneralizedSelfAdjointEigenSolver() [1/3]

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 64 of file GeneralizedSelfAdjointEigenSolver.h.

64 : Base() {}

◆ GeneralizedSelfAdjointEigenSolver() [2/3]

template<typename MatrixType_ >
Eigen::GeneralizedSelfAdjointEigenSolver< MatrixType_ >::GeneralizedSelfAdjointEigenSolver ( Index  size)
inlineexplicit

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 78 of file GeneralizedSelfAdjointEigenSolver.h.

79  : Base(size)
80  {}

◆ GeneralizedSelfAdjointEigenSolver() [3/3]

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:

MatrixXd A = X + X.transpose();
cout << "Here is a random symmetric matrix, A:" << endl << A << endl;
MatrixXd B = X * X.transpose();
cout << "and a random positive-definite matrix, B:" << endl << B << endl << endl;
GeneralizedSelfAdjointEigenSolver<MatrixXd> es(A,B);
cout << "The eigenvalues of the pencil (A,B) are:" << endl << es.eigenvalues() << endl;
cout << "The matrix of eigenvectors, V, is:" << endl << es.eigenvectors() << endl << endl;
double lambda = es.eigenvalues()[0];
cout << "Consider the first eigenvalue, lambda = " << lambda << endl;
VectorXd v = es.eigenvectors().col(0);
cout << "If v is the corresponding eigenvector, then A * v = " << endl << A * v << endl;
cout << "... and lambda * B * v = " << endl << lambda * B * v << endl << endl;
Array< int, Dynamic, 1 > v
MatrixXcf A
cout<< "The eigenvalues of A are:"<< endl<< ces.eigenvalues()<< endl;cout<< "The matrix of eigenvectors, V, is:"<< endl<< ces.eigenvectors()<< endl<< endl;complex< float > lambda
EigenSolver< MatrixXf > es
MatrixXf B
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< double, Dynamic, 1 > VectorXd
DynamicĂ—1 vector of type double.
Definition: Matrix.h:502
Matrix< double, Dynamic, Dynamic > MatrixXd
DynamicĂ—Dynamic matrix of type double.
Definition: Matrix.h:502

Output:

Here is a random symmetric matrix, A:
  1.36 -0.816  0.521   1.43 -0.144
-0.816 -0.659  0.794 -0.173 -0.406
 0.521  0.794 -0.541  0.461  0.179
  1.43 -0.173  0.461  -1.43  0.822
-0.144 -0.406  0.179  0.822  -1.37
and a random positive-definite matrix, B:
  0.132  0.0109 -0.0512  0.0674  -0.143
 0.0109    1.68    1.13   -1.12   0.916
-0.0512    1.13     2.3   -2.14    1.86
 0.0674   -1.12   -2.14    2.69   -2.01
 -0.143   0.916    1.86   -2.01    1.68

The eigenvalues of the pencil (A,B) are:
  -227
  -3.9
-0.837
 0.101
  54.2
The matrix of eigenvectors, V, is:
   14.2   -1.03  0.0766 -0.0273   -8.36
 0.0546  -0.115   0.729   0.478   0.374
  -9.23   0.624 -0.0165   0.499    3.01
   7.88     1.3   0.225   0.109   -3.85
   20.8   0.805  -0.567 -0.0828   -8.73

Consider the first eigenvalue, lambda = -227
If v is the corresponding eigenvector, then A * v = 
 22.8
-28.8
 19.8
 21.9
-25.9
... and lambda * B * v = 
 22.8
-28.8
 19.8
 21.9
-25.9

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

Definition at line 108 of file GeneralizedSelfAdjointEigenSolver.h.

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

Member Function Documentation

◆ compute()

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

According 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:

MatrixXd A = X * X.transpose();
MatrixXd B = X * X.transpose();
GeneralizedSelfAdjointEigenSolver<MatrixXd> es(A,B,EigenvaluesOnly);
cout << "The eigenvalues of the pencil (A,B) are:" << endl << es.eigenvalues() << endl;
es.compute(B,A,false);
cout << "The eigenvalues of the pencil (B,A) are:" << endl << es.eigenvalues() << endl;
@ EigenvaluesOnly
Definition: Constants.h:404

Output:

The eigenvalues of the pencil (A,B) are:
  0.0289
   0.299
    2.11
    8.64
2.08e+03
The eigenvalues of the pencil (B,A) are:
0.000481
   0.116
   0.473
    3.34
    34.6
See also
GeneralizedSelfAdjointEigenSolver(const MatrixType&, const MatrixType&, int)

Definition at line 164 of file GeneralizedSelfAdjointEigenSolver.h.

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

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