Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ > Class Template Reference

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

#include <CholmodSupport.h>

+ Inheritance diagram for Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >:

Public Types

typedef MatrixType_ MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef TriangularView< const MatrixType, Eigen::LowerMatrixL
 
typedef TriangularView< const typename MatrixType::AdjointReturnType, Eigen::UpperMatrixU
 
- Public Types inherited from Eigen::CholmodBase< MatrixType_, Lower, CholmodSimplicialLLT< MatrixType_, Lower > >
enum  
 
enum  
 
typedef MatrixType_ MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType CholMatrixType
 
typedef MatrixType::StorageIndex StorageIndex
 

Public Member Functions

 CholmodSimplicialLLT ()
 
 CholmodSimplicialLLT (const MatrixType &matrix)
 
 ~CholmodSimplicialLLT ()
 
MatrixL matrixL () const
 
MatrixU matrixU () const
 
- Public Member Functions inherited from Eigen::CholmodBase< MatrixType_, Lower, CholmodSimplicialLLT< MatrixType_, Lower > >
 CholmodBase ()
 
 CholmodBase (const MatrixType &matrix)
 
 ~CholmodBase ()
 
StorageIndex cols () const
 
StorageIndex rows () const
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
CholmodSimplicialLLT< MatrixType_, Lower > & compute (const MatrixType &matrix)
 
void analyzePattern (const MatrixType &matrix)
 
void factorize (const MatrixType &matrix)
 
cholmod_common & cholmod ()
 
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
 
void _solve_impl (const SparseMatrixBase< RhsDerived > &b, SparseMatrixBase< DestDerived > &dest) const
 
CholmodSimplicialLLT< MatrixType_, Lower > & setShift (const RealScalar &offset)
 
Scalar determinant () const
 
Scalar logDeterminant () const
 
void dumpMemory (Stream &)
 
- Public Member Functions inherited from Eigen::SparseSolverBase< Derived >
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Protected Member Functions

void init ()
 
- Protected Member Functions inherited from Eigen::CholmodBase< MatrixType_, Lower, CholmodSimplicialLLT< MatrixType_, Lower > >
CholmodSimplicialLLT< MatrixType_, Lower > & derived ()
 
const CholmodSimplicialLLT< MatrixType_, Lower > & derived () const
 

Private Types

typedef CholmodBase< MatrixType_, UpLo_, CholmodSimplicialLLTBase
 

Private Attributes

cholmod_common m_cholmod
 

Additional Inherited Members

- Protected Types inherited from Eigen::CholmodBase< MatrixType_, Lower, CholmodSimplicialLLT< MatrixType_, Lower > >
typedef SparseSolverBase< CholmodSimplicialLLT< MatrixType_, Lower > > Base
 
- Protected Attributes inherited from Eigen::CholmodBase< MatrixType_, Lower, CholmodSimplicialLLT< MatrixType_, Lower > >
cholmod_common m_cholmod
 
cholmod_factor * m_cholmodFactor
 
double m_shiftOffset [2]
 
ComputationInfo m_info
 
int m_factorizationIsOk
 
int m_analysisIsOk
 
bool m_isInitialized
 
- Protected Attributes inherited from Eigen::SparseSolverBase< Derived >
bool m_isInitialized
 

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.

\implsparsesolverconcept

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

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , int UpLo_ = Lower>
typedef CholmodBase<MatrixType_, UpLo_, CholmodSimplicialLLT> Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::Base
private

◆ MatrixL

template<typename MatrixType_ , int UpLo_ = Lower>
typedef TriangularView<const MatrixType, Eigen::Lower> Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::MatrixL

◆ MatrixType

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

◆ MatrixU

template<typename MatrixType_ , int UpLo_ = Lower>
typedef TriangularView<const typename MatrixType::AdjointReturnType, Eigen::Upper> Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::MatrixU

◆ RealScalar

template<typename MatrixType_ , int UpLo_ = Lower>
typedef MatrixType::RealScalar Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::RealScalar

◆ Scalar

template<typename MatrixType_ , int UpLo_ = Lower>
typedef MatrixType::Scalar Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::Scalar

◆ StorageIndex

template<typename MatrixType_ , int UpLo_ = Lower>
typedef MatrixType::StorageIndex Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::StorageIndex

Constructor & Destructor Documentation

◆ CholmodSimplicialLLT() [1/2]

template<typename MatrixType_ , int UpLo_ = Lower>
Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::CholmodSimplicialLLT ( )
inline
504 : Base() { init(); }
void init()
Definition: CholmodSupport.h:520
CholmodBase< MatrixType_, UpLo_, CholmodSimplicialLLT > Base
Definition: CholmodSupport.h:493

References Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::init().

◆ CholmodSimplicialLLT() [2/2]

template<typename MatrixType_ , int UpLo_ = Lower>
Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::CholmodSimplicialLLT ( const MatrixType matrix)
inline
506  : Base() {
507  init();
508  this->compute(matrix);
509  }
CholmodSimplicialLLT< MatrixType_, Lower > & compute(const MatrixType &matrix)
Definition: CholmodSupport.h:294
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85

References Eigen::CholmodBase< MatrixType_, Lower, CholmodSimplicialLLT< MatrixType_, Lower > >::compute(), and Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::init().

◆ ~CholmodSimplicialLLT()

template<typename MatrixType_ , int UpLo_ = Lower>
Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::~CholmodSimplicialLLT ( )
inline
511 {}

Member Function Documentation

◆ init()

template<typename MatrixType_ , int UpLo_ = Lower>
void Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::init ( )
inlineprotected
520  {
521  m_cholmod.final_asis = 0;
522  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
523  m_cholmod.final_ll = 1;
524  }
cholmod_common m_cholmod
Definition: CholmodSupport.h:460

References Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::m_cholmod.

Referenced by Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::CholmodSimplicialLLT().

◆ matrixL()

template<typename MatrixType_ , int UpLo_ = Lower>
MatrixL Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::matrixL ( ) const
inline
Returns
an expression of the factor L
514 { return viewAsEigen<Scalar, StorageIndex>(*Base::m_cholmodFactor); }
cholmod_factor * m_cholmodFactor
Definition: CholmodSupport.h:461

References Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmodFactor.

Referenced by Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::matrixU().

◆ matrixU()

template<typename MatrixType_ , int UpLo_ = Lower>
MatrixU Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::matrixU ( ) const
inline
Returns
an expression of the factor U (= L^*)
517 { return matrixL().adjoint(); }
MatrixL matrixL() const
Definition: CholmodSupport.h:514
EIGEN_DEVICE_FUNC const AdjointReturnType adjoint() const
Definition: TriangularMatrix.h:224

References Eigen::TriangularView< MatrixType_, Mode_ >::adjoint(), and Eigen::CholmodSimplicialLLT< MatrixType_, UpLo_ >::matrixL().

Member Data Documentation

◆ m_cholmod

template<typename MatrixType_ , int UpLo_ = Lower>
cholmod_common Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmod
mutableprivate

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