compbio
Public Types | Public Member Functions | List of all members
Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering > Class Template Reference

A direct sparse LLT Cholesky factorizations. More...

#include <SimplicialCholesky.h>

Inheritance diagram for Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >:
Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > > Eigen::SparseSolverBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > > Eigen::internal::noncopyable

Public Types

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

Public Member Functions

 SimplicialLLT ()
 Default constructor.
 
 SimplicialLLT (const MatrixType &matrix)
 Constructs and performs the LLT factorization of matrix.
 
const MatrixL matrixL () const
 
const MatrixU matrixU () const
 
SimplicialLLTcompute (const MatrixType &matrix)
 Computes the sparse Cholesky decomposition of matrix.
 
void analyzePattern (const MatrixType &a)
 Performs a symbolic decomposition on the sparcity of matrix. More...
 
void factorize (const MatrixType &a)
 Performs a numeric decomposition of matrix. More...
 
Scalar determinant () const
 
- Public Member Functions inherited from Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >
 SimplicialCholeskyBase ()
 Default constructor.
 
 SimplicialCholeskyBase (const MatrixType &matrix)
 
SimplicialLLT< _MatrixType, _UpLo, _Ordering > & derived ()
 
const SimplicialLLT< _MatrixType, _UpLo, _Ordering > & derived () const
 
Index cols () const
 
Index rows () 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
 
SimplicialLLT< _MatrixType, _UpLo, _Ordering > & setShift (const RealScalar &offset, const RealScalar &scale=1)
 Sets the shift parameters that will be used to adjust the diagonal coefficients during the numerical factorization. More...
 
void dumpMemory (Stream &s)
 
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
 
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
- Public Member Functions inherited from Eigen::SparseSolverBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >
 SparseSolverBase ()
 Default constructor.
 
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
 
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Additional Inherited Members

- Protected Member Functions inherited from Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >
void compute (const MatrixType &matrix)
 Computes the sparse Cholesky decomposition of matrix.
 
void factorize (const MatrixType &a)
 
void factorize_preordered (const CholMatrixType &a)
 
void analyzePattern (const MatrixType &a, bool doLDLT)
 
void analyzePattern_preordered (const CholMatrixType &a, bool doLDLT)
 
void ordering (const MatrixType &a, ConstCholMatrixPtr &pmat, CholMatrixType &ap)
 
- Protected Attributes inherited from Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >
ComputationInfo m_info
 
bool m_factorizationIsOk
 
bool m_analysisIsOk
 
CholMatrixType m_matrix
 
VectorType m_diag
 
VectorI m_parent
 
VectorI m_nonZerosPerCol
 
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_P
 
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_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
_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<>
See also
class SimplicialLDLT, class AMDOrdering, class NaturalOrdering

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

§ 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

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

§ 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

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

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