Eigen  3.2.91
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 
15 namespace Eigen {
16 
17 namespace internal {
18 
28 template <typename VectorV, typename VectorI>
29 Index QuickSplit(VectorV &row, VectorI &ind, Index ncut)
30 {
31  typedef typename VectorV::RealScalar RealScalar;
32  using std::swap;
33  using std::abs;
34  Index mid;
35  Index n = row.size(); /* length of the vector */
36  Index first, last ;
37 
38  ncut--; /* to fit the zero-based indices */
39  first = 0;
40  last = n-1;
41  if (ncut < first || ncut > last ) return 0;
42 
43  do {
44  mid = first;
45  RealScalar abskey = abs(row(mid));
46  for (Index j = first + 1; j <= last; j++) {
47  if ( abs(row(j)) > abskey) {
48  ++mid;
49  swap(row(mid), row(j));
50  swap(ind(mid), ind(j));
51  }
52  }
53  /* Interchange for the pivot element */
54  swap(row(mid), row(first));
55  swap(ind(mid), ind(first));
56 
57  if (mid > ncut) last = mid - 1;
58  else if (mid < ncut ) first = mid + 1;
59  } while (mid != ncut );
60 
61  return 0; /* mid is equal to ncut */
62 }
63 
64 }// end namespace internal
65 
96 template <typename _Scalar, typename _StorageIndex = int>
97 class IncompleteLUT : public SparseSolverBase<IncompleteLUT<_Scalar, _StorageIndex> >
98 {
99  protected:
101  using Base::m_isInitialized;
102  public:
103  typedef _Scalar Scalar;
104  typedef _StorageIndex StorageIndex;
105  typedef typename NumTraits<Scalar>::Real RealScalar;
109 
110  public:
111 
112  // this typedef is only to export the scalar type and compile-time dimensions to solve_retval
114 
115  IncompleteLUT()
116  : m_droptol(NumTraits<Scalar>::dummy_precision()), m_fillfactor(10),
117  m_analysisIsOk(false), m_factorizationIsOk(false)
118  {}
119 
120  template<typename MatrixType>
121  explicit IncompleteLUT(const MatrixType& mat, const RealScalar& droptol=NumTraits<Scalar>::dummy_precision(), int fillfactor = 10)
122  : m_droptol(droptol),m_fillfactor(fillfactor),
123  m_analysisIsOk(false),m_factorizationIsOk(false)
124  {
125  eigen_assert(fillfactor != 0);
126  compute(mat);
127  }
128 
129  Index rows() const { return m_lu.rows(); }
130 
131  Index cols() const { return m_lu.cols(); }
132 
139  {
140  eigen_assert(m_isInitialized && "IncompleteLUT is not initialized.");
141  return m_info;
142  }
143 
144  template<typename MatrixType>
145  void analyzePattern(const MatrixType& amat);
146 
147  template<typename MatrixType>
148  void factorize(const MatrixType& amat);
149 
155  template<typename MatrixType>
156  IncompleteLUT& compute(const MatrixType& amat)
157  {
158  analyzePattern(amat);
159  factorize(amat);
160  return *this;
161  }
162 
163  void setDroptol(const RealScalar& droptol);
164  void setFillfactor(int fillfactor);
165 
166  template<typename Rhs, typename Dest>
167  void _solve_impl(const Rhs& b, Dest& x) const
168  {
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:
176 
178  struct keep_diag {
179  inline bool operator() (const Index& row, const Index& col, const Scalar&) const
180  {
181  return row!=col;
182  }
183  };
184 
185 protected:
186 
187  FactorType m_lu;
188  RealScalar m_droptol;
189  int m_fillfactor;
190  bool m_analysisIsOk;
191  bool m_factorizationIsOk;
192  ComputationInfo m_info;
193  PermutationMatrix<Dynamic,Dynamic,StorageIndex> m_P; // Fill-reducing permutation
194  PermutationMatrix<Dynamic,Dynamic,StorageIndex> m_Pinv; // Inverse permutation
195 };
196 
201 template<typename Scalar, typename StorageIndex>
202 void IncompleteLUT<Scalar,StorageIndex>::setDroptol(const RealScalar& droptol)
203 {
204  this->m_droptol = droptol;
205 }
206 
211 template<typename Scalar, typename StorageIndex>
213 {
214  this->m_fillfactor = fillfactor;
215 }
216 
217 template <typename Scalar, typename StorageIndex>
218 template<typename _MatrixType>
219 void IncompleteLUT<Scalar,StorageIndex>::analyzePattern(const _MatrixType& amat)
220 {
221  // Compute the Fill-reducing permutation
223  SparseMatrix<Scalar,ColMajor, StorageIndex> mat2 = amat.transpose();
224  // Symmetrize the pattern
225  // FIXME for a matrix with nearly symmetric pattern, mat2+mat1 is the appropriate choice.
226  // on the other hand for a really non-symmetric pattern, mat2*mat1 should be prefered...
228  AtA.prune(keep_diag());
229  internal::minimum_degree_ordering<Scalar, StorageIndex>(AtA, m_P); // Then compute the AMD ordering...
230 
231  m_Pinv = m_P.inverse(); // ... and the inverse permutation
232 
233  m_analysisIsOk = true;
234  m_factorizationIsOk = false;
235  m_isInitialized = true;
236 }
237 
238 template <typename Scalar, typename StorageIndex>
239 template<typename _MatrixType>
240 void IncompleteLUT<Scalar,StorageIndex>::factorize(const _MatrixType& amat)
241 {
242  using std::sqrt;
243  using std::swap;
244  using std::abs;
245  using internal::convert_index;
246 
247  eigen_assert((amat.rows() == amat.cols()) && "The factorization should be done on a square matrix");
248  Index n = amat.cols(); // Size of the matrix
249  m_lu.resize(n,n);
250  // Declare Working vectors and variables
251  Vector u(n) ; // real values of the row -- maximum size is n --
252  VectorI ju(n); // column position of the values in u -- maximum size is n
253  VectorI jr(n); // Indicate the position of the nonzero elements in the vector u -- A zero location is indicated by -1
254 
255  // Apply the fill-reducing permutation
256  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
257  SparseMatrix<Scalar,RowMajor, StorageIndex> mat;
258  mat = amat.twistedBy(m_Pinv);
259 
260  // Initialization
261  jr.fill(-1);
262  ju.fill(0);
263  u.fill(0);
264 
265  // number of largest elements to keep in each row:
266  Index fill_in = (amat.nonZeros()*m_fillfactor)/n + 1;
267  if (fill_in > n) fill_in = n;
268 
269  // number of largest nonzero elements to keep in the L and the U part of the current row:
270  Index nnzL = fill_in/2;
271  Index nnzU = nnzL;
272  m_lu.reserve(n * (nnzL + nnzU + 1));
273 
274  // global loop over the rows of the sparse matrix
275  for (Index ii = 0; ii < n; ii++)
276  {
277  // 1 - copy the lower and the upper part of the row i of mat in the working vector u
278 
279  Index sizeu = 1; // number of nonzero elements in the upper part of the current row
280  Index sizel = 0; // number of nonzero elements in the lower part of the current row
281  ju(ii) = convert_index<StorageIndex>(ii);
282  u(ii) = 0;
283  jr(ii) = convert_index<StorageIndex>(ii);
284  RealScalar rownorm = 0;
285 
286  typename FactorType::InnerIterator j_it(mat, ii); // Iterate through the current row ii
287  for (; j_it; ++j_it)
288  {
289  Index k = j_it.index();
290  if (k < ii)
291  {
292  // copy the lower part
293  ju(sizel) = convert_index<StorageIndex>(k);
294  u(sizel) = j_it.value();
295  jr(k) = convert_index<StorageIndex>(sizel);
296  ++sizel;
297  }
298  else if (k == ii)
299  {
300  u(ii) = j_it.value();
301  }
302  else
303  {
304  // copy the upper part
305  Index jpos = ii + sizeu;
306  ju(jpos) = convert_index<StorageIndex>(k);
307  u(jpos) = j_it.value();
308  jr(k) = convert_index<StorageIndex>(jpos);
309  ++sizeu;
310  }
311  rownorm += numext::abs2(j_it.value());
312  }
313 
314  // 2 - detect possible zero row
315  if(rownorm==0)
316  {
317  m_info = NumericalIssue;
318  return;
319  }
320  // Take the 2-norm of the current row as a relative tolerance
321  rownorm = sqrt(rownorm);
322 
323  // 3 - eliminate the previous nonzero rows
324  Index jj = 0;
325  Index len = 0;
326  while (jj < sizel)
327  {
328  // In order to eliminate in the correct order,
329  // we must select first the smallest column index among ju(jj:sizel)
330  Index k;
331  Index minrow = ju.segment(jj,sizel-jj).minCoeff(&k); // k is relative to the segment
332  k += jj;
333  if (minrow != ju(jj))
334  {
335  // swap the two locations
336  Index j = ju(jj);
337  swap(ju(jj), ju(k));
338  jr(minrow) = convert_index<StorageIndex>(jj);
339  jr(j) = convert_index<StorageIndex>(k);
340  swap(u(jj), u(k));
341  }
342  // Reset this location
343  jr(minrow) = -1;
344 
345  // Start elimination
346  typename FactorType::InnerIterator ki_it(m_lu, minrow);
347  while (ki_it && ki_it.index() < minrow) ++ki_it;
348  eigen_internal_assert(ki_it && ki_it.col()==minrow);
349  Scalar fact = u(jj) / ki_it.value();
350 
351  // drop too small elements
352  if(abs(fact) <= m_droptol)
353  {
354  jj++;
355  continue;
356  }
357 
358  // linear combination of the current row ii and the row minrow
359  ++ki_it;
360  for (; ki_it; ++ki_it)
361  {
362  Scalar prod = fact * ki_it.value();
363  Index j = ki_it.index();
364  Index jpos = jr(j);
365  if (jpos == -1) // fill-in element
366  {
367  Index newpos;
368  if (j >= ii) // dealing with the upper part
369  {
370  newpos = ii + sizeu;
371  sizeu++;
372  eigen_internal_assert(sizeu<=n);
373  }
374  else // dealing with the lower part
375  {
376  newpos = sizel;
377  sizel++;
378  eigen_internal_assert(sizel<=ii);
379  }
380  ju(newpos) = convert_index<StorageIndex>(j);
381  u(newpos) = -prod;
382  jr(j) = convert_index<StorageIndex>(newpos);
383  }
384  else
385  u(jpos) -= prod;
386  }
387  // store the pivot element
388  u(len) = fact;
389  ju(len) = convert_index<StorageIndex>(minrow);
390  ++len;
391 
392  jj++;
393  } // end of the elimination on the row ii
394 
395  // reset the upper part of the pointer jr to zero
396  for(Index k = 0; k <sizeu; k++) jr(ju(ii+k)) = -1;
397 
398  // 4 - partially sort and insert the elements in the m_lu matrix
399 
400  // sort the L-part of the row
401  sizel = len;
402  len = (std::min)(sizel, nnzL);
403  typename Vector::SegmentReturnType ul(u.segment(0, sizel));
404  typename VectorI::SegmentReturnType jul(ju.segment(0, sizel));
405  internal::QuickSplit(ul, jul, len);
406 
407  // store the largest m_fill elements of the L part
408  m_lu.startVec(ii);
409  for(Index k = 0; k < len; k++)
410  m_lu.insertBackByOuterInnerUnordered(ii,ju(k)) = u(k);
411 
412  // store the diagonal element
413  // apply a shifting rule to avoid zero pivots (we are doing an incomplete factorization)
414  if (u(ii) == Scalar(0))
415  u(ii) = sqrt(m_droptol) * rownorm;
416  m_lu.insertBackByOuterInnerUnordered(ii, ii) = u(ii);
417 
418  // sort the U-part of the row
419  // apply the dropping rule first
420  len = 0;
421  for(Index k = 1; k < sizeu; k++)
422  {
423  if(abs(u(ii+k)) > m_droptol * rownorm )
424  {
425  ++len;
426  u(ii + len) = u(ii + k);
427  ju(ii + len) = ju(ii + k);
428  }
429  }
430  sizeu = len + 1; // +1 to take into account the diagonal element
431  len = (std::min)(sizeu, nnzU);
432  typename Vector::SegmentReturnType uu(u.segment(ii+1, sizeu-1));
433  typename VectorI::SegmentReturnType juu(ju.segment(ii+1, sizeu-1));
434  internal::QuickSplit(uu, juu, len);
435 
436  // store the largest elements of the U part
437  for(Index k = ii + 1; k < ii + len; k++)
438  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:202
Index cols() const
Definition: SparseMatrix.h:132
A base class for sparse solvers.
Definition: SparseSolverBase.h:53
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: IncompleteLUT.h:138
const Solve< IncompleteLUT< _Scalar, _StorageIndex >, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SparseSolverBase.h:74
Definition: LDLT.h:16
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:107
Definition: IncompleteLUT.h:178
void prune(const Scalar &reference, const RealScalar &epsilon=NumTraits< RealScalar >::dummy_precision())
Definition: SparseMatrix.h:495
Definition: Constants.h:426
Incomplete LU factorization with dual-threshold strategy.
Definition: IncompleteLUT.h:97
Definition: Constants.h:424
Definition: Eigen_Colamd.h:54
IncompleteLUT & compute(const MatrixType &amat)
Definition: IncompleteLUT.h:156
ComputationInfo
Definition: Constants.h:422
Index rows() const
Definition: SparseMatrix.h:130
void setFillfactor(int fillfactor)
Definition: IncompleteLUT.h:212