A direct sparse LDLT Cholesky factorizations without square root. More...
#include <SimplicialCholesky.h>
Public Types | |
enum | { UpLo = _UpLo } |
typedef SimplicialCholeskyBase< SimplicialLDLT > | Base |
typedef SparseMatrix< Scalar, ColMajor, Index > | CholMatrixType |
typedef MatrixType::Index | Index |
typedef Traits::MatrixL | MatrixL |
typedef _MatrixType | MatrixType |
typedef Traits::MatrixU | MatrixU |
typedef MatrixType::RealScalar | RealScalar |
typedef MatrixType::Scalar | Scalar |
typedef internal::traits< SimplicialLDLT > | Traits |
typedef Matrix< Scalar, Dynamic, 1 > | VectorType |
Public Types inherited from Eigen::SimplicialCholeskyBase< SimplicialLDLT< _MatrixType, _UpLo > > | |
enum | |
typedef SparseMatrix< Scalar, ColMajor, Index > | CholMatrixType |
typedef MatrixType::Index | Index |
typedef internal::traits< SimplicialLDLT< _MatrixType, _UpLo > >::MatrixType | MatrixType |
typedef MatrixType::RealScalar | RealScalar |
typedef MatrixType::Scalar | Scalar |
typedef Matrix< Scalar, Dynamic, 1 > | VectorType |
Additional Inherited Members | |
Protected Member Functions inherited from Eigen::SimplicialCholeskyBase< SimplicialLDLT< _MatrixType, _UpLo > > | |
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, CholMatrixType &ap) |
Protected Attributes inherited from Eigen::SimplicialCholeskyBase< SimplicialLDLT< _MatrixType, _UpLo > > | |
bool | m_analysisIsOk |
VectorType | m_diag |
bool | m_factorizationIsOk |
ComputationInfo | m_info |
bool | m_isInitialized |
CholMatrixType | m_matrix |
VectorXi | m_nonZerosPerCol |
PermutationMatrix< Dynamic, Dynamic, Index > | m_P |
VectorXi | m_parent |
PermutationMatrix< Dynamic, Dynamic, Index > | m_Pinv |
RealScalar | m_shiftOffset |
RealScalar | m_shiftScale |
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.
_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. |
Definition at line 244 of file SimplicialCholesky.h.
typedef SimplicialCholeskyBase<SimplicialLDLT> Eigen::SimplicialLDLT< _MatrixType, _UpLo >::Base |
Definition at line 394 of file SimplicialCholesky.h.
typedef SparseMatrix<Scalar,ColMajor,Index> Eigen::SimplicialLDLT< _MatrixType, _UpLo >::CholMatrixType |
Definition at line 398 of file SimplicialCholesky.h.
typedef MatrixType::Index Eigen::SimplicialLDLT< _MatrixType, _UpLo >::Index |
Definition at line 397 of file SimplicialCholesky.h.
typedef Traits::MatrixL Eigen::SimplicialLDLT< _MatrixType, _UpLo >::MatrixL |
Definition at line 401 of file SimplicialCholesky.h.
typedef _MatrixType Eigen::SimplicialLDLT< _MatrixType, _UpLo >::MatrixType |
Definition at line 392 of file SimplicialCholesky.h.
typedef Traits::MatrixU Eigen::SimplicialLDLT< _MatrixType, _UpLo >::MatrixU |
Definition at line 402 of file SimplicialCholesky.h.
typedef MatrixType::RealScalar Eigen::SimplicialLDLT< _MatrixType, _UpLo >::RealScalar |
Definition at line 396 of file SimplicialCholesky.h.
typedef MatrixType::Scalar Eigen::SimplicialLDLT< _MatrixType, _UpLo >::Scalar |
Definition at line 395 of file SimplicialCholesky.h.
typedef internal::traits<SimplicialLDLT> Eigen::SimplicialLDLT< _MatrixType, _UpLo >::Traits |
Definition at line 400 of file SimplicialCholesky.h.
typedef Matrix<Scalar,Dynamic,1> Eigen::SimplicialLDLT< _MatrixType, _UpLo >::VectorType |
Definition at line 399 of file SimplicialCholesky.h.
anonymous enum |
Enumerator | |
---|---|
UpLo |
Definition at line 393 of file SimplicialCholesky.h.
|
inline |
Default constructor
Definition at line 405 of file SimplicialCholesky.h.
|
inline |
Constructs and performs the LLT factorization of matrix
Definition at line 408 of file SimplicialCholesky.h.
|
inline |
Performs a symbolic decomposition on the sparcity of matrix.
This function is particularly useful when solving for several problems having the same structure.
Definition at line 441 of file SimplicialCholesky.h.
|
inline |
Computes the sparse Cholesky decomposition of matrix
Definition at line 429 of file SimplicialCholesky.h.
|
inline |
Definition at line 458 of file SimplicialCholesky.h.
|
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.
Definition at line 452 of file SimplicialCholesky.h.
|
inline |
Definition at line 417 of file SimplicialCholesky.h.
|
inline |
Definition at line 423 of file SimplicialCholesky.h.
|
inline |
Definition at line 412 of file SimplicialCholesky.h.