template<typename _MatrixType>
GeneralizedSelfAdjointEigenSolver class
Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.
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. |
Contents
This is defined in the Eigenvalues module. #include <Eigen/Eigenvalues>
This class solves the generalized eigenvalue problem . In this case, the matrix should be selfadjoint and the matrix 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.
Base classes
-
template<typename _MatrixType>class SelfAdjointEigenSolver
- Computes eigenvalues and eigenvectors of selfadjoint matrices.
Constructors, destructors, conversion operators
- GeneralizedSelfAdjointEigenSolver()
- Default constructor for fixed-size matrices.
- GeneralizedSelfAdjointEigenSolver(Index size) explicit
- Constructor, pre-allocates memory for dynamic-size matrices.
-
GeneralizedSelfAdjointEigenSolver(const MatrixType& matA,
const MatrixType& matB,
int options = ComputeEigenvectors|Ax_
lBx) - Constructor; computes generalized eigendecomposition of given matrix pencil.
Public functions
-
auto compute(const MatrixType& matA,
const MatrixType& matB,
int options = ComputeEigenvectors|Ax_
lBx) -> GeneralizedSelfAdjointEigenSolver& - Computes generalized eigendecomposition of given matrix pencil.
Function documentation
template<typename _MatrixType>
Eigen:: GeneralizedSelfAdjointEigenSolver<_MatrixType>:: GeneralizedSelfAdjointEigenSolver()
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.
template<typename _MatrixType>
Eigen:: GeneralizedSelfAdjointEigenSolver<_MatrixType>:: GeneralizedSelfAdjointEigenSolver(Index size) explicit
Constructor, pre-allocates memory for dynamic-size matrices.
Parameters | |
---|---|
size in | Positive 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.
template<typename _MatrixType>
Eigen:: GeneralizedSelfAdjointEigenSolver<_MatrixType>:: GeneralizedSelfAdjointEigenSolver(const MatrixType& matA,
const MatrixType& matB,
int options = ComputeEigenvectors|Ax_ lBx)
Constructor; computes generalized eigendecomposition of given matrix pencil.
Parameters | |
---|---|
matA in | Selfadjoint matrix in matrix pencil. Only the lower triangular part of the matrix is referenced. |
matB in | Positive-definite matrix in matrix pencil. Only the lower triangular part of the matrix is referenced. |
options in | A or-ed set of flags {ComputeEigenvectors,EigenvaluesOnly} | {Ax_ |
This constructor calls compute(const MatrixType&, const MatrixType&, int) to compute the eigenvalues and (if requested) the eigenvectors of the generalized eigenproblem with matA the selfadjoint matrix and matB the positive definite matrix . Each eigenvector satisfies the property . The eigenvectors are computed if options contains ComputeEigenvectors.
In addition, the two following variants can be solved via options:
ABx_lx:
BAx_lx:
Example:
MatrixXd X = MatrixXd::Random(5,5); MatrixXd A = X + X.transpose(); cout << "Here is a random symmetric matrix, A:" << endl << A << endl; X = MatrixXd::Random(5,5); MatrixXd B = X * X.transpose(); cout << "and a random postive-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;
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 postive-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
template<typename _MatrixType>
GeneralizedSelfAdjointEigenSolver& Eigen:: GeneralizedSelfAdjointEigenSolver<_MatrixType>:: compute(const MatrixType& matA,
const MatrixType& matB,
int options = ComputeEigenvectors|Ax_ lBx)
Computes generalized eigendecomposition of given matrix pencil.
Parameters | |
---|---|
matA in | Selfadjoint matrix in matrix pencil. Only the lower triangular part of the matrix is referenced. |
matB in | Positive-definite matrix in matrix pencil. Only the lower triangular part of the matrix is referenced. |
options in | A or-ed set of flags {ComputeEigenvectors,EigenvaluesOnly} | {Ax_ |
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:
ABx_lx:
BAx_lx:
with matA the selfadjoint matrix and matB the positive definite matrix . In addition, each eigenvector satisfies the property .
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 and computes the classical eigendecomposition of the selfadjoint matrix if options
contains Ax_lBx and of otherwise. This solves the generalized eigenproblem, because any solution of the generalized eigenproblem corresponds to a solution of the eigenproblem for . Similar statements can be made for the two other variants.
Example:
MatrixXd X = MatrixXd::Random(5,5); MatrixXd A = X * X.transpose(); X = MatrixXd::Random(5,5); 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;
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