Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ > Class Template Reference

A direct sparse LLT Cholesky factorizations. More...

+ Inheritance diagram for Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >:

Public Types

enum  { UpLo }
 
typedef SimplicialCholeskyBase< SimplicialLLTBase
 
typedef SparseMatrix< Scalar, ColMajor, IndexCholMatrixType
 
typedef Traits::MatrixL MatrixL
 
typedef MatrixType_ MatrixType
 
typedef Traits::MatrixU MatrixU
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef internal::traits< SimplicialLLTTraits
 
typedef Matrix< Scalar, Dynamic, 1 > VectorType
 
- Public Types inherited from Eigen::SimplicialCholeskyBase< SimplicialLLT< MatrixType_, UpLo_, Ordering_ > >
enum  
 
enum  
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexCholMatrixType
 
typedef CholMatrixType const * ConstCholMatrixPtr
 
typedef internal::traits< SimplicialLLT< MatrixType_, UpLo_, Ordering_ > >::MatrixType MatrixType
 
typedef internal::traits< SimplicialLLT< MatrixType_, UpLo_, Ordering_ > >::OrderingType OrderingType
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef Matrix< StorageIndex, Dynamic, 1 > VectorI
 
typedef Matrix< Scalar, Dynamic, 1 > VectorType
 

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_ > >
Index cols () const
 
SimplicialLLT< MatrixType_, UpLo_, Ordering_ > & derived ()
 
SimplicialLLT< MatrixType_, UpLo_, Ordering_ > & derived ()
 
const SimplicialLLT< MatrixType_, UpLo_, Ordering_ > & derived () const
 
const SimplicialLLT< MatrixType_, UpLo_, Ordering_ > & derived () const
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationP () const
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationPinv () const
 
Index rows () const
 
SimplicialLLT< MatrixType_, UpLo_, Ordering_ > & setShift (const RealScalar &offset, const RealScalar &scale=1)
 
 SimplicialCholeskyBase ()
 
 SimplicialCholeskyBase (const MatrixType &matrix)
 
 ~SimplicialCholeskyBase ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< SimplicialLLT< MatrixType_, UpLo_, Ordering_ > >
SimplicialLLT< MatrixType_, UpLo_, Ordering_ > & derived ()
 
const SimplicialLLT< MatrixType_, UpLo_, Ordering_ > & derived () const
 
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 ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 

Additional Inherited Members

- Protected Member Functions inherited from Eigen::SimplicialCholeskyBase< SimplicialLLT< MatrixType_, UpLo_, Ordering_ > >
void analyzePattern (const MatrixType &a, bool doLDLT)
 
void analyzePattern_preordered (const CholMatrixType &a, bool doLDLT)
 
void compute (const MatrixType &matrix)
 
void factorize (const MatrixType &a)
 
void factorize_preordered (const CholMatrixType &a)
 
void ordering (const MatrixType &a, ConstCholMatrixPtr &pmat, CholMatrixType &ap)
 
- Protected Attributes inherited from Eigen::SimplicialCholeskyBase< SimplicialLLT< MatrixType_, UpLo_, Ordering_ > >
bool m_analysisIsOk
 
VectorType m_diag
 
bool m_factorizationIsOk
 
ComputationInfo m_info
 
CholMatrixType m_matrix
 
VectorI m_nonZerosPerCol
 
PermutationMatrix< Dynamic, Dynamic, StorageIndexm_P
 
VectorI m_parent
 
PermutationMatrix< Dynamic, Dynamic, StorageIndexm_Pinv
 
RealScalar m_shiftOffset
 
RealScalar m_shiftScale
 
- Protected Attributes inherited from Eigen::SparseSolverBase< SimplicialLLT< MatrixType_, UpLo_, Ordering_ > >
bool m_isInitialized
 

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

Definition at line 340 of file SimplicialCholesky.h.

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef SimplicialCholeskyBase<SimplicialLLT> Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::Base

Definition at line 345 of file SimplicialCholesky.h.

◆ CholMatrixType

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef SparseMatrix<Scalar,ColMajor,Index> Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::CholMatrixType

Definition at line 349 of file SimplicialCholesky.h.

◆ MatrixL

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef Traits::MatrixL Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::MatrixL

Definition at line 352 of file SimplicialCholesky.h.

◆ MatrixType

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

Definition at line 343 of file SimplicialCholesky.h.

◆ MatrixU

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef Traits::MatrixU Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::MatrixU

Definition at line 353 of file SimplicialCholesky.h.

◆ RealScalar

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef MatrixType::RealScalar Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::RealScalar

Definition at line 347 of file SimplicialCholesky.h.

◆ Scalar

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef MatrixType::Scalar Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::Scalar

Definition at line 346 of file SimplicialCholesky.h.

◆ StorageIndex

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef MatrixType::StorageIndex Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::StorageIndex

Definition at line 348 of file SimplicialCholesky.h.

◆ Traits

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef internal::traits<SimplicialLLT> Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::Traits

Definition at line 351 of file SimplicialCholesky.h.

◆ VectorType

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
typedef Matrix<Scalar,Dynamic,1> Eigen::SimplicialLLT< MatrixType_, UpLo_, Ordering_ >::VectorType

Definition at line 350 of file SimplicialCholesky.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , int UpLo_, typename Ordering_ >
anonymous enum
Enumerator
UpLo 

Definition at line 344 of file SimplicialCholesky.h.

344 { UpLo = UpLo_ };

Constructor & Destructor Documentation

◆ SimplicialLLT() [1/2]

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

Default constructor

Definition at line 356 of file SimplicialCholesky.h.

356 : Base() {}
SimplicialCholeskyBase< SimplicialLLT > Base

◆ 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

Definition at line 358 of file SimplicialCholesky.h.

359  : Base(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()

Definition at line 386 of file SimplicialCholesky.h.

387  {
388  Base::analyzePattern(a, false);
389  }
void analyzePattern(const MatrixType &a, bool doLDLT)

◆ 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

Definition at line 374 of file SimplicialCholesky.h.

375  {
376  Base::template compute<false>(matrix);
377  return *this;
378  }

◆ 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

Definition at line 403 of file SimplicialCholesky.h.

404  {
405  Scalar detL = Base::m_matrix.diagonal().prod();
406  return numext::abs2(detL);
407  }
MatrixType::Scalar Scalar
const ConstDiagonalReturnType diagonal() const
Definition: SparseMatrix.h:764
bool abs2(bool x)

◆ 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()

Definition at line 397 of file SimplicialCholesky.h.

398  {
399  Base::template factorize<false>(a);
400  }

◆ 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

Definition at line 362 of file SimplicialCholesky.h.

362  {
363  eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
364  return Traits::getL(Base::m_matrix);
365  }
#define eigen_assert(x)
Definition: Macros.h:902

◆ 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^*)

Definition at line 368 of file SimplicialCholesky.h.

368  {
369  eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
370  return Traits::getU(Base::m_matrix);
371  }

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