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

A direct sparse LDLT Cholesky factorizations without square root. More...

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

Public Types

enum  { UpLo }
 
typedef SimplicialCholeskyBase< SimplicialLDLTBase
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexCholMatrixType
 
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< SimplicialLDLTTraits
 
typedef Matrix< Scalar, Dynamic, 1 > VectorType
 
- Public Types inherited from Eigen::SimplicialCholeskyBase< SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > >
enum  
 
enum  
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexCholMatrixType
 
typedef CholMatrixType const * ConstCholMatrixPtr
 
typedef internal::traits< SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > >::MatrixType MatrixType
 
typedef internal::traits< SimplicialLDLT< 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)
 
SimplicialLDLTcompute (const MatrixType &matrix)
 
Scalar determinant () const
 
void factorize (const MatrixType &a)
 
const MatrixL matrixL () const
 
const MatrixU matrixU () const
 
 SimplicialLDLT ()
 
 SimplicialLDLT (const MatrixType &matrix)
 
const VectorType vectorD () const
 
- Public Member Functions inherited from Eigen::SimplicialCholeskyBase< SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > >
Index cols () const
 
SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > & derived ()
 
SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > & derived ()
 
const SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > & derived () const
 
const SimplicialLDLT< 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
 
SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > & setShift (const RealScalar &offset, const RealScalar &scale=1)
 
 SimplicialCholeskyBase ()
 
 SimplicialCholeskyBase (const MatrixType &matrix)
 
 ~SimplicialCholeskyBase ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > >
SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > & derived ()
 
const SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > & derived () const
 
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
 
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 

Additional Inherited Members

- Protected Member Functions inherited from Eigen::SimplicialCholeskyBase< SimplicialLDLT< 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< SimplicialLDLT< 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< SimplicialLDLT< MatrixType_, UpLo_, Ordering_ > >
bool m_isInitialized
 

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
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 SimplicialLLT, class AMDOrdering, class NaturalOrdering

Definition at line 431 of file SimplicialCholesky.h.

Member Typedef Documentation

◆ Base

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

Definition at line 436 of file SimplicialCholesky.h.

◆ CholMatrixType

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

Definition at line 440 of file SimplicialCholesky.h.

◆ MatrixL

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

Definition at line 443 of file SimplicialCholesky.h.

◆ MatrixType

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

Definition at line 434 of file SimplicialCholesky.h.

◆ MatrixU

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

Definition at line 444 of file SimplicialCholesky.h.

◆ RealScalar

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

Definition at line 438 of file SimplicialCholesky.h.

◆ Scalar

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

Definition at line 437 of file SimplicialCholesky.h.

◆ StorageIndex

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

Definition at line 439 of file SimplicialCholesky.h.

◆ Traits

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

Definition at line 442 of file SimplicialCholesky.h.

◆ VectorType

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

Definition at line 441 of file SimplicialCholesky.h.

Member Enumeration Documentation

◆ anonymous enum

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

Definition at line 435 of file SimplicialCholesky.h.

435 { UpLo = UpLo_ };

Constructor & Destructor Documentation

◆ SimplicialLDLT() [1/2]

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

Default constructor

Definition at line 447 of file SimplicialCholesky.h.

447 : Base() {}
SimplicialCholeskyBase< SimplicialLDLT > Base

◆ SimplicialLDLT() [2/2]

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

Definition at line 450 of file SimplicialCholesky.h.

451  : Base(matrix) {}

Member Function Documentation

◆ analyzePattern()

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

Definition at line 483 of file SimplicialCholesky.h.

484  {
485  Base::analyzePattern(a, true);
486  }
void analyzePattern(const MatrixType &a, bool doLDLT)

◆ compute()

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

Definition at line 471 of file SimplicialCholesky.h.

472  {
473  Base::template compute<true>(matrix);
474  return *this;
475  }

◆ determinant()

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

Definition at line 500 of file SimplicialCholesky.h.

501  {
502  return Base::m_diag.prod();
503  }

◆ factorize()

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

Definition at line 494 of file SimplicialCholesky.h.

495  {
496  Base::template factorize<true>(a);
497  }

◆ matrixL()

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

Definition at line 459 of file SimplicialCholesky.h.

459  {
460  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
461  return Traits::getL(Base::m_matrix);
462  }
#define eigen_assert(x)
Definition: Macros.h:902

◆ matrixU()

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

Definition at line 465 of file SimplicialCholesky.h.

465  {
466  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
467  return Traits::getU(Base::m_matrix);
468  }

◆ vectorD()

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

Definition at line 454 of file SimplicialCholesky.h.

454  {
455  eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
456  return Base::m_diag;
457  }

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