Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ > Class Template Reference

Detailed Description

template<typename MatrixType_, int UpLo_ = Lower>
class Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >

A simplicial direct Cholesky (LLT) factorization and solver based on Cholmod.

This class allows to solve for A.X = B sparse linear problems via a simplicial LL^T Cholesky factorization using the Cholmod library. This simplicial variant is equivalent to Eigen's built-in SimplicialLLT class. Therefore, it has little practical interest. The sparse matrix A must be selfadjoint and positive definite. The vectors or matrices X and B can be either dense or sparse.

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.

This class follows the sparse solver concept .

This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.

Warning
Only double precision real and complex scalar types are supported by Cholmod.
See also
Sparse solver concept, class CholmodSupernodalLLT, class SimplicialLLT
+ Inheritance diagram for Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >:

Public Member Functions

MatrixL matrixL () const
 
MatrixU matrixU () const
 
- Public Member Functions inherited from Eigen::CholmodBase< MatrixType_, Lower, CholmodSimplicialLLT< MatrixType_, Lower > >
void analyzePattern (const MatrixType &matrix)
 
cholmod_common & cholmod ()
 
CholmodSimplicialLLT< MatrixType_, Lower > & compute (const MatrixType &matrix)
 
Scalar determinant () const
 
void factorize (const MatrixType &matrix)
 
ComputationInfo info () const
 Reports whether previous computation was successful.
 
Scalar logDeterminant () const
 
CholmodSimplicialLLT< MatrixType_, Lower > & setShift (const RealScalar &offset)
 
- Public Member Functions inherited from Eigen::SparseSolverBase< CholmodSimplicialLLT< MatrixType_, Lower > >
const Solve< CholmodSimplicialLLT< MatrixType_, Lower >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< CholmodSimplicialLLT< MatrixType_, Lower >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 

Member Function Documentation

◆ matrixL()

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

◆ matrixU()

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

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