$treeview $search $mathjax
Eigen
3.2.5
$projectbrief
|
$projectbrief
|
$searchbox |
A direct sparse LLT Cholesky factorizations. More...
Public Member Functions | |
void | analyzePattern (const MatrixType &a) |
SimplicialLLT & | compute (const MatrixType &matrix) |
Scalar | determinant () const |
void | factorize (const MatrixType &a) |
ComputationInfo | info () const |
Reports whether previous computation was successful. | |
const MatrixL | matrixL () const |
const MatrixU | matrixU () const |
const PermutationMatrix < Dynamic, Dynamic, Index > & | permutationP () const |
const PermutationMatrix < Dynamic, Dynamic, Index > & | permutationPinv () const |
SimplicialLLT< _MatrixType, _UpLo, _Ordering > & | setShift (const RealScalar &offset, const RealScalar &scale=1) |
SimplicialLLT (const MatrixType &matrix) | |
SimplicialLLT () | |
const internal::sparse_solve_retval < SimplicialCholeskyBase, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
const internal::solve_retval < SimplicialCholeskyBase, Rhs > | solve (const MatrixBase< Rhs > &b) const |
Protected Member Functions | |
void | compute (const MatrixType &matrix) |
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.
_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<> |
SimplicialLLT | ( | ) | [inline] |
Default constructor
SimplicialLLT | ( | const MatrixType & | matrix | ) | [inline] |
Constructs and performs the LLT factorization of matrix
void 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.
void compute | ( | const MatrixType & | matrix | ) | [inline, protected, inherited] |
Computes the sparse Cholesky decomposition of matrix
SimplicialLLT& compute | ( | const MatrixType & | matrix | ) | [inline] |
Computes the sparse Cholesky decomposition of matrix
Scalar determinant | ( | ) | const [inline] |
void 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.
ComputationInfo info | ( | ) | const [inline, inherited] |
Reports whether previous computation was successful.
Success
if computation was succesful, NumericalIssue
if the matrix.appears to be negative. const MatrixL matrixL | ( | ) | const [inline] |
const MatrixU matrixU | ( | ) | const [inline] |
const PermutationMatrix<Dynamic,Dynamic,Index>& permutationP | ( | ) | const [inline, inherited] |
const PermutationMatrix<Dynamic,Dynamic,Index>& permutationPinv | ( | ) | const [inline, inherited] |
SimplicialLLT< _MatrixType, _UpLo, _Ordering > & setShift | ( | const RealScalar & | offset, | |
const RealScalar & | scale = 1 | |||
) | [inline, inherited] |
Sets the shift parameters that will be used to adjust the diagonal coefficients during the numerical factorization.
During the numerical factorization, the diagonal coefficients are transformed by the following linear model:
d_ii
= offset + scale * d_ii
The default is the identity transformation with offset=0, and scale=1.
*this
. const internal::sparse_solve_retval<SimplicialCholeskyBase, Rhs> solve | ( | const SparseMatrixBase< Rhs > & | b | ) | const [inline, inherited] |
const internal::solve_retval<SimplicialCholeskyBase, Rhs> solve | ( | const MatrixBase< Rhs > & | b | ) | const [inline, inherited] |