11#ifndef EIGEN_INCOMPLETE_CHOlESKY_H
12#define EIGEN_INCOMPLETE_CHOlESKY_H
44template <
typename Scalar,
int _UpLo = Lower,
typename _OrderingType = AMDOrdering<
int> >
49 using Base::m_isInitialized;
51 typedef typename NumTraits<Scalar>::Real RealScalar;
52 typedef _OrderingType OrderingType;
53 typedef typename OrderingType::PermutationType PermutationType;
54 typedef typename PermutationType::StorageIndex StorageIndex;
59 typedef std::vector<std::list<StorageIndex> > VectorList;
60 enum { UpLo = _UpLo };
77 template<
typename MatrixType>
78 IncompleteCholesky(
const MatrixType& matrix) : m_initialShift(1e-3),m_analysisIsOk(false),m_factorizationIsOk(false)
84 EIGEN_CONSTEXPR
Index rows() const EIGEN_NOEXCEPT {
return m_L.
rows(); }
87 EIGEN_CONSTEXPR
Index cols() const EIGEN_NOEXCEPT {
return m_L.
cols(); }
100 eigen_assert(m_isInitialized &&
"IncompleteCholesky is not initialized.");
110 template<
typename MatrixType>
114 PermutationType pinv;
115 ord(mat.template selfadjointView<UpLo>(), pinv);
116 if(pinv.size()>0) m_perm = pinv.inverse();
117 else m_perm.resize(0);
118 m_L.
resize(mat.rows(), mat.cols());
119 m_analysisIsOk =
true;
120 m_isInitialized =
true;
131 template<
typename MatrixType>
140 template<
typename MatrixType>
148 template<
typename Rhs,
typename Dest>
149 void _solve_impl(
const Rhs& b, Dest& x)
const
151 eigen_assert(m_factorizationIsOk &&
"factorize() should be called first");
152 if (m_perm.rows() == b.rows()) x = m_perm * b;
154 x = m_scale.asDiagonal() * x;
155 x = m_L.template triangularView<Lower>().solve(x);
156 x = m_L.adjoint().template triangularView<Upper>().solve(x);
157 x = m_scale.asDiagonal() * x;
158 if (m_perm.rows() == b.rows())
159 x = m_perm.inverse() * x;
169 const PermutationType&
permutationP()
const { eigen_assert(
"m_analysisIsOk");
return m_perm; }
174 RealScalar m_initialShift;
176 bool m_factorizationIsOk;
178 PermutationType m_perm;
188template<
typename Scalar,
int _UpLo,
typename OrderingType>
189template<
typename _MatrixType>
190void IncompleteCholesky<Scalar,_UpLo, OrderingType>::factorize(
const _MatrixType& mat)
193 eigen_assert(m_analysisIsOk &&
"analyzePattern() should be called first");
198 if (m_perm.rows() == mat.rows() )
201 FactorType tmp(mat.rows(), mat.cols());
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();
212 Map<VectorSx> vals(m_L.valuePtr(), nnz);
213 Map<VectorIx> rowIdx(m_L.innerIndexPtr(), nnz);
214 Map<VectorIx> colPtr( m_L.outerIndexPtr(), n+1);
215 VectorIx firstElt(n-1);
216 VectorList listCol(n);
217 VectorSx col_vals(n);
218 VectorIx col_irow(n);
219 VectorIx col_pattern(n);
220 col_pattern.fill(-1);
221 StorageIndex col_nnz;
227 for (Index j = 0; j < n; j++)
228 for (Index k = colPtr[j]; k < colPtr[j+1]; k++)
230 m_scale(j) += numext::abs2(vals(k));
232 m_scale(rowIdx[k]) += numext::abs2(vals(k));
235 m_scale = m_scale.cwiseSqrt().cwiseSqrt();
237 for (Index j = 0; j < n; ++j)
238 if(m_scale(j)>(std::numeric_limits<RealScalar>::min)())
239 m_scale(j) = RealScalar(1)/m_scale(j);
246 RealScalar mindiag = NumTraits<RealScalar>::highest();
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");
252 mindiag = numext::mini(numext::real(vals[colPtr[j]]), mindiag);
255 FactorType L_save = m_L;
257 RealScalar shift = 0;
258 if(mindiag <= RealScalar(0.))
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++)
281 StorageIndex l = rowIdx[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);
293 eigen_internal_assert(rowIdx[jk]==j);
294 Scalar v_j_jk = numext::conj(vals[jk]);
297 for (Index i = jk; i < colPtr[*k+1]; i++)
299 StorageIndex l = rowIdx[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);
315 if(numext::real(diag) <= 0)
321 shift = numext::maxi(m_initialShift,RealScalar(2)*shift);
323 vals = Map<const VectorSx>(L_save.valuePtr(), nnz);
324 rowIdx = Map<const VectorIx>(L_save.innerIndexPtr(), nnz);
325 colPtr = Map<const VectorIx>(L_save.outerIndexPtr(), n+1);
326 col_pattern.fill(-1);
327 for(Index i=0; i<n; ++i)
333 RealScalar rdiag = sqrt(numext::real(diag));
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;
341 vals[colPtr[i]] -= numext::abs2(col_vals(k));
345 Index p = colPtr[j+1] - colPtr[j] - 1 ;
346 Ref<VectorSx> cvals = col_vals.head(col_nnz);
347 Ref<VectorIx> cirow = col_irow.head(col_nnz);
348 internal::QuickSplit(cvals,cirow, p);
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;
369 }
while(m_info!=Success);
372template<
typename Scalar,
int _UpLo,
typename OrderingType>
373inline void IncompleteCholesky<Scalar,_UpLo, OrderingType>::updateList(Ref<const VectorIx> colPtr, Ref<VectorIx> rowIdx, Ref<VectorSx> vals,
const Index& col,
const Index& jk, VectorIx& firstElt, VectorList& listCol)
375 if (jk < colPtr(col+1) )
377 Index p = colPtr(col+1) - jk;
379 rowIdx.segment(jk,p).minCoeff(&minpos);
381 if (rowIdx(minpos) != rowIdx(jk))
384 std::swap(rowIdx(jk),rowIdx(minpos));
385 std::swap(vals(jk),vals(minpos));
387 firstElt(col) = internal::convert_index<StorageIndex,Index>(jk);
388 listCol[rowIdx(jk)].push_back(internal::convert_index<StorageIndex,Index>(col));
Modified Incomplete Cholesky with dual threshold.
Definition: IncompleteCholesky.h:46
void setInitialShift(RealScalar shift)
Set the initial shift parameter .
Definition: IncompleteCholesky.h:106
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)
Definition: IncompleteCholesky.h:78
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: IncompleteCholesky.h:84
const VectorRx & scalingS() const
Definition: IncompleteCholesky.h:166
void compute(const MatrixType &mat)
Definition: IncompleteCholesky.h:141
const PermutationType & permutationP() const
Definition: IncompleteCholesky.h:169
void factorize(const MatrixType &mat)
Performs the numerical factorization of the input matrix mat.
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: IncompleteCholesky.h:87
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: IncompleteCholesky.h:98
IncompleteCholesky()
Definition: IncompleteCholesky.h:73
const FactorType & matrixL() const
Definition: IncompleteCholesky.h:163
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:180
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:283
A versatible sparse matrix representation.
Definition: SparseMatrix.h:98
Index rows() const
Definition: SparseMatrix.h:138
Index cols() const
Definition: SparseMatrix.h:140
void resize(Index rows, Index cols)
Definition: SparseMatrix.h:626
A base class for sparse solvers.
Definition: SparseSolverBase.h:68
ComputationInfo
Definition: Constants.h:440
@ NumericalIssue
Definition: Constants.h:444
@ Success
Definition: Constants.h:442
Namespace containing all symbols from the Eigen library.
Definition: Core:141
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
const int Dynamic
Definition: Constants.h:22