IncompleteLUT.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) 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;
114 
116 
117  public:
119  : m_droptol(NumTraits<Scalar>::dummy_precision()),
120  m_fillfactor(10),
121  m_analysisIsOk(false),
122  m_factorizationIsOk(false) {}
123 
124  template <typename MatrixType>
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 
135 
142  eigen_assert(m_isInitialized && "IncompleteLUT is not initialized.");
143  return m_info;
144  }
145 
146  template <typename MatrixType>
147  void analyzePattern(const MatrixType& amat);
148 
149  template <typename MatrixType>
150  void factorize(const MatrixType& amat);
151 
157  template <typename MatrixType>
159  analyzePattern(amat);
160  factorize(amat);
161  return *this;
162  }
163 
164  void setDroptol(const RealScalar& droptol);
165  void setFillfactor(int fillfactor);
166 
167  template <typename Rhs, typename Dest>
168  void _solve_impl(const Rhs& b, Dest& x) const {
169  x = m_Pinv * b;
170  x = m_lu.template triangularView<UnitLower>().solve(x);
171  x = m_lu.template triangularView<Upper>().solve(x);
172  x = m_P * x;
173  }
174 
175  protected:
177  struct keep_diag {
178  inline bool operator()(const Index& row, const Index& col, const Scalar&) const { return row != col; }
179  };
180 
181  protected:
190 };
191 
196 template <typename Scalar, typename StorageIndex>
198  this->m_droptol = droptol;
199 }
200 
205 template <typename Scalar, typename StorageIndex>
207  this->m_fillfactor = fillfactor;
208 }
209 
210 template <typename Scalar, typename StorageIndex>
211 template <typename MatrixType_>
213  // Compute the Fill-reducing permutation
214  // Since ILUT does not perform any numerical pivoting,
215  // it is highly preferable to keep the diagonal through symmetric permutations.
216  // To this end, let's symmetrize the pattern and perform AMD on it.
219  // FIXME for a matrix with nearly symmetric pattern, mat2+mat1 is the appropriate choice.
220  // on the other hand for a really non-symmetric pattern, mat2*mat1 should be preferred...
222  AMDOrdering<StorageIndex> ordering;
223  ordering(AtA, m_P);
224  m_Pinv = m_P.inverse(); // cache the inverse permutation
225  m_analysisIsOk = true;
226  m_factorizationIsOk = false;
227  m_isInitialized = true;
228 }
229 
230 template <typename Scalar, typename StorageIndex>
231 template <typename MatrixType_>
232 void IncompleteLUT<Scalar, StorageIndex>::factorize(const MatrixType_& amat) {
234  using std::abs;
235  using std::sqrt;
236  using std::swap;
237 
238  eigen_assert((amat.rows() == amat.cols()) && "The factorization should be done on a square matrix");
239  Index n = amat.cols(); // Size of the matrix
240  m_lu.resize(n, n);
241  // Declare Working vectors and variables
242  Vector u(n); // real values of the row -- maximum size is n --
243  VectorI ju(n); // column position of the values in u -- maximum size is n
244  VectorI jr(n); // Indicate the position of the nonzero elements in the vector u -- A zero location is indicated by -1
245 
246  // Apply the fill-reducing permutation
247  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
249  mat = amat.twistedBy(m_Pinv);
250 
251  // Initialization
252  jr.fill(-1);
253  ju.fill(0);
254  u.fill(0);
255 
256  // number of largest elements to keep in each row:
257  Index fill_in = (amat.nonZeros() * m_fillfactor) / n + 1;
258  if (fill_in > n) fill_in = n;
259 
260  // number of largest nonzero elements to keep in the L and the U part of the current row:
261  Index nnzL = fill_in / 2;
262  Index nnzU = nnzL;
263  m_lu.reserve(n * (nnzL + nnzU + 1));
264 
265  // global loop over the rows of the sparse matrix
266  for (Index ii = 0; ii < n; ii++) {
267  // 1 - copy the lower and the upper part of the row i of mat in the working vector u
268 
269  Index sizeu = 1; // number of nonzero elements in the upper part of the current row
270  Index sizel = 0; // number of nonzero elements in the lower part of the current row
271  ju(ii) = convert_index<StorageIndex>(ii);
272  u(ii) = 0;
273  jr(ii) = convert_index<StorageIndex>(ii);
274  RealScalar rownorm = 0;
275 
276  typename FactorType::InnerIterator j_it(mat, ii); // Iterate through the current row ii
277  for (; j_it; ++j_it) {
278  Index k = j_it.index();
279  if (k < ii) {
280  // copy the lower part
281  ju(sizel) = convert_index<StorageIndex>(k);
282  u(sizel) = j_it.value();
283  jr(k) = convert_index<StorageIndex>(sizel);
284  ++sizel;
285  } else if (k == ii) {
286  u(ii) = j_it.value();
287  } else {
288  // copy the upper part
289  Index jpos = ii + sizeu;
290  ju(jpos) = convert_index<StorageIndex>(k);
291  u(jpos) = j_it.value();
292  jr(k) = convert_index<StorageIndex>(jpos);
293  ++sizeu;
294  }
295  rownorm += numext::abs2(j_it.value());
296  }
297 
298  // 2 - detect possible zero row
299  if (rownorm == 0) {
300  m_info = NumericalIssue;
301  return;
302  }
303  // Take the 2-norm of the current row as a relative tolerance
304  rownorm = sqrt(rownorm);
305 
306  // 3 - eliminate the previous nonzero rows
307  Index jj = 0;
308  Index len = 0;
309  while (jj < sizel) {
310  // In order to eliminate in the correct order,
311  // we must select first the smallest column index among ju(jj:sizel)
312  Index k;
313  Index minrow = ju.segment(jj, sizel - jj).minCoeff(&k); // k is relative to the segment
314  k += jj;
315  if (minrow != ju(jj)) {
316  // swap the two locations
317  Index j = ju(jj);
318  swap(ju(jj), ju(k));
319  jr(minrow) = convert_index<StorageIndex>(jj);
320  jr(j) = convert_index<StorageIndex>(k);
321  swap(u(jj), u(k));
322  }
323  // Reset this location
324  jr(minrow) = -1;
325 
326  // Start elimination
327  typename FactorType::InnerIterator ki_it(m_lu, minrow);
328  while (ki_it && ki_it.index() < minrow) ++ki_it;
329  eigen_internal_assert(ki_it && ki_it.col() == minrow);
330  Scalar fact = u(jj) / ki_it.value();
331 
332  // drop too small elements
333  if (abs(fact) <= m_droptol) {
334  jj++;
335  continue;
336  }
337 
338  // linear combination of the current row ii and the row minrow
339  ++ki_it;
340  for (; ki_it; ++ki_it) {
341  Scalar prod = fact * ki_it.value();
342  Index j = ki_it.index();
343  Index jpos = jr(j);
344  if (jpos == -1) // fill-in element
345  {
346  Index newpos;
347  if (j >= ii) // dealing with the upper part
348  {
349  newpos = ii + sizeu;
350  sizeu++;
351  eigen_internal_assert(sizeu <= n);
352  } else // dealing with the lower part
353  {
354  newpos = sizel;
355  sizel++;
356  eigen_internal_assert(sizel <= ii);
357  }
358  ju(newpos) = convert_index<StorageIndex>(j);
359  u(newpos) = -prod;
360  jr(j) = convert_index<StorageIndex>(newpos);
361  } else
362  u(jpos) -= prod;
363  }
364  // store the pivot element
365  u(len) = fact;
366  ju(len) = convert_index<StorageIndex>(minrow);
367  ++len;
368 
369  jj++;
370  } // end of the elimination on the row ii
371 
372  // reset the upper part of the pointer jr to zero
373  for (Index k = 0; k < sizeu; k++) jr(ju(ii + k)) = -1;
374 
375  // 4 - partially sort and insert the elements in the m_lu matrix
376 
377  // sort the L-part of the row
378  sizel = len;
379  len = (std::min)(sizel, nnzL);
380  typename Vector::SegmentReturnType ul(u.segment(0, sizel));
381  typename VectorI::SegmentReturnType jul(ju.segment(0, sizel));
382  internal::QuickSplit(ul, jul, len);
383 
384  // store the largest m_fill elements of the L part
385  m_lu.startVec(ii);
386  for (Index k = 0; k < len; k++) m_lu.insertBackByOuterInnerUnordered(ii, ju(k)) = u(k);
387 
388  // store the diagonal element
389  // apply a shifting rule to avoid zero pivots (we are doing an incomplete factorization)
390  if (u(ii) == Scalar(0)) u(ii) = sqrt(m_droptol) * rownorm;
391  m_lu.insertBackByOuterInnerUnordered(ii, ii) = u(ii);
392 
393  // sort the U-part of the row
394  // apply the dropping rule first
395  len = 0;
396  for (Index k = 1; k < sizeu; k++) {
397  if (abs(u(ii + k)) > m_droptol * rownorm) {
398  ++len;
399  u(ii + len) = u(ii + k);
400  ju(ii + len) = ju(ii + k);
401  }
402  }
403  sizeu = len + 1; // +1 to take into account the diagonal element
404  len = (std::min)(sizeu, nnzU);
405  typename Vector::SegmentReturnType uu(u.segment(ii + 1, sizeu - 1));
406  typename VectorI::SegmentReturnType juu(ju.segment(ii + 1, sizeu - 1));
407  internal::QuickSplit(uu, juu, len);
408 
409  // store the largest elements of the U part
410  for (Index k = ii + 1; k < ii + len; k++) m_lu.insertBackByOuterInnerUnordered(ii, ju(k)) = u(k);
411  }
412  m_lu.finalize();
413  m_lu.makeCompressed();
414 
415  m_factorizationIsOk = true;
416  m_info = Success;
417 }
418 
419 } // end namespace Eigen
420 
421 #endif // EIGEN_INCOMPLETE_LUT_H
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Eigen::SparseMatrix< double > mat
Definition: EigenUnitTest.cpp:10
#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)
m row(1)
std::vector< int > ind
Definition: Slicing_stdvector_cxx11.cpp:1
MatrixXd mat1(size, size)
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
Definition: Ordering.h:50
Incomplete LU factorization with dual-threshold strategy.
Definition: IncompleteLUT.h:102
ComputationInfo m_info
Definition: IncompleteLUT.h:187
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_Pinv
Definition: IncompleteLUT.h:189
RealScalar m_droptol
Definition: IncompleteLUT.h:183
void analyzePattern(const MatrixType &amat)
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: IncompleteLUT.h:132
void setFillfactor(int fillfactor)
Definition: IncompleteLUT.h:206
void factorize(const MatrixType &amat)
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: IncompleteLUT.h:134
IncompleteLUT()
Definition: IncompleteLUT.h:118
FactorType m_lu
Definition: IncompleteLUT.h:182
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_P
Definition: IncompleteLUT.h:188
NumTraits< Scalar >::Real RealScalar
Definition: IncompleteLUT.h:110
@ ColsAtCompileTime
Definition: IncompleteLUT.h:115
@ MaxColsAtCompileTime
Definition: IncompleteLUT.h:115
IncompleteLUT & compute(const MatrixType &amat)
Definition: IncompleteLUT.h:158
void _solve_impl(const Rhs &b, Dest &x) const
Definition: IncompleteLUT.h:168
int m_fillfactor
Definition: IncompleteLUT.h:184
void setDroptol(const RealScalar &droptol)
Definition: IncompleteLUT.h:197
bool m_factorizationIsOk
Definition: IncompleteLUT.h:186
SparseSolverBase< IncompleteLUT > Base
Definition: IncompleteLUT.h:104
SparseMatrix< Scalar, RowMajor, StorageIndex > FactorType
Definition: IncompleteLUT.h:113
StorageIndex_ StorageIndex
Definition: IncompleteLUT.h:109
Matrix< StorageIndex, Dynamic, 1 > VectorI
Definition: IncompleteLUT.h:112
bool m_analysisIsOk
Definition: IncompleteLUT.h:185
bool m_isInitialized
Definition: SparseSolverBase.h:110
Matrix< Scalar, Dynamic, 1 > Vector
Definition: IncompleteLUT.h:111
IncompleteLUT(const MatrixType &mat, const RealScalar &droptol=NumTraits< Scalar >::dummy_precision(), int fillfactor=10)
Definition: IncompleteLUT.h:125
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: IncompleteLUT.h:141
Scalar_ Scalar
Definition: IncompleteLUT.h:108
SparseSymmetricPermutationProduct< Derived, Upper|Lower > twistedBy(const PermutationMatrix< Dynamic, Dynamic, StorageIndex > &perm) const
Definition: SparseMatrixBase.h:325
TransposeReturnType transpose()
Definition: SparseMatrixBase.h:358
Index cols() const
Definition: SparseMatrix.h:161
Index rows() const
Definition: SparseMatrix.h:159
Base::InnerIterator InnerIterator
Definition: SparseMatrix.h:138
A base class for sparse solvers.
Definition: SparseSolverBase.h:67
bool m_isInitialized
Definition: SparseSolverBase.h:110
#define min(a, b)
Definition: datatypes.h:22
static constexpr const last_t last
Definition: IndexedViewHelper.h:48
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
char char char int int * k
Definition: level2_impl.h:374
@ Rhs
Definition: TensorContractionMapper.h:20
EIGEN_DEVICE_FUNC IndexDest convert_index(const IndexSrc &idx)
Definition: XprHelper.h:63
Index QuickSplit(VectorV &row, VectorI &ind, Index ncut)
Definition: IncompleteLUT.h:31
void swap(scoped_array< T > &a, scoped_array< T > &b)
Definition: Memory.h:734
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
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE 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:655
const Product< Lhs, Rhs > prod(const Lhs &lhs, const Rhs &rhs)
Definition: evaluators.cpp:7
Definition: Eigen_Colamd.h:49
list x
Definition: plotDoE.py:28
Definition: IncompleteLUT.h:177
bool operator()(const Index &row, const Index &col, const Scalar &) const
Definition: IncompleteLUT.h:178
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