Eigen  3.2.92
Eigen::SimplicialLDLT< _MatrixType, _UpLo, _Ordering > Class Template Reference

Detailed Description

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

A direct sparse LDLT Cholesky factorizations without square root.

This class provides a LDL^T Cholesky factorizations without square root of sparse matrices that are selfadjoint and positive definite. The factorization allows for solving A.X = B where X and B can be either dense or sparse.

In order to reduce the fill-in, a symmetric permutation P is applied prior to the factorization such that the factorized matrix is P A P^-1.

Template Parameters
_MatrixTypethe type of the sparse matrix A, it must be a SparseMatrix<>
_UpLothe triangular part that will be used for the computations. It can be Lower or Upper. Default is Lower.
_OrderingThe ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<>

This class follows the sparse solver concept .

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

Public Member Functions

void analyzePattern (const MatrixType &a)
 
SimplicialLDLTcompute (const MatrixType &matrix)
 
Scalar determinant () const
 
void factorize (const MatrixType &a)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
const MatrixL matrixL () const
 
const MatrixU matrixU () const
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationP () const
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationPinv () const
 
SimplicialLDLT< _MatrixType, _UpLo, _Ordering > & setShift (const RealScalar &offset, const RealScalar &scale=1)
 
 SimplicialLDLT ()
 
 SimplicialLDLT (const MatrixType &matrix)
 
const Solve< SimplicialLDLT< _MatrixType, _UpLo, _Ordering >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< SimplicialLDLT< _MatrixType, _UpLo, _Ordering >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
const VectorType vectorD () const
 

Constructor & Destructor Documentation

template<typename _MatrixType , int _UpLo, typename _Ordering >
Eigen::SimplicialLDLT< _MatrixType, _UpLo, _Ordering >::SimplicialLDLT ( )
inline

Default constructor

template<typename _MatrixType , int _UpLo, typename _Ordering >
Eigen::SimplicialLDLT< _MatrixType, _UpLo, _Ordering >::SimplicialLDLT ( const MatrixType &  matrix)
inlineexplicit

Constructs and performs the LLT factorization of matrix

Member Function Documentation

template<typename _MatrixType , int _UpLo, typename _Ordering >
void Eigen::SimplicialLDLT< _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 >
SimplicialLDLT& Eigen::SimplicialLDLT< _MatrixType, _UpLo, _Ordering >::compute ( const MatrixType &  matrix)
inline

Computes the sparse Cholesky decomposition of matrix

template<typename _MatrixType , int _UpLo, typename _Ordering >
Scalar Eigen::SimplicialLDLT< _MatrixType, _UpLo, _Ordering >::determinant ( ) const
inline
Returns
the determinant of the underlying matrix from the current factorization
template<typename _MatrixType , int _UpLo, typename _Ordering >
void Eigen::SimplicialLDLT< _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< SimplicialLDLT< _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.
template<typename _MatrixType , int _UpLo, typename _Ordering >
const MatrixL Eigen::SimplicialLDLT< _MatrixType, _UpLo, _Ordering >::matrixL ( ) const
inline
Returns
an expression of the factor L
template<typename _MatrixType , int _UpLo, typename _Ordering >
const MatrixU Eigen::SimplicialLDLT< _MatrixType, _UpLo, _Ordering >::matrixU ( ) const
inline
Returns
an expression of the factor U (= L^*)
const PermutationMatrix<Dynamic,Dynamic,StorageIndex>& Eigen::SimplicialCholeskyBase< SimplicialLDLT< _MatrixType, _UpLo, _Ordering > >::permutationP ( ) const
inlineinherited
Returns
the permutation P
See also
permutationPinv()
const PermutationMatrix<Dynamic,Dynamic,StorageIndex>& Eigen::SimplicialCholeskyBase< SimplicialLDLT< _MatrixType, _UpLo, _Ordering > >::permutationPinv ( ) const
inlineinherited
Returns
the inverse P^-1 of the permutation P
See also
permutationP()
SimplicialLDLT< _MatrixType, _UpLo, _Ordering > & Eigen::SimplicialCholeskyBase< SimplicialLDLT< _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<SimplicialLDLT< _MatrixType, _UpLo, _Ordering > , Rhs> Eigen::SparseSolverBase< SimplicialLDLT< _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<SimplicialLDLT< _MatrixType, _UpLo, _Ordering > , Rhs> Eigen::SparseSolverBase< SimplicialLDLT< _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()
template<typename _MatrixType , int _UpLo, typename _Ordering >
const VectorType Eigen::SimplicialLDLT< _MatrixType, _UpLo, _Ordering >::vectorD ( ) const
inline
Returns
a vector expression of the diagonal D

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