template<typename _MatrixType, int _UpLo, typename _Ordering>
Eigen::SimplicialLDLT class

A direct sparse LDLT Cholesky factorizations without square root.

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 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.

This class follows the sparse solver concept.

Base classes

template<typename Derived>
class SimplicialCholeskyBase
A base class for direct sparse Cholesky factorizations.

Constructors, destructors, conversion operators

SimplicialLDLT()
SimplicialLDLT(const MatrixType& matrix) explicit

Public functions

void analyzePattern(const MatrixType& a)
auto compute(const MatrixType& matrix) -> SimplicialLDLT&
auto determinant() const -> Scalar
void factorize(const MatrixType& a)
auto matrixL() const -> const MatrixL
auto matrixU() const -> const MatrixU
auto vectorD() const -> const VectorType

Function documentation

template<typename _MatrixType, int _UpLo, typename _Ordering>
Eigen::SimplicialLDLT<_MatrixType, _UpLo, _Ordering>::SimplicialLDLT()

Default constructor

template<typename _MatrixType, int _UpLo, typename _Ordering>
Eigen::SimplicialLDLT<_MatrixType, _UpLo, _Ordering>::SimplicialLDLT(const MatrixType& matrix) explicit

Constructs and performs the LLT factorization of matrix

template<typename _MatrixType, int _UpLo, typename _Ordering>
void Eigen::SimplicialLDLT<_MatrixType, _UpLo, _Ordering>::analyzePattern(const MatrixType& a)

Performs a symbolic decomposition on the sparcity of matrix.

This function is particularly useful when solving for several problems having the same structure.

template<typename _MatrixType, int _UpLo, typename _Ordering>
SimplicialLDLT& Eigen::SimplicialLDLT<_MatrixType, _UpLo, _Ordering>::compute(const MatrixType& matrix)

Computes the sparse Cholesky decomposition of matrix

template<typename _MatrixType, int _UpLo, typename _Ordering>
Scalar Eigen::SimplicialLDLT<_MatrixType, _UpLo, _Ordering>::determinant() const

Returns the determinant of the underlying matrix from the current factorization

template<typename _MatrixType, int _UpLo, typename _Ordering>
void Eigen::SimplicialLDLT<_MatrixType, _UpLo, _Ordering>::factorize(const MatrixType& a)

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.

template<typename _MatrixType, int _UpLo, typename _Ordering>
const MatrixL Eigen::SimplicialLDLT<_MatrixType, _UpLo, _Ordering>::matrixL() const

Returns an expression of the factor L

template<typename _MatrixType, int _UpLo, typename _Ordering>
const MatrixU Eigen::SimplicialLDLT<_MatrixType, _UpLo, _Ordering>::matrixU() const

Returns an expression of the factor U (= L^*)

template<typename _MatrixType, int _UpLo, typename _Ordering>
const VectorType Eigen::SimplicialLDLT<_MatrixType, _UpLo, _Ordering>::vectorD() const

Returns a vector expression of the diagonal D