Eigen  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ > Class Template Reference

Detailed Description

template<typename MatrixType_, int UpLo_, typename Ordering_>
class Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >

A direct sparse LLT Cholesky factorizations.

This class provides a LL^T Cholesky factorizations 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
MatrixType_the type of the sparse matrix A, it must be a SparseMatrix<>
UpLo_the triangular part that will be used for the computations. It can be Lower or Upper. Default is Lower.
Ordering_The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<>

This class follows the sparse solver concept .

See also
class SimplicialLDLT, class AMDOrdering, class NaturalOrdering
+ Inheritance diagram for Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >:

Public Member Functions

void analyzePattern (const MatrixType &a)
 
SimplicialLLTcompute (const MatrixType &matrix)
 
Scalar determinant () const
 
void factorize (const MatrixType &a)
 
const MatrixL matrixL () const
 
const MatrixU matrixU () const
 
 SimplicialLLT ()
 
 SimplicialLLT (const MatrixType &matrix)
 
- Public Member Functions inherited from Eigen::SimplicialCholeskyBase< SimplicialLLT< MatrixType_, UpLo_, Ordering_ > >
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationP () const
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationPinv () const
 
SimplicialLLT< MatrixType_, UpLo_, Ordering_ > & setShift (const RealScalar &offset, const RealScalar &scale=1)
 
 SimplicialCholeskyBase ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< SimplicialLLT< MatrixType_, UpLo_, Ordering_ > >
const Solve< SimplicialLLT< MatrixType_, UpLo_, Ordering_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< SimplicialLLT< MatrixType_, UpLo_, Ordering_ >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 

Additional Inherited Members

- Protected Member Functions inherited from Eigen::SimplicialCholeskyBase< SimplicialLLT< MatrixType_, UpLo_, Ordering_ > >
void compute (const MatrixType &matrix)
 

Constructor & Destructor Documentation

◆ SimplicialLLT() [1/2]

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::SimplicialLLT ( )
inline

Default constructor

◆ SimplicialLLT() [2/2]

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::SimplicialLLT ( const MatrixType &  matrix)
inlineexplicit

Constructs and performs the LLT factorization of matrix

Member Function Documentation

◆ analyzePattern()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
void Eigen::SimplicialLLT< 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()

◆ compute()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
SimplicialLLT& Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::compute ( const MatrixType &  matrix)
inline

Computes the sparse Cholesky decomposition of matrix

◆ determinant()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
Scalar Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::determinant ( ) const
inline
Returns
the determinant of the underlying matrix from the current factorization

◆ factorize()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
void Eigen::SimplicialLLT< 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()

◆ matrixL()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
const MatrixL Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::matrixL ( ) const
inline
Returns
an expression of the factor L

◆ matrixU()

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
const MatrixU Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::matrixU ( ) const
inline
Returns
an expression of the factor U (= L^*)

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