Eigen  3.2.91
Eigen::SimplicialCholesky< _MatrixType, _UpLo, _Ordering > Class Template Reference

Detailed Description

template<typename _MatrixType, int _UpLo, typename _Ordering>
class Eigen::SimplicialCholesky< _MatrixType, _UpLo, _Ordering >

Deprecated:

use SimplicialLDLT or class SimplicialLLT

See also
class SimplicialLDLT, class SimplicialLLT
+ Inheritance diagram for Eigen::SimplicialCholesky< _MatrixType, _UpLo, _Ordering >:

Public Member Functions

void analyzePattern (const MatrixType &a)
 
SimplicialCholeskycompute (const MatrixType &matrix)
 
void factorize (const MatrixType &a)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationP () const
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationPinv () const
 
SimplicialCholesky< _MatrixType, _UpLo, _Ordering > & setShift (const RealScalar &offset, const RealScalar &scale=1)
 
const Solve< SimplicialCholesky< _MatrixType, _UpLo, _Ordering >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< SimplicialCholesky< _MatrixType, _UpLo, _Ordering >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 

Member Function Documentation

template<typename _MatrixType, int _UpLo, typename _Ordering>
void Eigen::SimplicialCholesky< _MatrixType, _UpLo, _Ordering >::analyzePattern ( const MatrixType &  a)
inline

Performs a symbolic decomposition on the sparcity of matrix.

This function is particularly useful when solving for several problems having the same structure.

See also
factorize()
template<typename _MatrixType, int _UpLo, typename _Ordering>
SimplicialCholesky& Eigen::SimplicialCholesky< _MatrixType, _UpLo, _Ordering >::compute ( const MatrixType &  matrix)
inline

Computes the sparse Cholesky decomposition of matrix

template<typename _MatrixType, int _UpLo, typename _Ordering>
void Eigen::SimplicialCholesky< _MatrixType, _UpLo, _Ordering >::factorize ( const MatrixType &  a)
inline

Performs a numeric decomposition of matrix

The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.

See also
analyzePattern()
ComputationInfo Eigen::SimplicialCholeskyBase< SimplicialCholesky< _MatrixType, _UpLo, _Ordering > >::info ( ) const
inlineinherited

Reports whether previous computation was successful.

Returns
Success if computation was succesful, NumericalIssue if the matrix.appears to be negative.
const PermutationMatrix<Dynamic,Dynamic,StorageIndex>& Eigen::SimplicialCholeskyBase< SimplicialCholesky< _MatrixType, _UpLo, _Ordering > >::permutationP ( ) const
inlineinherited
Returns
the permutation P
See also
permutationPinv()
const PermutationMatrix<Dynamic,Dynamic,StorageIndex>& Eigen::SimplicialCholeskyBase< SimplicialCholesky< _MatrixType, _UpLo, _Ordering > >::permutationPinv ( ) const
inlineinherited
Returns
the inverse P^-1 of the permutation P
See also
permutationP()
SimplicialCholesky< _MatrixType, _UpLo, _Ordering > & Eigen::SimplicialCholeskyBase< SimplicialCholesky< _MatrixType, _UpLo, _Ordering > >::setShift ( const RealScalar &  offset,
const RealScalar &  scale = 1 
)
inlineinherited

Sets the shift parameters that will be used to adjust the diagonal coefficients during the numerical factorization.

During the numerical factorization, the diagonal coefficients are transformed by the following linear model:
d_ii = offset + scale * d_ii

The default is the identity transformation with offset=0, and scale=1.

Returns
a reference to *this.
const Solve<SimplicialCholesky< _MatrixType, _UpLo, _Ordering > , Rhs> Eigen::SparseSolverBase< SimplicialCholesky< _MatrixType, _UpLo, _Ordering > >::solve ( const MatrixBase< Rhs > &  b) const
inlineinherited
Returns
an expression of the solution x of $ A x = b $ using the current decomposition of A.
See also
compute()
const Solve<SimplicialCholesky< _MatrixType, _UpLo, _Ordering > , Rhs> Eigen::SparseSolverBase< SimplicialCholesky< _MatrixType, _UpLo, _Ordering > >::solve ( const SparseMatrixBase< Rhs > &  b) const
inlineinherited
Returns
an expression of the solution x of $ A x = b $ using the current decomposition of A.
See also
compute()

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