IncompleteCholesky.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
5 // Copyright (C) 2015 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_INCOMPLETE_CHOlESKY_H
12 #define EIGEN_INCOMPLETE_CHOlESKY_H
13 
14 #include <vector>
15 #include <list>
16 
17 // IWYU pragma: private
18 #include "./InternalHeaderCheck.h"
19 
20 namespace Eigen {
50 template <typename Scalar, int UpLo_ = Lower, typename OrderingType_ = AMDOrdering<int> >
51 class IncompleteCholesky : public SparseSolverBase<IncompleteCholesky<Scalar, UpLo_, OrderingType_> > {
52  protected:
55 
56  public:
58  typedef OrderingType_ OrderingType;
59  typedef typename OrderingType::PermutationType PermutationType;
60  typedef typename PermutationType::StorageIndex StorageIndex;
65  typedef std::vector<std::list<StorageIndex> > VectorList;
66  enum { UpLo = UpLo_ };
68 
69  public:
77 
80  template <typename MatrixType>
82  : m_initialShift(1e-3), m_analysisIsOk(false), m_factorizationIsOk(false) {
83  compute(matrix);
84  }
85 
88 
91 
101  eigen_assert(m_isInitialized && "IncompleteCholesky is not initialized.");
102  return m_info;
103  }
104 
108 
111  template <typename MatrixType>
113  OrderingType ord;
115  ord(mat.template selfadjointView<UpLo>(), pinv);
116  if (pinv.size() > 0)
117  m_perm = pinv.inverse();
118  else
119  m_perm.resize(0);
120  m_L.resize(mat.rows(), mat.cols());
121  m_analysisIsOk = true;
122  m_isInitialized = true;
123  m_info = Success;
124  }
125 
133  template <typename MatrixType>
134  void factorize(const MatrixType& mat);
135 
142  template <typename MatrixType>
143  void compute(const MatrixType& mat) {
145  factorize(mat);
146  }
147 
148  // internal
149  template <typename Rhs, typename Dest>
150  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())
153  x = m_perm * b;
154  else
155  x = b;
156  x = m_scale.asDiagonal() * x;
157  x = m_L.template triangularView<Lower>().solve(x);
158  x = m_L.adjoint().template triangularView<Upper>().solve(x);
159  x = m_scale.asDiagonal() * x;
160  if (m_perm.rows() == b.rows()) x = m_perm.inverse() * x;
161  }
162 
164  const FactorType& matrixL() const {
165  eigen_assert(m_factorizationIsOk && "factorize() should be called first");
166  return m_L;
167  }
168 
170  const VectorRx& scalingS() const {
171  eigen_assert(m_factorizationIsOk && "factorize() should be called first");
172  return m_scale;
173  }
174 
176  const PermutationType& permutationP() const {
177  eigen_assert(m_analysisIsOk && "analyzePattern() should be called first");
178  return m_perm;
179  }
180 
182  RealScalar shift() const { return m_shift; }
183 
184  protected:
185  FactorType m_L; // The lower part stored in CSC
186  VectorRx m_scale; // The vector for scaling the matrix
187  RealScalar m_initialShift; // The initial shift parameter
192  RealScalar m_shift; // The final shift parameter.
193 
194  private:
195  inline void updateList(Ref<const VectorIx> colPtr, Ref<VectorIx> rowIdx, Ref<VectorSx> vals, const Index& col,
196  const Index& jk, VectorIx& firstElt, VectorList& listCol);
197 };
198 
199 // Based on the following paper:
200 // C-J. Lin and J. J. Moré, Incomplete Cholesky Factorizations with
201 // Limited memory, SIAM J. Sci. Comput. 21(1), pp. 24-45, 1999
202 // http://ftp.mcs.anl.gov/pub/tech_reports/reports/P682.pdf
203 template <typename Scalar, int UpLo_, typename OrderingType>
204 template <typename MatrixType_>
206  using std::sqrt;
207  eigen_assert(m_analysisIsOk && "analyzePattern() should be called first");
208 
209  // Dropping strategy : Keep only the p largest elements per column, where p is the number of elements in the column of
210  // the original matrix. Other strategies will be added
211 
212  // Apply the fill-reducing permutation computed in analyzePattern()
213  if (m_perm.rows() == mat.rows()) // To detect the null permutation
214  {
215  // The temporary is needed to make sure that the diagonal entry is properly sorted
216  FactorType tmp(mat.rows(), mat.cols());
217  tmp = mat.template selfadjointView<UpLo_>().twistedBy(m_perm);
218  m_L.template selfadjointView<Lower>() = tmp.template selfadjointView<Lower>();
219  } else {
220  m_L.template selfadjointView<Lower>() = mat.template selfadjointView<UpLo_>();
221  }
222 
223  // The algorithm will insert increasingly large shifts on the diagonal until
224  // factorization succeeds. Therefore we have to make sure that there is a
225  // space in the datastructure to store such values, even if the original
226  // matrix has a zero on the diagonal.
227  bool modified = false;
228  for (Index i = 0; i < mat.cols(); ++i) {
229  bool inserted = false;
230  m_L.findOrInsertCoeff(i, i, &inserted);
231  if (inserted) {
232  modified = true;
233  }
234  }
235  if (modified) m_L.makeCompressed();
236 
237  Index n = m_L.cols();
238  Index nnz = m_L.nonZeros();
239  Map<VectorSx> vals(m_L.valuePtr(), nnz); // values
240  Map<VectorIx> rowIdx(m_L.innerIndexPtr(), nnz); // Row indices
241  Map<VectorIx> colPtr(m_L.outerIndexPtr(), n + 1); // Pointer to the beginning of each row
242  VectorIx firstElt(n - 1); // for each j, points to the next entry in vals that will be used in the factorization
243  VectorList listCol(n); // listCol(j) is a linked list of columns to update column j
244  VectorSx col_vals(n); // Store a nonzero values in each column
245  VectorIx col_irow(n); // Row indices of nonzero elements in each column
246  VectorIx col_pattern(n);
247  col_pattern.fill(-1);
248  StorageIndex col_nnz;
249 
250  // Computes the scaling factors
251  m_scale.resize(n);
252  m_scale.setZero();
253  for (Index j = 0; j < n; j++)
254  for (Index k = colPtr[j]; k < colPtr[j + 1]; k++) {
255  m_scale(j) += numext::abs2(vals(k));
256  if (rowIdx[k] != j) m_scale(rowIdx[k]) += numext::abs2(vals(k));
257  }
258 
259  m_scale = m_scale.cwiseSqrt().cwiseSqrt();
260 
261  for (Index j = 0; j < n; ++j)
262  if (m_scale(j) > (std::numeric_limits<RealScalar>::min)())
263  m_scale(j) = RealScalar(1) / m_scale(j);
264  else
265  m_scale(j) = 1;
266 
267  // TODO disable scaling if not needed, i.e., if it is roughly uniform? (this will make solve() faster)
268 
269  // Scale and compute the shift for the matrix
271  for (Index j = 0; j < n; j++) {
272  for (Index k = colPtr[j]; k < colPtr[j + 1]; k++) vals[k] *= (m_scale(j) * m_scale(rowIdx[k]));
273  eigen_internal_assert(rowIdx[colPtr[j]] == j &&
274  "IncompleteCholesky: only the lower triangular part must be stored");
275  mindiag = numext::mini(numext::real(vals[colPtr[j]]), mindiag);
276  }
277 
278  FactorType L_save = m_L;
279 
280  m_shift = RealScalar(0);
281  if (mindiag <= RealScalar(0.)) m_shift = m_initialShift - mindiag;
282 
283  m_info = NumericalIssue;
284 
285  // Try to perform the incomplete factorization using the current shift
286  int iter = 0;
287  do {
288  // Apply the shift to the diagonal elements of the matrix
289  for (Index j = 0; j < n; j++) vals[colPtr[j]] += m_shift;
290 
291  // jki version of the Cholesky factorization
292  Index j = 0;
293  for (; j < n; ++j) {
294  // Left-looking factorization of the j-th column
295  // First, load the j-th column into col_vals
296  Scalar diag = vals[colPtr[j]]; // It is assumed that only the lower part is stored
297  col_nnz = 0;
298  for (Index i = colPtr[j] + 1; i < colPtr[j + 1]; i++) {
299  StorageIndex l = rowIdx[i];
300  col_vals(col_nnz) = vals[i];
301  col_irow(col_nnz) = l;
302  col_pattern(l) = col_nnz;
303  col_nnz++;
304  }
305  {
306  typename std::list<StorageIndex>::iterator k;
307  // Browse all previous columns that will update column j
308  for (k = listCol[j].begin(); k != listCol[j].end(); k++) {
309  Index jk = firstElt(*k); // First element to use in the column
310  eigen_internal_assert(rowIdx[jk] == j);
311  Scalar v_j_jk = numext::conj(vals[jk]);
312 
313  jk += 1;
314  for (Index i = jk; i < colPtr[*k + 1]; i++) {
315  StorageIndex l = rowIdx[i];
316  if (col_pattern[l] < 0) {
317  col_vals(col_nnz) = vals[i] * v_j_jk;
318  col_irow[col_nnz] = l;
319  col_pattern(l) = col_nnz;
320  col_nnz++;
321  } else
322  col_vals(col_pattern[l]) -= vals[i] * v_j_jk;
323  }
324  updateList(colPtr, rowIdx, vals, *k, jk, firstElt, listCol);
325  }
326  }
327 
328  // Scale the current column
329  if (numext::real(diag) <= 0) {
330  if (++iter >= 10) return;
331 
332  // increase shift
333  m_shift = numext::maxi(m_initialShift, RealScalar(2) * m_shift);
334  // restore m_L, col_pattern, and listCol
335  vals = Map<const VectorSx>(L_save.valuePtr(), nnz);
336  rowIdx = Map<const VectorIx>(L_save.innerIndexPtr(), nnz);
337  colPtr = Map<const VectorIx>(L_save.outerIndexPtr(), n + 1);
338  col_pattern.fill(-1);
339  for (Index i = 0; i < n; ++i) listCol[i].clear();
340 
341  break;
342  }
343 
344  RealScalar rdiag = sqrt(numext::real(diag));
345  vals[colPtr[j]] = rdiag;
346  for (Index k = 0; k < col_nnz; ++k) {
347  Index i = col_irow[k];
348  // Scale
349  col_vals(k) /= rdiag;
350  // Update the remaining diagonals with col_vals
351  vals[colPtr[i]] -= numext::abs2(col_vals(k));
352  }
353  // Select the largest p elements
354  // p is the original number of elements in the column (without the diagonal)
355  Index p = colPtr[j + 1] - colPtr[j] - 1;
356  Ref<VectorSx> cvals = col_vals.head(col_nnz);
357  Ref<VectorIx> cirow = col_irow.head(col_nnz);
358  internal::QuickSplit(cvals, cirow, p);
359  // Insert the largest p elements in the matrix
360  Index cpt = 0;
361  for (Index i = colPtr[j] + 1; i < colPtr[j + 1]; i++) {
362  vals[i] = col_vals(cpt);
363  rowIdx[i] = col_irow(cpt);
364  // restore col_pattern:
365  col_pattern(col_irow(cpt)) = -1;
366  cpt++;
367  }
368  // Get the first smallest row index and put it after the diagonal element
369  Index jk = colPtr(j) + 1;
370  updateList(colPtr, rowIdx, vals, j, jk, firstElt, listCol);
371  }
372 
373  if (j == n) {
374  m_factorizationIsOk = true;
375  m_info = Success;
376  }
377  } while (m_info != Success);
378 }
379 
380 template <typename Scalar, int UpLo_, typename OrderingType>
382  Ref<VectorIx> rowIdx, Ref<VectorSx> vals,
383  const Index& col, const Index& jk,
384  VectorIx& firstElt, VectorList& listCol) {
385  if (jk < colPtr(col + 1)) {
386  Index p = colPtr(col + 1) - jk;
387  Index minpos;
388  rowIdx.segment(jk, p).minCoeff(&minpos);
389  minpos += jk;
390  if (rowIdx(minpos) != rowIdx(jk)) {
391  // Swap
392  std::swap(rowIdx(jk), rowIdx(minpos));
393  std::swap(vals(jk), vals(minpos));
394  }
395  firstElt(col) = internal::convert_index<StorageIndex, Index>(jk);
396  listCol[rowIdx(jk)].push_back(internal::convert_index<StorageIndex, Index>(col));
397  }
398 }
399 
400 } // end namespace Eigen
401 
402 #endif
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
#define eigen_internal_assert(x)
Definition: Macros.h:916
#define EIGEN_NOEXCEPT
Definition: Macros.h:1267
#define EIGEN_CONSTEXPR
Definition: Macros.h:758
#define eigen_assert(x)
Definition: Macros.h:910
m col(1)
float * p
Definition: Tutorial_Map_using.cpp:9
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
Modified Incomplete Cholesky with dual threshold.
Definition: IncompleteCholesky.h:51
PermutationType m_perm
Definition: IncompleteCholesky.h:191
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: IncompleteCholesky.h:100
void updateList(Ref< const VectorIx > colPtr, Ref< VectorIx > rowIdx, Ref< VectorSx > vals, const Index &col, const Index &jk, VectorIx &firstElt, VectorList &listCol)
Definition: IncompleteCholesky.h:381
IncompleteCholesky(const MatrixType &matrix)
Definition: IncompleteCholesky.h:81
VectorRx m_scale
Definition: IncompleteCholesky.h:186
Matrix< StorageIndex, Dynamic, 1 > VectorIx
Definition: IncompleteCholesky.h:64
std::vector< std::list< StorageIndex > > VectorList
Definition: IncompleteCholesky.h:65
@ ColsAtCompileTime
Definition: IncompleteCholesky.h:67
@ MaxColsAtCompileTime
Definition: IncompleteCholesky.h:67
Matrix< Scalar, Dynamic, 1 > VectorSx
Definition: IncompleteCholesky.h:62
ComputationInfo m_info
Definition: IncompleteCholesky.h:190
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: IncompleteCholesky.h:87
@ UpLo
Definition: IncompleteCholesky.h:66
void factorize(const MatrixType &mat)
Performs the numerical factorization of the input matrix mat.
const VectorRx & scalingS() const
Definition: IncompleteCholesky.h:170
OrderingType_ OrderingType
Definition: IncompleteCholesky.h:58
IncompleteCholesky()
Definition: IncompleteCholesky.h:76
const PermutationType & permutationP() const
Definition: IncompleteCholesky.h:176
bool m_factorizationIsOk
Definition: IncompleteCholesky.h:189
void _solve_impl(const Rhs &b, Dest &x) const
Definition: IncompleteCholesky.h:150
Matrix< RealScalar, Dynamic, 1 > VectorRx
Definition: IncompleteCholesky.h:63
OrderingType::PermutationType PermutationType
Definition: IncompleteCholesky.h:59
bool m_analysisIsOk
Definition: IncompleteCholesky.h:188
FactorType m_L
Definition: IncompleteCholesky.h:185
SparseSolverBase< IncompleteCholesky< Scalar, UpLo_, OrderingType_ > > Base
Definition: IncompleteCholesky.h:53
bool m_isInitialized
Definition: SparseSolverBase.h:110
void compute(const MatrixType &mat)
Definition: IncompleteCholesky.h:143
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: IncompleteCholesky.h:90
RealScalar m_shift
Definition: IncompleteCholesky.h:192
const FactorType & matrixL() const
Definition: IncompleteCholesky.h:164
void analyzePattern(const MatrixType &mat)
Computes the fill reducing permutation vector using the sparsity pattern of mat.
Definition: IncompleteCholesky.h:112
RealScalar shift() const
Definition: IncompleteCholesky.h:182
RealScalar m_initialShift
Definition: IncompleteCholesky.h:187
PermutationType::StorageIndex StorageIndex
Definition: IncompleteCholesky.h:60
SparseMatrix< Scalar, ColMajor, StorageIndex > FactorType
Definition: IncompleteCholesky.h:61
NumTraits< Scalar >::Real RealScalar
Definition: IncompleteCholesky.h:57
void setInitialShift(RealScalar shift)
Set the initial shift parameter .
Definition: IncompleteCholesky.h:107
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:264
SparseSymmetricPermutationProduct< Derived, Upper|Lower > twistedBy(const PermutationMatrix< Dynamic, Dynamic, StorageIndex > &perm) const
Definition: SparseMatrixBase.h:325
const AdjointReturnType adjoint() const
Definition: SparseMatrixBase.h:360
Index cols() const
Definition: SparseMatrix.h:161
const Scalar * valuePtr() const
Definition: SparseMatrix.h:171
void resize(Index rows, Index cols)
Definition: SparseMatrix.h:734
Index rows() const
Definition: SparseMatrix.h:159
const StorageIndex * outerIndexPtr() const
Definition: SparseMatrix.h:189
const StorageIndex * innerIndexPtr() const
Definition: SparseMatrix.h:180
A base class for sparse solvers.
Definition: SparseSolverBase.h:67
bool m_isInitialized
Definition: SparseSolverBase.h:110
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85
float real
Definition: datatypes.h:10
#define min(a, b)
Definition: datatypes.h:22
ComputationInfo
Definition: Constants.h:438
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440
EIGEN_BLAS_FUNC() swap(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
Definition: level1_impl.h:117
const char const char const char * diag
Definition: level2_impl.h:86
char char char int int * k
Definition: level2_impl.h:374
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
@ Rhs
Definition: TensorContractionMapper.h:20
Index QuickSplit(VectorV &row, VectorI &ind, Index ncut)
Definition: IncompleteLUT.h:31
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Definition: MathFunctions.h:926
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
Definition: MathFunctions.h:920
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition: MathFunctions.h:1102
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
const int Dynamic
Definition: Constants.h:25
bool pinv
Definition: jeffery_hamel.cc:86
list x
Definition: plotDoE.py:28
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2