11#ifndef EIGEN_INCOMPLETE_CHOlESKY_H
12#define EIGEN_INCOMPLETE_CHOlESKY_H
44template <
typename Scalar,
int _UpLo = Lower,
typename _OrderingType = AMDOrdering<
int> >
59 typedef std::vector<std::list<StorageIndex> >
VectorList;
77 template<
typename MatrixType>
110 template<
typename MatrixType>
115 ord(mat.template selfadjointView<UpLo>(), pinv);
116 if(pinv.size()>0)
m_perm = pinv.inverse();
131 template<
typename MatrixType>
140 template<
typename MatrixType>
148 template<
typename Rhs,
typename Dest>
155 x =
m_L.template triangularView<Lower>().solve(x);
156 x =
m_L.
adjoint().template triangularView<Upper>().solve(x);
188template<
typename Scalar,
int _UpLo,
typename OrderingType>
189template<
typename _MatrixType>
193 eigen_assert(m_analysisIsOk &&
"analyzePattern() should be called first");
198 if (m_perm.rows() == mat.rows() )
202 tmp = mat.template selfadjointView<_UpLo>().
twistedBy(m_perm);
203 m_L.template selfadjointView<Lower>() = tmp.template selfadjointView<Lower>();
207 m_L.template selfadjointView<Lower>() = mat.template selfadjointView<_UpLo>();
210 Index n = m_L.cols();
211 Index nnz = m_L.nonZeros();
220 col_pattern.fill(-1);
227 for (
Index j = 0; j < n; j++)
228 for (
Index k = colPtr[j]; k < colPtr[j+1]; k++)
235 m_scale = m_scale.cwiseSqrt().cwiseSqrt();
237 for (
Index j = 0; j < n; ++j)
247 for (
Index j = 0; j < n; j++)
249 for (
Index k = colPtr[j]; k < colPtr[j+1]; k++)
250 vals[k] *= (m_scale(j)*m_scale(rowIdx[k]));
251 eigen_internal_assert(rowIdx[colPtr[j]]==j &&
"IncompleteCholesky: only the lower triangular part must be stored");
259 shift = m_initialShift - mindiag;
268 for (
Index j = 0; j < n; j++)
269 vals[colPtr[j]] += shift;
277 Scalar diag = vals[colPtr[j]];
279 for (
Index i = colPtr[j] + 1; i < colPtr[j+1]; i++)
282 col_vals(col_nnz) = vals[i];
283 col_irow(col_nnz) = l;
284 col_pattern(l) = col_nnz;
288 typename std::list<StorageIndex>::iterator k;
290 for(k = listCol[j].
begin(); k != listCol[j].end(); k++)
292 Index jk = firstElt(*k);
294 Scalar v_j_jk = numext::conj(vals[jk]);
297 for (
Index i = jk; i < colPtr[*k+1]; i++)
302 col_vals(col_nnz) = vals[i] * v_j_jk;
303 col_irow[col_nnz] = l;
304 col_pattern(l) = col_nnz;
308 col_vals(col_pattern[l]) -= vals[i] * v_j_jk;
310 updateList(colPtr,rowIdx,vals, *k, jk, firstElt, listCol);
326 col_pattern.fill(-1);
327 for(
Index i=0; i<n; ++i)
334 vals[colPtr[j]] = rdiag;
335 for (
Index k = 0; k<col_nnz; ++k)
337 Index i = col_irow[k];
339 col_vals(k) /= rdiag;
345 Index p = colPtr[j+1] - colPtr[j] - 1 ;
351 for (
Index i = colPtr[j]+1; i < colPtr[j+1]; i++)
353 vals[i] = col_vals(cpt);
354 rowIdx[i] = col_irow(cpt);
356 col_pattern(col_irow(cpt)) = -1;
360 Index jk = colPtr(j)+1;
361 updateList(colPtr,rowIdx,vals,j,jk,firstElt,listCol);
366 m_factorizationIsOk =
true;
372template<
typename Scalar,
int _UpLo,
typename OrderingType>
375 if (jk < colPtr(
col+1) )
379 rowIdx.segment(jk,p).minCoeff(&minpos);
381 if (rowIdx(minpos) != rowIdx(jk))
387 firstElt(
col) = internal::convert_index<StorageIndex,Index>(jk);
388 listCol[rowIdx(jk)].push_back(internal::convert_index<StorageIndex,Index>(
col));
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ColXpr col(Index i)
This is the const version of col().
Definition: BlockMethods.h:1097
EIGEN_DEVICE_FUNC RealReturnType real() const
Definition: CommonCwiseUnaryOps.h:100
#define eigen_internal_assert(x)
Definition: Macros.h:1053
#define EIGEN_NOEXCEPT
Definition: Macros.h:1428
#define EIGEN_CONSTEXPR
Definition: Macros.h:797
#define eigen_assert(x)
Definition: Macros.h:1047
Modified Incomplete Cholesky with dual threshold.
Definition: IncompleteCholesky.h:46
@ MaxColsAtCompileTime
Definition: IncompleteCholesky.h:63
@ ColsAtCompileTime
Definition: IncompleteCholesky.h:62
PermutationType::StorageIndex StorageIndex
Definition: IncompleteCholesky.h:54
Matrix< StorageIndex, Dynamic, 1 > VectorIx
Definition: IncompleteCholesky.h:58
SparseMatrix< Scalar, ColMajor, StorageIndex > FactorType
Definition: IncompleteCholesky.h:55
RealScalar m_initialShift
Definition: IncompleteCholesky.h:174
void setInitialShift(RealScalar shift)
Set the initial shift parameter .
Definition: IncompleteCholesky.h:106
Matrix< Scalar, Dynamic, 1 > VectorSx
Definition: IncompleteCholesky.h:56
bool m_analysisIsOk
Definition: IncompleteCholesky.h:175
void analyzePattern(const MatrixType &mat)
Computes the fill reducing permutation vector using the sparsity pattern of mat.
Definition: IncompleteCholesky.h:111
IncompleteCholesky(const MatrixType &matrix)
Constructor computing the incomplete factorization for the given matrix matrix.
Definition: IncompleteCholesky.h:78
@ UpLo
Definition: IncompleteCholesky.h:60
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: IncompleteCholesky.h:84
const VectorRx & scalingS() const
Definition: IncompleteCholesky.h:166
void compute(const MatrixType &mat)
Computes or re-computes the incomplete Cholesky factorization of the input matrix mat.
Definition: IncompleteCholesky.h:141
void _solve_impl(const Rhs &b, Dest &x) const
Definition: IncompleteCholesky.h:149
SparseSolverBase< IncompleteCholesky< Scalar, _UpLo, _OrderingType > > Base
Definition: IncompleteCholesky.h:48
VectorRx m_scale
Definition: IncompleteCholesky.h:173
std::vector< std::list< StorageIndex > > VectorList
Definition: IncompleteCholesky.h:59
OrderingType::PermutationType PermutationType
Definition: IncompleteCholesky.h:53
bool m_factorizationIsOk
Definition: IncompleteCholesky.h:176
const PermutationType & permutationP() const
Definition: IncompleteCholesky.h:169
ComputationInfo m_info
Definition: IncompleteCholesky.h:177
void factorize(const MatrixType &mat)
Performs the numerical factorization of the input matrix mat.
bool m_isInitialized
Definition: SparseSolverBase.h:119
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: IncompleteCholesky.h:87
_OrderingType OrderingType
Definition: IncompleteCholesky.h:52
FactorType m_L
Definition: IncompleteCholesky.h:172
PermutationType m_perm
Definition: IncompleteCholesky.h:178
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: IncompleteCholesky.h:98
IncompleteCholesky()
Default constructor leaving the object in a partly non-initialized stage.
Definition: IncompleteCholesky.h:73
const FactorType & matrixL() const
Definition: IncompleteCholesky.h:163
NumTraits< Scalar >::Real RealScalar
Definition: IncompleteCholesky.h:51
Matrix< RealScalar, Dynamic, 1 > VectorRx
Definition: IncompleteCholesky.h:57
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Resizes *this to a rows x cols matrix.
Definition: PlainObjectBase.h:271
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:283
SparseSymmetricPermutationProduct< Derived, Upper|Lower > twistedBy(const PermutationMatrix< Dynamic, Dynamic, StorageIndex > &perm) const
Definition: SparseMatrixBase.h:329
const AdjointReturnType adjoint() const
Definition: SparseMatrixBase.h:356
Index rows() const
Definition: SparseMatrix.h:138
const StorageIndex * innerIndexPtr() const
Definition: SparseMatrix.h:159
const StorageIndex * outerIndexPtr() const
Definition: SparseMatrix.h:168
Index cols() const
Definition: SparseMatrix.h:140
const Scalar * valuePtr() const
Definition: SparseMatrix.h:150
void resize(Index rows, Index cols)
Resizes the matrix to a rows x cols matrix and initializes it to zero.
Definition: SparseMatrix.h:626
A base class for sparse solvers.
Definition: SparseSolverBase.h:68
bool m_isInitialized
Definition: SparseSolverBase.h:119
auto sqrt(const UnitType &value) noexcept -> unit_t< square_root< typename units::traits::unit_t_traits< UnitType >::unit_type >, typename units::traits::unit_t_traits< UnitType >::underlying_type, linear_scale >
computes the square root of value
Definition: math.h:483
ComputationInfo
Enum for reporting the status of a computation.
Definition: Constants.h:440
@ NumericalIssue
The provided data did not satisfy the prerequisites.
Definition: Constants.h:444
@ Success
Computation was successful.
Definition: Constants.h:442
constexpr common_t< T1, T2 > min(const T1 x, const T2 y) noexcept
Compile-time pairwise minimum function.
Definition: min.hpp:35
Index QuickSplit(VectorV &row, VectorI &ind, Index ncut)
Definition: IncompleteLUT.h:29
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Definition: MathFunctions.h:1091
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
Definition: MathFunctions.h:1083
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition: MathFunctions.h:1292
Namespace containing all symbols from the Eigen library.
Definition: MatrixExponential.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
const int Dynamic
This value means that a positive quantity (e.g., a size) is not known at compile-time,...
Definition: Constants.h:22
GHC_FS_API directory_iterator begin(directory_iterator iter) noexcept
Definition: filesystem.hpp:5746
void swap(wpi::SmallVectorImpl< T > &LHS, wpi::SmallVectorImpl< T > &RHS)
Implement std::swap in terms of SmallVector swap.
Definition: SmallVector.h:1299
static constexpr const charge::coulomb_t e(1.6021766208e-19)
elementary charge.
Holds information about the various numeric (i.e.
Definition: NumTraits.h:233