Eigen  3.3.3
Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering > Class Template Reference

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:
_MatrixTypethe type of the sparse matrix A, it must be a SparseMatrix<>
_UpLothe triangular part that will be used for the computations. It can be Lower or Upper. Default is Lower.
_OrderingThe 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
+ Inheritance diagram for Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >:

List of all members.

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)

Constructor & Destructor Documentation

template<typename _MatrixType , int _UpLo, typename _Ordering >
Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::SimplicialLLT ( ) [inline]

Default constructor

template<typename _MatrixType , int _UpLo, typename _Ordering >
Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::SimplicialLLT ( const MatrixType &  matrix) [inline, explicit]

Constructs and performs the LLT factorization of matrix


Member Function Documentation

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

Reimplemented from Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >.

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

Reimplemented from Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >.

template<typename _MatrixType , int _UpLo, typename _Ordering >
const MatrixL Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::matrixL ( ) const [inline]
Returns:
an expression of the factor L
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^*)

The documentation for this class was generated from the following file:
 All Classes Functions Variables Typedefs Enumerations Enumerator Friends