Public Types | |
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 MatrixType::Scalar | Scalar |
Scalar type for matrices of type MatrixType . More... | |
Public Member Functions | |
ArpackGeneralizedSelfAdjointEigenSolver () | |
Default constructor. More... | |
ArpackGeneralizedSelfAdjointEigenSolver (const MatrixType &A, const MatrixType &B, Index nbrEigenvalues, std::string eigs_sigma="LM", int options=ComputeEigenvectors, RealScalar tol=0.0) | |
Constructor; computes generalized eigenvalues of given matrix with respect to another matrix. More... | |
ArpackGeneralizedSelfAdjointEigenSolver (const MatrixType &A, Index nbrEigenvalues, std::string eigs_sigma="LM", int options=ComputeEigenvectors, RealScalar tol=0.0) | |
Constructor; computes eigenvalues of given matrix. More... | |
ArpackGeneralizedSelfAdjointEigenSolver & | compute (const MatrixType &A, const MatrixType &B, Index nbrEigenvalues, std::string eigs_sigma="LM", int options=ComputeEigenvectors, RealScalar tol=0.0) |
Computes generalized eigenvalues / eigenvectors of given matrix using the external ARPACK library. More... | |
ArpackGeneralizedSelfAdjointEigenSolver & | compute (const MatrixType &A, Index nbrEigenvalues, std::string eigs_sigma="LM", int options=ComputeEigenvectors, RealScalar tol=0.0) |
Computes eigenvalues / eigenvectors of given matrix using the external ARPACK library. More... | |
const Matrix< Scalar, Dynamic, 1 > & | eigenvalues () const |
Returns the eigenvalues of given matrix. More... | |
const Matrix< Scalar, Dynamic, Dynamic > & | eigenvectors () const |
Returns the eigenvectors of given matrix. More... | |
size_t | getNbrConvergedEigenValues () const |
size_t | getNbrIterations () const |
ComputationInfo | info () const |
Reports whether previous computation was successful. More... | |
Matrix< Scalar, Dynamic, Dynamic > | operatorInverseSqrt () const |
Computes the inverse square root of the matrix. More... | |
Matrix< Scalar, Dynamic, Dynamic > | operatorSqrt () const |
Computes the positive-definite square root of the matrix. More... | |
Protected Attributes | |
bool | m_eigenvectorsOk |
Matrix< Scalar, Dynamic, 1 > | m_eivalues |
Matrix< Scalar, Dynamic, Dynamic > | m_eivec |
ComputationInfo | m_info |
bool | m_isInitialized |
size_t | m_nbrConverged |
size_t | m_nbrIterations |
Definition at line 27 of file ArpackSelfAdjointEigenSolver.h.
typedef MatrixType::Index Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::Index |
Definition at line 34 of file ArpackSelfAdjointEigenSolver.h.
typedef NumTraits<Scalar>::Real Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::RealScalar |
Real scalar type for MatrixType
.
This is just Scalar
if Scalar is real (e.g., float
or Scalar
), and the type of the real part of Scalar
if Scalar is complex.
Definition at line 42 of file ArpackSelfAdjointEigenSolver.h.
typedef internal::plain_col_type<MatrixType, RealScalar>::type Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::RealVectorType |
Type for vector of eigenvalues as returned by eigenvalues().
This is a column vector with entries of type RealScalar. The length of the vector is the size of nbrEigenvalues
.
Definition at line 49 of file ArpackSelfAdjointEigenSolver.h.
typedef MatrixType::Scalar Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::Scalar |
Scalar type for matrices of type MatrixType
.
Definition at line 33 of file ArpackSelfAdjointEigenSolver.h.
|
inline |
Default constructor.
The default constructor is for cases in which the user intends to perform decompositions via compute().
Definition at line 57 of file ArpackSelfAdjointEigenSolver.h.
|
inline |
Constructor; computes generalized eigenvalues of given matrix with respect to another matrix.
[in] | A | Self-adjoint matrix whose eigenvalues / eigenvectors will computed. By default, the upper triangular part is used, but can be changed through the template parameter. |
[in] | B | Self-adjoint matrix for the generalized eigenvalue problem. |
[in] | nbrEigenvalues | The number of eigenvalues / eigenvectors to compute. Must be less than the size of the input matrix, or an error is returned. |
[in] | eigs_sigma | String containing either "LM", "SM", "LA", or "SA", with respective meanings to find the largest magnitude , smallest magnitude, largest algebraic, or smallest algebraic eigenvalues. Alternatively, this value can contain floating point value in string form, in which case the eigenvalues closest to this value will be found. |
[in] | options | Can be ComputeEigenvectors (default) or EigenvaluesOnly. |
[in] | tol | What tolerance to find the eigenvalues to. Default is 0, which means machine precision. |
This constructor calls compute(const MatrixType&, const MatrixType&, Index, string, int, RealScalar) to compute the eigenvalues of the matrix A
with respect to B
. The eigenvectors are computed if options
equals ComputeEigenvectors.
Definition at line 88 of file ArpackSelfAdjointEigenSolver.h.
|
inline |
Constructor; computes eigenvalues of given matrix.
[in] | A | Self-adjoint matrix whose eigenvalues / eigenvectors will computed. By default, the upper triangular part is used, but can be changed through the template parameter. |
[in] | nbrEigenvalues | The number of eigenvalues / eigenvectors to compute. Must be less than the size of the input matrix, or an error is returned. |
[in] | eigs_sigma | String containing either "LM", "SM", "LA", or "SA", with respective meanings to find the largest magnitude , smallest magnitude, largest algebraic, or smallest algebraic eigenvalues. Alternatively, this value can contain floating point value in string form, in which case the eigenvalues closest to this value will be found. |
[in] | options | Can be ComputeEigenvectors (default) or EigenvaluesOnly. |
[in] | tol | What tolerance to find the eigenvalues to. Default is 0, which means machine precision. |
This constructor calls compute(const MatrixType&, Index, string, int, RealScalar) to compute the eigenvalues of the matrix A
. The eigenvectors are computed if options
equals ComputeEigenvectors.
Definition at line 123 of file ArpackSelfAdjointEigenSolver.h.
ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD > & Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::compute | ( | const MatrixType & | A, |
const MatrixType & | B, | ||
Index | nbrEigenvalues, | ||
std::string | eigs_sigma = "LM" , |
||
int | options = ComputeEigenvectors , |
||
RealScalar | tol = 0.0 |
||
) |
Computes generalized eigenvalues / eigenvectors of given matrix using the external ARPACK library.
[in] | A | Selfadjoint matrix whose eigendecomposition is to be computed. |
[in] | B | Selfadjoint matrix for generalized eigenvalues. |
[in] | nbrEigenvalues | The number of eigenvalues / eigenvectors to compute. Must be less than the size of the input matrix, or an error is returned. |
[in] | eigs_sigma | String containing either "LM", "SM", "LA", or "SA", with respective meanings to find the largest magnitude , smallest magnitude, largest algebraic, or smallest algebraic eigenvalues. Alternatively, this value can contain floating point value in string form, in which case the eigenvalues closest to this value will be found. |
[in] | options | Can be ComputeEigenvectors (default) or EigenvaluesOnly. |
[in] | tol | What tolerance to find the eigenvalues to. Default is 0, which means machine precision. |
*this
This function computes the generalized eigenvalues of A
with respect to B
using ARPACK. The eigenvalues() function can be used to retrieve them. If options
equals ComputeEigenvectors, then the eigenvectors are also computed and can be retrieved by calling eigenvectors().
Definition at line 335 of file ArpackSelfAdjointEigenSolver.h.
ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD > & Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::compute | ( | const MatrixType & | A, |
Index | nbrEigenvalues, | ||
std::string | eigs_sigma = "LM" , |
||
int | options = ComputeEigenvectors , |
||
RealScalar | tol = 0.0 |
||
) |
Computes eigenvalues / eigenvectors of given matrix using the external ARPACK library.
[in] | A | Selfadjoint matrix whose eigendecomposition is to be computed. |
[in] | nbrEigenvalues | The number of eigenvalues / eigenvectors to compute. Must be less than the size of the input matrix, or an error is returned. |
[in] | eigs_sigma | String containing either "LM", "SM", "LA", or "SA", with respective meanings to find the largest magnitude , smallest magnitude, largest algebraic, or smallest algebraic eigenvalues. Alternatively, this value can contain floating point value in string form, in which case the eigenvalues closest to this value will be found. |
[in] | options | Can be ComputeEigenvectors (default) or EigenvaluesOnly. |
[in] | tol | What tolerance to find the eigenvalues to. Default is 0, which means machine precision. |
*this
This function computes the eigenvalues of A
using ARPACK. The eigenvalues() function can be used to retrieve them. If options
equals ComputeEigenvectors, then the eigenvectors are also computed and can be retrieved by calling eigenvectors().
Definition at line 322 of file ArpackSelfAdjointEigenSolver.h.
|
inline |
Returns the eigenvalues of given matrix.
The eigenvalues are repeated according to their algebraic multiplicity, so there are as many eigenvalues as rows in the matrix. The eigenvalues are sorted in increasing order.
Example:
Output:
The eigenvalues of the 3x3 matrix of ones are: -3.09e-16 0 3
Definition at line 232 of file ArpackSelfAdjointEigenSolver.h.
|
inline |
Returns the eigenvectors of given matrix.
Column \( k \) of the returned matrix is an eigenvector corresponding to eigenvalue number \( k \) as returned by eigenvalues(). The eigenvectors are normalized to have (Euclidean) norm equal to one. If this object was used to solve the eigenproblem for the selfadjoint matrix \( A \), then the matrix returned by this function is the matrix \( V \) in the eigendecomposition \( A V = D V \). For the generalized eigenproblem, the matrix returned is the solution \( A V = D B V \)
Example:
Output:
The first eigenvector of the 3x3 matrix of ones is: -0.816 0.408 0.408
Definition at line 210 of file ArpackSelfAdjointEigenSolver.h.
|
inline |
Definition at line 298 of file ArpackSelfAdjointEigenSolver.h.
|
inline |
Definition at line 301 of file ArpackSelfAdjointEigenSolver.h.
|
inline |
Reports whether previous computation was successful.
Success
if computation was successful, NoConvergence
otherwise. Definition at line 292 of file ArpackSelfAdjointEigenSolver.h.
|
inline |
Computes the inverse square root of the matrix.
This function uses the eigendecomposition \( A = V D V^{-1} \) to compute the inverse square root as \( V D^{-1/2} V^{-1} \). This is cheaper than first computing the square root with operatorSqrt() and then its inverse with MatrixBase::inverse().
Example:
Output:
Here is a random positive-definite matrix, A: 1.41 -0.697 -0.111 0.508 -0.697 0.423 0.0991 -0.4 -0.111 0.0991 1.25 0.902 0.508 -0.4 0.902 1.4 The inverse square root of A is: 1.88 2.78 -0.546 0.605 2.78 8.61 -2.3 2.74 -0.546 -2.3 1.92 -1.36 0.605 2.74 -1.36 2.18 We can also compute it with operatorSqrt() and inverse(). That yields: 1.88 2.78 -0.546 0.605 2.78 8.61 -2.3 2.74 -0.546 -2.3 1.92 -1.36 0.605 2.74 -1.36 2.18
Definition at line 281 of file ArpackSelfAdjointEigenSolver.h.
|
inline |
Computes the positive-definite square root of the matrix.
The square root of a positive-definite matrix \( A \) is the positive-definite matrix whose square equals \( A \). This function uses the eigendecomposition \( A = V D V^{-1} \) to compute the square root as \( A^{1/2} = V D^{1/2} V^{-1} \).
Example:
Output:
Here is a random positive-definite matrix, A: 1.41 -0.697 -0.111 0.508 -0.697 0.423 0.0991 -0.4 -0.111 0.0991 1.25 0.902 0.508 -0.4 0.902 1.4 The square root of A is: 1.09 -0.432 -0.0685 0.2 -0.432 0.379 0.141 -0.269 -0.0685 0.141 1 0.468 0.2 -0.269 0.468 1.04 If we square this, we get: 1.41 -0.697 -0.111 0.508 -0.697 0.423 0.0991 -0.4 -0.111 0.0991 1.25 0.902 0.508 -0.4 0.902 1.4
Definition at line 256 of file ArpackSelfAdjointEigenSolver.h.
|
protected |
Definition at line 309 of file ArpackSelfAdjointEigenSolver.h.
|
protected |
Definition at line 306 of file ArpackSelfAdjointEigenSolver.h.
|
protected |
Definition at line 305 of file ArpackSelfAdjointEigenSolver.h.
|
protected |
Definition at line 307 of file ArpackSelfAdjointEigenSolver.h.
|
protected |
Definition at line 308 of file ArpackSelfAdjointEigenSolver.h.
|
protected |
Definition at line 311 of file ArpackSelfAdjointEigenSolver.h.
|
protected |
Definition at line 312 of file ArpackSelfAdjointEigenSolver.h.