Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::GeneralizedSelfAdjointEigenSolver< _MatrixType > Class Template Reference

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

#include <src/eigen/Eigen/src/Eigenvalues/GeneralizedSelfAdjointEigenSolver.h>

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

Public Types

typedef _MatrixType MatrixType
 
enum  { Size = MatrixType::RowsAtCompileTime , ColsAtCompileTime = MatrixType::ColsAtCompileTime , Options = MatrixType::Options , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef MatrixType::Scalar Scalar
 Scalar type for matrices of type _MatrixType.
 
typedef Eigen::Index Index
 
typedef Matrix< Scalar, Size, Size, ColMajor, MaxColsAtCompileTime, MaxColsAtCompileTimeEigenvectorsType
 
typedef NumTraits< Scalar >::Real RealScalar
 Real scalar type for _MatrixType.
 
typedef internal::plain_col_type< MatrixType, RealScalar >::type RealVectorType
 Type for vector of eigenvalues as returned by eigenvalues().
 
typedef Tridiagonalization< MatrixTypeTridiagonalizationType
 
typedef TridiagonalizationType::SubDiagonalType SubDiagonalType
 

Public Member Functions

 GeneralizedSelfAdjointEigenSolver ()
 Default constructor for fixed-size matrices.
 
 GeneralizedSelfAdjointEigenSolver (Index size)
 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.
 
GeneralizedSelfAdjointEigenSolvercompute (const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
 Computes generalized eigendecomposition of given matrix pencil.
 
template<typename InputType >
EIGEN_DEVICE_FUNC SelfAdjointEigenSolvercompute (const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
 Computes eigendecomposition of given matrix.
 
EIGEN_DEVICE_FUNC SelfAdjointEigenSolvercomputeDirect (const MatrixType &matrix, int options=ComputeEigenvectors)
 Computes eigendecomposition of given matrix using a closed-form algorithm.
 
SelfAdjointEigenSolvercomputeFromTridiagonal (const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
 Computes the eigen decomposition from a tridiagonal symmetric matrix.
 
EIGEN_DEVICE_FUNC const EigenvectorsTypeeigenvectors () const
 Returns the eigenvectors of given matrix.
 
EIGEN_DEVICE_FUNC const RealVectorTypeeigenvalues () const
 Returns the eigenvalues of given matrix.
 
EIGEN_DEVICE_FUNC MatrixType operatorSqrt () const
 Computes the positive-definite square root of the matrix.
 
EIGEN_DEVICE_FUNC MatrixType operatorInverseSqrt () const
 Computes the inverse square root of the matrix.
 
EIGEN_DEVICE_FUNC ComputationInfo info () const
 Reports whether previous computation was successful.
 

Static Public Attributes

static const int m_maxIterations = 30
 Maximum number of iterations.
 

Static Protected Member Functions

static void check_template_parameters ()
 

Protected Attributes

EigenvectorsType m_eivec
 
RealVectorType m_eivalues
 
TridiagonalizationType::SubDiagonalType m_subdiag
 
ComputationInfo m_info
 
bool m_isInitialized
 
bool m_eigenvectorsOk
 

Private Types

typedef SelfAdjointEigenSolver< _MatrixType > Base
 

Detailed Description

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

Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.

\eigenvalues_module

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

Member Typedef Documentation

◆ Base

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

◆ EigenvectorsType

template<typename _MatrixType >
typedef Matrix<Scalar,Size,Size,ColMajor,MaxColsAtCompileTime,MaxColsAtCompileTime> Eigen::SelfAdjointEigenSolver< _MatrixType >::EigenvectorsType
inherited

◆ Index

template<typename _MatrixType >
typedef Eigen::Index Eigen::SelfAdjointEigenSolver< _MatrixType >::Index
inherited

◆ MatrixType

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

◆ RealScalar

template<typename _MatrixType >
typedef NumTraits<Scalar>::Real Eigen::SelfAdjointEigenSolver< _MatrixType >::RealScalar
inherited

Real scalar type for _MatrixType.

This is just Scalar if Scalar is real (e.g., float or double), and the type of the real part of Scalar if Scalar is complex.

◆ RealVectorType

template<typename _MatrixType >
typedef internal::plain_col_type<MatrixType,RealScalar>::type Eigen::SelfAdjointEigenSolver< _MatrixType >::RealVectorType
inherited

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 _MatrixType.

◆ Scalar

template<typename _MatrixType >
typedef MatrixType::Scalar Eigen::SelfAdjointEigenSolver< _MatrixType >::Scalar
inherited

Scalar type for matrices of type _MatrixType.

◆ SubDiagonalType

template<typename _MatrixType >
typedef TridiagonalizationType::SubDiagonalType Eigen::SelfAdjointEigenSolver< _MatrixType >::SubDiagonalType
inherited

◆ TridiagonalizationType

template<typename _MatrixType >
typedef Tridiagonalization<MatrixType> Eigen::SelfAdjointEigenSolver< _MatrixType >::TridiagonalizationType
inherited

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType >
anonymous enum
inherited
Enumerator
Size 
ColsAtCompileTime 
Options 
MaxColsAtCompileTime 
75 {
76 Size = MatrixType::RowsAtCompileTime,
77 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
78 Options = MatrixType::Options,
79 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
80 };
@ Size
Definition SelfAdjointEigenSolver.h:76
@ Options
Definition SelfAdjointEigenSolver.h:78
@ ColsAtCompileTime
Definition SelfAdjointEigenSolver.h:77
@ MaxColsAtCompileTime
Definition SelfAdjointEigenSolver.h:79

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.

62: Base() {}
SelfAdjointEigenSolver< _MatrixType > Base
Definition GeneralizedSelfAdjointEigenSolver.h:50

◆ 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
77 : Base(size)
78 {}

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

Output:

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

References Eigen::GeneralizedSelfAdjointEigenSolver< _MatrixType >::compute().

+ Here is the call graph for this function:

Member Function Documentation

◆ check_template_parameters()

template<typename _MatrixType >
static void Eigen::SelfAdjointEigenSolver< _MatrixType >::check_template_parameters ( )
inlinestaticprotectedinherited
358 {
360 }
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition StaticAssert.h:184
MatrixType::Scalar Scalar
Scalar type for matrices of type _MatrixType.
Definition SelfAdjointEigenSolver.h:83

References EIGEN_STATIC_ASSERT_NON_INTEGER.

◆ compute() [1/2]

template<typename _MatrixType >
template<typename InputType >
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver & Eigen::SelfAdjointEigenSolver< _MatrixType >::compute ( const EigenBase< InputType > &  matrix,
int  options = ComputeEigenvectors 
)
inherited

Computes eigendecomposition of given matrix.

Parameters
[in]matrixSelfadjoint matrix whose eigendecomposition is to be computed. Only the lower triangular part of the matrix is referenced.
[in]optionsCan be ComputeEigenvectors (default) or EigenvaluesOnly.
Returns
Reference to *this

This function computes the eigenvalues of matrix. 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().

This implementation uses a symmetric QR algorithm. The matrix is first reduced to tridiagonal form using the Tridiagonalization class. The tridiagonal matrix is then brought to diagonal form with implicit symmetric QR steps with Wilkinson shift. Details can be found in Section 8.3 of Golub & Van Loan, Matrix Computations.

The cost of the computation is about $ 9n^3 $ if the eigenvectors are required and $ 4n^3/3 $ if they are not required.

This method reuses the memory in the SelfAdjointEigenSolver object that was allocated when the object was constructed, if the size of the matrix does not change.

Example:

Output:

See also
SelfAdjointEigenSolver(const MatrixType&, int)

Referenced by Eigen::SelfAdjointEigenSolver< _MatrixType >::SelfAdjointEigenSolver().

+ Here is the caller graph for this function:

◆ compute() [2/2]

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

References Eigen::ABx_lx, Eigen::Ax_lBx, Eigen::BAx_lx, Eigen::ComputeEigenvectors, eigen_assert, Eigen::EigenvaluesOnly, Eigen::EigVecMask, Eigen::GenEigMask, Eigen::LLT< _MatrixType, _UpLo >::matrixL(), and Eigen::LLT< _MatrixType, _UpLo >::matrixU().

Referenced by Eigen::GeneralizedSelfAdjointEigenSolver< _MatrixType >::GeneralizedSelfAdjointEigenSolver().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ computeDirect()

template<typename MatrixType >
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver< MatrixType > & Eigen::SelfAdjointEigenSolver< MatrixType >::computeDirect ( const MatrixType matrix,
int  options = ComputeEigenvectors 
)
inherited

Computes eigendecomposition of given matrix using a closed-form algorithm.

This is a variant of compute(const MatrixType&, int options) which directly solves the underlying polynomial equation.

Currently only 2x2 and 3x3 matrices for which the sizes are known at compile time are supported (e.g., Matrix3d).

This method is usually significantly faster than the QR iterative algorithm but it might also be less accurate. It is also worth noting that for 3x3 matrices it involves trigonometric operations which are not necessarily available for all scalar types.

For the 3x3 case, we observed the following worst case relative error regarding the eigenvalues:

  • double: 1e-8
  • float: 1e-3
See also
compute(const MatrixType&, int options)
800{
801 internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>::run(*this,matrix,options);
802 return *this;
803}

Referenced by igl::polar_dec().

+ Here is the caller graph for this function:

◆ computeFromTridiagonal()

template<typename MatrixType >
SelfAdjointEigenSolver< MatrixType > & Eigen::SelfAdjointEigenSolver< MatrixType >::computeFromTridiagonal ( const RealVectorType diag,
const SubDiagonalType subdiag,
int  options = ComputeEigenvectors 
)
inherited

Computes the eigen decomposition from a tridiagonal symmetric matrix.

Parameters
[in]diagThe vector containing the diagonal of the matrix.
[in]subdiagThe subdiagonal of the matrix.
[in]optionsCan be ComputeEigenvectors (default) or EigenvaluesOnly.
Returns
Reference to *this

This function assumes that the matrix has been reduced to tridiagonal form.

See also
compute(const MatrixType&, int) for more information
452{
453 //TODO : Add an option to scale the values beforehand
454 bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
455
457 m_subdiag = subdiag;
458 if (computeEigenvectors)
459 {
460 m_eivec.setIdentity(diag.size(), diag.size());
461 }
463
464 m_isInitialized = true;
465 m_eigenvectorsOk = computeEigenvectors;
466 return *this;
467}
RealVectorType m_eivalues
Definition SelfAdjointEigenSolver.h:363
bool m_eigenvectorsOk
Definition SelfAdjointEigenSolver.h:367
TridiagonalizationType::SubDiagonalType m_subdiag
Definition SelfAdjointEigenSolver.h:364
bool m_isInitialized
Definition SelfAdjointEigenSolver.h:366
ComputationInfo m_info
Definition SelfAdjointEigenSolver.h:365
static const int m_maxIterations
Maximum number of iterations.
Definition SelfAdjointEigenSolver.h:354
ComputationInfo computeFromTridiagonal_impl(DiagType &diag, SubDiagType &subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType &eivec)
Definition SelfAdjointEigenSolver.h:482
IGL_INLINE void diag(const Eigen::SparseMatrix< T > &X, Eigen::SparseVector< T > &V)
Definition diag.cpp:17

References Eigen::ComputeEigenvectors, and Eigen::internal::computeFromTridiagonal_impl().

+ Here is the call graph for this function:

◆ eigenvalues()

template<typename _MatrixType >
EIGEN_DEVICE_FUNC const RealVectorType & Eigen::SelfAdjointEigenSolver< _MatrixType >::eigenvalues ( ) const
inlineinherited

Returns the eigenvalues of given matrix.

Returns
A const reference to the column vector containing the eigenvalues.
Precondition
The eigenvalues have been computed before.

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:

See also
eigenvectors(), MatrixBase::eigenvalues()
283 {
284 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
285 return m_eivalues;
286 }

References eigen_assert, Eigen::SelfAdjointEigenSolver< _MatrixType >::m_eivalues, and Eigen::SelfAdjointEigenSolver< _MatrixType >::m_isInitialized.

Referenced by Eigen::SelfAdjointView< _MatrixType, UpLo >::eigenvalues(), CurvatureCalculator::finalEigenStuff(), and igl::polar_dec().

+ Here is the caller graph for this function:

◆ eigenvectors()

template<typename _MatrixType >
EIGEN_DEVICE_FUNC const EigenvectorsType & Eigen::SelfAdjointEigenSolver< _MatrixType >::eigenvectors ( ) const
inlineinherited

Returns the eigenvectors of given matrix.

Returns
A const reference to the matrix whose columns are the eigenvectors.
Precondition
The eigenvectors have been computed before.

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^{-1} $.

Example:

Output:

See also
eigenvalues()
260 {
261 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
262 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
263 return m_eivec;
264 }

References eigen_assert, Eigen::SelfAdjointEigenSolver< _MatrixType >::m_eigenvectorsOk, Eigen::SelfAdjointEigenSolver< _MatrixType >::m_eivec, and Eigen::SelfAdjointEigenSolver< _MatrixType >::m_isInitialized.

Referenced by CurvatureCalculator::finalEigenStuff(), igl::fit_plane(), and igl::polar_dec().

+ Here is the caller graph for this function:

◆ info()

template<typename _MatrixType >
EIGEN_DEVICE_FUNC ComputationInfo Eigen::SelfAdjointEigenSolver< _MatrixType >::info ( ) const
inlineinherited

Reports whether previous computation was successful.

Returns
Success if computation was succesful, NoConvergence otherwise.
344 {
345 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
346 return m_info;
347 }

References eigen_assert, Eigen::SelfAdjointEigenSolver< _MatrixType >::m_info, and Eigen::SelfAdjointEigenSolver< _MatrixType >::m_isInitialized.

◆ operatorInverseSqrt()

template<typename _MatrixType >
EIGEN_DEVICE_FUNC MatrixType Eigen::SelfAdjointEigenSolver< _MatrixType >::operatorInverseSqrt ( ) const
inlineinherited

Computes the inverse square root of the matrix.

Returns
the inverse positive-definite square root of the matrix
Precondition
The eigenvalues and eigenvectors of a positive-definite matrix have been computed before.

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:

See also
operatorSqrt(), MatrixBase::inverse(), MatrixFunctions Module
332 {
333 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
334 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
335 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
336 }

References eigen_assert, Eigen::SelfAdjointEigenSolver< _MatrixType >::m_eigenvectorsOk, Eigen::SelfAdjointEigenSolver< _MatrixType >::m_eivalues, Eigen::SelfAdjointEigenSolver< _MatrixType >::m_eivec, and Eigen::SelfAdjointEigenSolver< _MatrixType >::m_isInitialized.

◆ operatorSqrt()

template<typename _MatrixType >
EIGEN_DEVICE_FUNC MatrixType Eigen::SelfAdjointEigenSolver< _MatrixType >::operatorSqrt ( ) const
inlineinherited

Computes the positive-definite square root of the matrix.

Returns
the positive-definite square root of the matrix
Precondition
The eigenvalues and eigenvectors of a positive-definite matrix have been computed before.

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:

See also
operatorInverseSqrt(), MatrixFunctions Module
307 {
308 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
309 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
310 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
311 }

References eigen_assert, Eigen::SelfAdjointEigenSolver< _MatrixType >::m_eigenvectorsOk, Eigen::SelfAdjointEigenSolver< _MatrixType >::m_eivalues, Eigen::SelfAdjointEigenSolver< _MatrixType >::m_eivec, and Eigen::SelfAdjointEigenSolver< _MatrixType >::m_isInitialized.

Member Data Documentation

◆ m_eigenvectorsOk

◆ m_eivalues

◆ m_eivec

◆ m_info

template<typename _MatrixType >
ComputationInfo Eigen::SelfAdjointEigenSolver< _MatrixType >::m_info
protectedinherited

◆ m_isInitialized

◆ m_maxIterations

template<typename _MatrixType >
const int Eigen::SelfAdjointEigenSolver< _MatrixType >::m_maxIterations = 30
staticinherited

Maximum number of iterations.

The algorithm terminates if it does not converge within m_maxIterations * n iterations, where n denotes the size of the matrix. This value is currently set to 30 (copied from LAPACK).

◆ m_subdiag

template<typename _MatrixType >
TridiagonalizationType::SubDiagonalType Eigen::SelfAdjointEigenSolver< _MatrixType >::m_subdiag
protectedinherited

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