![]() |
|
Modified Incomplete Cholesky with dual threshold. More...
#include <IncompleteCholesky.h>
Public Types | |
enum | { UpLo = UpLo_ } |
enum | { ColsAtCompileTime = Dynamic , MaxColsAtCompileTime = Dynamic } |
typedef NumTraits< Scalar >::Real | RealScalar |
typedef OrderingType_ | OrderingType |
typedef OrderingType::PermutationType | PermutationType |
typedef PermutationType::StorageIndex | StorageIndex |
typedef SparseMatrix< Scalar, ColMajor, StorageIndex > | FactorType |
typedef Matrix< Scalar, Dynamic, 1 > | VectorSx |
typedef Matrix< RealScalar, Dynamic, 1 > | VectorRx |
typedef Matrix< StorageIndex, Dynamic, 1 > | VectorIx |
typedef std::vector< std::list< StorageIndex > > | VectorList |
Public Member Functions | |
IncompleteCholesky () | |
template<typename MatrixType > | |
IncompleteCholesky (const MatrixType &matrix) | |
EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
ComputationInfo | info () const |
Reports whether previous computation was successful. More... | |
void | setInitialShift (RealScalar shift) |
Set the initial shift parameter \( \sigma \). More... | |
template<typename MatrixType > | |
void | analyzePattern (const MatrixType &mat) |
Computes the fill reducing permutation vector using the sparsity pattern of mat. More... | |
template<typename MatrixType > | |
void | factorize (const MatrixType &mat) |
Performs the numerical factorization of the input matrix mat. More... | |
template<typename MatrixType > | |
void | compute (const MatrixType &mat) |
template<typename Rhs , typename Dest > | |
void | _solve_impl (const Rhs &b, Dest &x) const |
const FactorType & | matrixL () const |
const VectorRx & | scalingS () const |
const PermutationType & | permutationP () const |
RealScalar | shift () const |
template<typename MatrixType_ > | |
void | factorize (const MatrixType_ &mat) |
![]() | |
SparseSolverBase () | |
SparseSolverBase (SparseSolverBase &&other) | |
~SparseSolverBase () | |
IncompleteCholesky< Scalar, Lower, AMDOrdering< int > > & | derived () |
const IncompleteCholesky< Scalar, Lower, AMDOrdering< int > > & | derived () const |
const Solve< IncompleteCholesky< Scalar, Lower, AMDOrdering< int > >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
const Solve< IncompleteCholesky< Scalar, Lower, AMDOrdering< int > >, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
void | _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const |
Protected Types | |
typedef SparseSolverBase< IncompleteCholesky< Scalar, UpLo_, OrderingType_ > > | Base |
Protected Attributes | |
FactorType | m_L |
VectorRx | m_scale |
RealScalar | m_initialShift |
bool | m_analysisIsOk |
bool | m_factorizationIsOk |
ComputationInfo | m_info |
PermutationType | m_perm |
RealScalar | m_shift |
bool | m_isInitialized |
![]() | |
bool | m_isInitialized |
Private Member Functions | |
void | updateList (Ref< const VectorIx > colPtr, Ref< VectorIx > rowIdx, Ref< VectorSx > vals, const Index &col, const Index &jk, VectorIx &firstElt, VectorList &listCol) |
Modified Incomplete Cholesky with dual threshold.
References : C-J. Lin and J. J. Moré, Incomplete Cholesky Factorizations with Limited memory, SIAM J. Sci. Comput. 21(1), pp. 24-45, 1999
Scalar | the scalar type of the input matrices |
UpLo_ | The triangular part that will be used for the computations. It can be Lower or Upper. Default is Lower. |
OrderingType_ | The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<int>. |
\implsparsesolverconcept
It performs the following incomplete factorization: \( S P A P' S + \sigma I \approx L L' \) where L is a lower triangular factor, S is a diagonal scaling matrix, P is a fill-in reducing permutation as computed by the ordering method, and \( \sigma \) is a shift for ensuring the decomposed matrix is positive definite.
Shifting strategy: Let \( B = S P A P' S \) be the scaled matrix on which the factorization is carried out, and \( \beta \) be the minimum value of the diagonal. If \( \beta > 0 \) then, the factorization is directly performed on the matrix B, and \sigma = 0. Otherwise, the factorization is performed on the shifted matrix \( B + \sigma I \) for a shifting factor \( \sigma \). We start with \( \sigma = \sigma_0 - \beta \), where \( \sigma_0 \) is the initial shift value as returned and set by setInitialShift() method. The default value is \( \sigma_0 = 10^{-3} \). If the factorization fails, then the shift in doubled until it succeed or a maximum of ten attempts. If it still fails, as returned by the info() method, then you can either increase the initial shift, or better use another preconditioning technique.
|
protected |
typedef SparseMatrix<Scalar, ColMajor, StorageIndex> Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::FactorType |
typedef OrderingType_ Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::OrderingType |
typedef OrderingType::PermutationType Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::PermutationType |
typedef NumTraits<Scalar>::Real Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::RealScalar |
typedef PermutationType::StorageIndex Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::StorageIndex |
typedef Matrix<StorageIndex, Dynamic, 1> Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::VectorIx |
typedef std::vector<std::list<StorageIndex> > Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::VectorList |
typedef Matrix<RealScalar, Dynamic, 1> Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::VectorRx |
typedef Matrix<Scalar, Dynamic, 1> Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::VectorSx |
anonymous enum |
anonymous enum |
|
inline |
Default constructor leaving the object in a partly non-initialized stage.
You must call compute() or the pair analyzePattern()/factorize() to make it valid.
|
inline |
Constructor computing the incomplete factorization for the given matrix matrix.
References Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::compute(), and matrix().
|
inline |
References Eigen::SparseMatrixBase< Derived >::adjoint(), b, eigen_assert, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_factorizationIsOk, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_L, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_perm, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_scale, and plotDoE::x.
|
inline |
Computes the fill reducing permutation vector using the sparsity pattern of mat.
References Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_analysisIsOk, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_info, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_isInitialized, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_L, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_perm, Global_Physical_Variables::pinv, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::resize(), Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::rows(), and Eigen::Success.
Referenced by Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::compute().
|
inline |
References Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_L.
Referenced by gdb.printers._MatrixEntryIterator::__next__(), gdb.printers.EigenMatrixPrinter::children(), gdb.printers.EigenSparseMatrixPrinter::children(), gdb.printers.EigenMatrixPrinter::to_string(), and gdb.printers.EigenSparseMatrixPrinter::to_string().
|
inline |
Computes or re-computes the incomplete Cholesky factorization of the input matrix mat
It is a shortcut for a sequential call to the analyzePattern() and factorize() methods.
References Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::analyzePattern(), and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::factorize().
Referenced by Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::IncompleteCholesky().
void Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::factorize | ( | const MatrixType & | mat | ) |
Performs the numerical factorization of the input matrix mat.
The method analyzePattern() or compute() must have been called beforehand with a matrix having the same pattern.
Referenced by Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::compute().
void Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::factorize | ( | const MatrixType_ & | mat | ) |
References Eigen::numext::abs2(), Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), conj(), diag, eigen_assert, eigen_internal_assert, i, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::innerIndexPtr(), j, k, Eigen::numext::maxi(), min, Eigen::numext::mini(), n, Eigen::NumericalIssue, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::outerIndexPtr(), p, Eigen::internal::QuickSplit(), Eigen::PlainObjectBase< Derived >::resize(), Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::rows(), sqrt(), Eigen::Success, tmp, Eigen::SparseMatrixBase< Derived >::twistedBy(), and Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::valuePtr().
|
inline |
Reports whether previous computation was successful.
It triggers an assertion if *this
has not been initialized through the respective constructor, or a call to compute() or analyzePattern().
Success
if computation was successful, NumericalIssue
if the matrix appears to be negative. References eigen_assert, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_info, and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_isInitialized.
|
inline |
References eigen_assert, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_factorizationIsOk, and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_L.
|
inline |
References eigen_assert, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_analysisIsOk, and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_perm.
|
inline |
References Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_L, and Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::rows().
Referenced by gdb.printers._MatrixEntryIterator::__next__(), gdb.printers.EigenMatrixPrinter::children(), gdb.printers.EigenSparseMatrixPrinter::children(), gdb.printers.EigenMatrixPrinter::to_string(), and gdb.printers.EigenSparseMatrixPrinter::to_string().
|
inline |
References eigen_assert, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_factorizationIsOk, and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_scale.
|
inline |
Set the initial shift parameter \( \sigma \).
References Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_initialShift, and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::shift().
|
inline |
References Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_shift.
Referenced by Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::setInitialShift().
|
inlineprivate |
|
protected |
|
protected |
|
protected |
|
protected |
|
mutableprotected |
|
protected |
Referenced by Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::_solve_impl(), Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::analyzePattern(), Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::cols(), Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::matrixL(), and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::rows().
|
protected |
|
protected |
|
protected |