$darkmode
Eigen  5.0.1-dev
IncompleteLUT.h
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) 2014 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_LUT_H
12 #define EIGEN_INCOMPLETE_LUT_H
13 
14 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 namespace internal {
20 
30 template <typename VectorV, typename VectorI>
31 Index QuickSplit(VectorV& row, VectorI& ind, Index ncut) {
32  typedef typename VectorV::RealScalar RealScalar;
33  using std::abs;
34  using std::swap;
35  Index mid;
36  Index n = row.size(); /* length of the vector */
37  Index first, last;
38 
39  ncut--; /* to fit the zero-based indices */
40  first = 0;
41  last = n - 1;
42  if (ncut < first || ncut > last) return 0;
43 
44  do {
45  mid = first;
46  RealScalar abskey = abs(row(mid));
47  for (Index j = first + 1; j <= last; j++) {
48  if (abs(row(j)) > abskey) {
49  ++mid;
50  swap(row(mid), row(j));
51  swap(ind(mid), ind(j));
52  }
53  }
54  /* Interchange for the pivot element */
55  swap(row(mid), row(first));
56  swap(ind(mid), ind(first));
57 
58  if (mid > ncut)
59  last = mid - 1;
60  else if (mid < ncut)
61  first = mid + 1;
62  } while (mid != ncut);
63 
64  return 0; /* mid is equal to ncut */
65 }
66 
67 } // end namespace internal
68 
101 template <typename Scalar_, typename StorageIndex_ = int>
102 class IncompleteLUT : public SparseSolverBase<IncompleteLUT<Scalar_, StorageIndex_> > {
103  protected:
105  using Base::m_isInitialized;
106 
107  public:
108  typedef Scalar_ Scalar;
109  typedef StorageIndex_ StorageIndex;
110  typedef typename NumTraits<Scalar>::Real RealScalar;
114 
115  enum { ColsAtCompileTime = Dynamic, MaxColsAtCompileTime = Dynamic };
116 
117  public:
118  IncompleteLUT()
119  : m_droptol(NumTraits<Scalar>::dummy_precision()),
120  m_fillfactor(10),
121  m_analysisIsOk(false),
122  m_factorizationIsOk(false) {}
123 
124  template <typename MatrixType>
125  explicit IncompleteLUT(const MatrixType& mat, const RealScalar& droptol = NumTraits<Scalar>::dummy_precision(),
126  int fillfactor = 10)
127  : m_droptol(droptol), m_fillfactor(fillfactor), m_analysisIsOk(false), m_factorizationIsOk(false) {
128  eigen_assert(fillfactor != 0);
129  compute(mat);
130  }
131 
133  const FactorType matrixL() const;
134 
136  const FactorType matrixU() const;
137 
138  constexpr Index rows() const noexcept { return m_lu.rows(); }
139 
140  constexpr Index cols() const noexcept { return m_lu.cols(); }
141 
148  eigen_assert(m_isInitialized && "IncompleteLUT is not initialized.");
149  return m_info;
150  }
151 
152  template <typename MatrixType>
153  void analyzePattern(const MatrixType& amat);
154 
155  template <typename MatrixType>
156  void factorize(const MatrixType& amat);
157 
163  template <typename MatrixType>
164  IncompleteLUT& compute(const MatrixType& amat) {
165  analyzePattern(amat);
166  factorize(amat);
167  return *this;
168  }
169 
170  void setDroptol(const RealScalar& droptol);
171  void setFillfactor(int fillfactor);
172 
173  template <typename Rhs, typename Dest>
174  void _solve_impl(const Rhs& b, Dest& x) const {
175  x = m_Pinv * b;
176  x = m_lu.template triangularView<UnitLower>().solve(x);
177  x = m_lu.template triangularView<Upper>().solve(x);
178  x = m_P * x;
179  }
180 
181  protected:
183  struct keep_diag {
184  inline bool operator()(const Index& row, const Index& col, const Scalar&) const { return row != col; }
185  };
186 
187  protected:
188  FactorType m_lu;
189  RealScalar m_droptol;
190  int m_fillfactor;
191  bool m_analysisIsOk;
192  bool m_factorizationIsOk;
193  ComputationInfo m_info;
194  PermutationMatrix<Dynamic, Dynamic, StorageIndex> m_P; // Fill-reducing permutation
195  PermutationMatrix<Dynamic, Dynamic, StorageIndex> m_Pinv; // Inverse permutation
196 };
197 
202 template <typename Scalar, typename StorageIndex>
203 void IncompleteLUT<Scalar, StorageIndex>::setDroptol(const RealScalar& droptol) {
204  this->m_droptol = droptol;
205 }
206 
211 template <typename Scalar, typename StorageIndex>
213  this->m_fillfactor = fillfactor;
214 }
215 
221 template <typename Scalar, typename StorageIndex>
223  eigen_assert(m_factorizationIsOk && "factorize() should be called first");
224  return m_lu.template triangularView<UnitLower>();
225 }
226 
232 template <typename Scalar, typename StorageIndex>
234  eigen_assert(m_factorizationIsOk && "Factorization must be computed first.");
235  return m_lu.template triangularView<Upper>();
236 }
237 
238 template <typename Scalar, typename StorageIndex>
239 template <typename MatrixType_>
240 void IncompleteLUT<Scalar, StorageIndex>::analyzePattern(const MatrixType_& amat) {
241  // Compute the Fill-reducing permutation
242  // Since ILUT does not perform any numerical pivoting,
243  // it is highly preferable to keep the diagonal through symmetric permutations.
244  // To this end, let's symmetrize the pattern and perform AMD on it.
246  SparseMatrix<Scalar, ColMajor, StorageIndex> mat2 = amat.transpose();
247  // FIXME for a matrix with nearly symmetric pattern, mat2+mat1 is the appropriate choice.
248  // on the other hand for a really non-symmetric pattern, mat2*mat1 should be preferred...
250  AMDOrdering<StorageIndex> ordering;
251  ordering(AtA, m_P);
252  m_Pinv = m_P.inverse(); // cache the inverse permutation
253  m_analysisIsOk = true;
254  m_factorizationIsOk = false;
255  m_isInitialized = true;
256 }
257 
258 template <typename Scalar, typename StorageIndex>
259 template <typename MatrixType_>
260 void IncompleteLUT<Scalar, StorageIndex>::factorize(const MatrixType_& amat) {
261  using internal::convert_index;
262  using std::abs;
263  using std::sqrt;
264  using std::swap;
265 
266  eigen_assert((amat.rows() == amat.cols()) && "The factorization should be done on a square matrix");
267  Index n = amat.cols(); // Size of the matrix
268  m_lu.resize(n, n);
269  // Declare Working vectors and variables
270  Vector u(n); // real values of the row -- maximum size is n --
271  VectorI ju(n); // column position of the values in u -- maximum size is n
272  VectorI jr(n); // Indicate the position of the nonzero elements in the vector u -- A zero location is indicated by -1
273 
274  // Apply the fill-reducing permutation
275  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
276  SparseMatrix<Scalar, RowMajor, StorageIndex> mat;
277  mat = amat.twistedBy(m_Pinv);
278 
279  // Initialization
280  jr.fill(-1);
281  ju.fill(0);
282  u.fill(0);
283 
284  // number of largest elements to keep in each row:
285  Index fill_in = (amat.nonZeros() * m_fillfactor) / n + 1;
286  if (fill_in > n) fill_in = n;
287 
288  // number of largest nonzero elements to keep in the L and the U part of the current row:
289  Index nnzL = fill_in / 2;
290  Index nnzU = nnzL;
291  m_lu.reserve(n * (nnzL + nnzU + 1));
292 
293  // global loop over the rows of the sparse matrix
294  for (Index ii = 0; ii < n; ii++) {
295  // 1 - copy the lower and the upper part of the row i of mat in the working vector u
296 
297  Index sizeu = 1; // number of nonzero elements in the upper part of the current row
298  Index sizel = 0; // number of nonzero elements in the lower part of the current row
299  ju(ii) = convert_index<StorageIndex>(ii);
300  u(ii) = 0;
301  jr(ii) = convert_index<StorageIndex>(ii);
302  RealScalar rownorm = 0;
303 
304  typename FactorType::InnerIterator j_it(mat, ii); // Iterate through the current row ii
305  for (; j_it; ++j_it) {
306  Index k = j_it.index();
307  if (k < ii) {
308  // copy the lower part
309  ju(sizel) = convert_index<StorageIndex>(k);
310  u(sizel) = j_it.value();
311  jr(k) = convert_index<StorageIndex>(sizel);
312  ++sizel;
313  } else if (k == ii) {
314  u(ii) = j_it.value();
315  } else {
316  // copy the upper part
317  Index jpos = ii + sizeu;
318  ju(jpos) = convert_index<StorageIndex>(k);
319  u(jpos) = j_it.value();
320  jr(k) = convert_index<StorageIndex>(jpos);
321  ++sizeu;
322  }
323  rownorm += numext::abs2(j_it.value());
324  }
325 
326  // 2 - detect possible zero row
327  if (rownorm == 0) {
328  m_info = NumericalIssue;
329  return;
330  }
331  // Take the 2-norm of the current row as a relative tolerance
332  rownorm = sqrt(rownorm);
333 
334  // 3 - eliminate the previous nonzero rows
335  Index jj = 0;
336  Index len = 0;
337  while (jj < sizel) {
338  // In order to eliminate in the correct order,
339  // we must select first the smallest column index among ju(jj:sizel)
340  Index k;
341  Index minrow = ju.segment(jj, sizel - jj).minCoeff(&k); // k is relative to the segment
342  k += jj;
343  if (minrow != ju(jj)) {
344  // swap the two locations
345  Index j = ju(jj);
346  swap(ju(jj), ju(k));
347  jr(minrow) = convert_index<StorageIndex>(jj);
348  jr(j) = convert_index<StorageIndex>(k);
349  swap(u(jj), u(k));
350  }
351  // Reset this location
352  jr(minrow) = -1;
353 
354  // Start elimination
355  typename FactorType::InnerIterator ki_it(m_lu, minrow);
356  while (ki_it && ki_it.index() < minrow) ++ki_it;
357  eigen_internal_assert(ki_it && ki_it.col() == minrow);
358  Scalar fact = u(jj) / ki_it.value();
359 
360  // drop too small elements
361  if (abs(fact) <= m_droptol) {
362  jj++;
363  continue;
364  }
365 
366  // linear combination of the current row ii and the row minrow
367  ++ki_it;
368  for (; ki_it; ++ki_it) {
369  Scalar prod = fact * ki_it.value();
370  Index j = ki_it.index();
371  Index jpos = jr(j);
372  if (jpos == -1) // fill-in element
373  {
374  Index newpos;
375  if (j >= ii) // dealing with the upper part
376  {
377  newpos = ii + sizeu;
378  sizeu++;
379  eigen_internal_assert(sizeu <= n);
380  } else // dealing with the lower part
381  {
382  newpos = sizel;
383  sizel++;
384  eigen_internal_assert(sizel <= ii);
385  }
386  ju(newpos) = convert_index<StorageIndex>(j);
387  u(newpos) = -prod;
388  jr(j) = convert_index<StorageIndex>(newpos);
389  } else
390  u(jpos) -= prod;
391  }
392  // store the pivot element
393  u(len) = fact;
394  ju(len) = convert_index<StorageIndex>(minrow);
395  ++len;
396 
397  jj++;
398  } // end of the elimination on the row ii
399 
400  // reset the upper part of the pointer jr to zero
401  for (Index k = 0; k < sizeu; k++) jr(ju(ii + k)) = -1;
402 
403  // 4 - partially sort and insert the elements in the m_lu matrix
404 
405  // sort the L-part of the row
406  sizel = len;
407  len = (std::min)(sizel, nnzL);
408  typename Vector::SegmentReturnType ul(u.segment(0, sizel));
409  typename VectorI::SegmentReturnType jul(ju.segment(0, sizel));
410  internal::QuickSplit(ul, jul, len);
411 
412  // store the largest m_fill elements of the L part
413  m_lu.startVec(ii);
414  for (Index k = 0; k < len; k++) m_lu.insertBackByOuterInnerUnordered(ii, ju(k)) = u(k);
415 
416  // store the diagonal element
417  // apply a shifting rule to avoid zero pivots (we are doing an incomplete factorization)
418  if (u(ii) == Scalar(0)) u(ii) = sqrt(m_droptol) * rownorm;
419  m_lu.insertBackByOuterInnerUnordered(ii, ii) = u(ii);
420 
421  // sort the U-part of the row
422  // apply the dropping rule first
423  len = 0;
424  for (Index k = 1; k < sizeu; k++) {
425  if (abs(u(ii + k)) > m_droptol * rownorm) {
426  ++len;
427  u(ii + len) = u(ii + k);
428  ju(ii + len) = ju(ii + k);
429  }
430  }
431  sizeu = len + 1; // +1 to take into account the diagonal element
432  len = (std::min)(sizeu, nnzU);
433  typename Vector::SegmentReturnType uu(u.segment(ii + 1, sizeu - 1));
434  typename VectorI::SegmentReturnType juu(ju.segment(ii + 1, sizeu - 1));
435  internal::QuickSplit(uu, juu, len);
436 
437  // store the largest elements of the U part
438  for (Index k = ii + 1; k < ii + len; k++) m_lu.insertBackByOuterInnerUnordered(ii, ju(k)) = u(k);
439  }
440  m_lu.finalize();
441  m_lu.makeCompressed();
442 
443  m_factorizationIsOk = true;
444  m_info = Success;
445 }
446 
447 } // end namespace Eigen
448 
449 #endif // EIGEN_INCOMPLETE_LUT_H
void setDroptol(const RealScalar &droptol)
Definition: IncompleteLUT.h:203
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Index rows() const
Definition: SparseMatrix.h:159
A base class for sparse solvers.
Definition: SparseSolverBase.h:67
const Solve< IncompleteLUT< Scalar_, StorageIndex_ >, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SparseSolverBase.h:84
Definition: Ordering.h:48
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
const FactorType matrixU() const
Extraction Method for U-Factor.
Definition: IncompleteLUT.h:233
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
Matrix< Type, Size, 1 > Vector
[c++11] Size×1 vector of type Type.
Definition: Matrix.h:522
Definition: IncompleteLUT.h:183
Definition: Constants.h:442
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
std::enable_if_t< std::is_base_of< DenseBase< std::decay_t< DerivedA > >, std::decay_t< DerivedA > >::value &&std::is_base_of< DenseBase< std::decay_t< DerivedB > >, std::decay_t< DerivedB > >::value, void > swap(DerivedA &&a, DerivedB &&b)
Definition: DenseBase.h:667
Incomplete LU factorization with dual-threshold strategy.
Definition: IncompleteLUT.h:102
Definition: Constants.h:440
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
Index cols() const
Definition: SparseMatrix.h:161
IncompleteLUT & compute(const MatrixType &amat)
Definition: IncompleteLUT.h:164
const int Dynamic
Definition: Constants.h:25
ComputationInfo
Definition: Constants.h:438
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: IncompleteLUT.h:147
static constexpr const last_t last
Definition: IndexedViewHelper.h:48
const FactorType matrixL() const
Extraction Method for L-Factor.
Definition: IncompleteLUT.h:222
void setFillfactor(int fillfactor)
Definition: IncompleteLUT.h:212