blob: cdcf709eb63e0c5785d314db3855c813b8f4e3f3 [file] [log] [blame]
Brian Silverman72890c22015-09-19 14:37:37 -04001// 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>
Austin Schuh189376f2018-12-20 22:11:15 +11005// Copyright (C) 2014 Gael Guennebaud <gael.guennebaud@inria.fr>
Brian Silverman72890c22015-09-19 14:37:37 -04006//
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
Austin Schuhc55b0172022-02-20 17:52:35 -080015namespace Eigen {
Brian Silverman72890c22015-09-19 14:37:37 -040016
17namespace internal {
Austin Schuhc55b0172022-02-20 17:52:35 -080018
Brian Silverman72890c22015-09-19 14:37:37 -040019/** \internal
Austin Schuhc55b0172022-02-20 17:52:35 -080020 * Compute a quick-sort split of a vector
Brian Silverman72890c22015-09-19 14:37:37 -040021 * On output, the vector row is permuted such that its elements satisfy
22 * abs(row(i)) >= abs(row(ncut)) if i<ncut
Austin Schuhc55b0172022-02-20 17:52:35 -080023 * abs(row(i)) <= abs(row(ncut)) if i>ncut
Brian Silverman72890c22015-09-19 14:37:37 -040024 * \param row The vector of values
25 * \param ind The array of index for the elements in @p row
26 * \param ncut The number of largest elements to keep
Austin Schuhc55b0172022-02-20 17:52:35 -080027 **/
Austin Schuh189376f2018-12-20 22:11:15 +110028template <typename VectorV, typename VectorI>
Brian Silverman72890c22015-09-19 14:37:37 -040029Index 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 ;
Austin Schuhc55b0172022-02-20 17:52:35 -080037
Brian Silverman72890c22015-09-19 14:37:37 -040038 ncut--; /* to fit the zero-based indices */
Austin Schuhc55b0172022-02-20 17:52:35 -080039 first = 0;
40 last = n-1;
Brian Silverman72890c22015-09-19 14:37:37 -040041 if (ncut < first || ncut > last ) return 0;
Austin Schuhc55b0172022-02-20 17:52:35 -080042
Brian Silverman72890c22015-09-19 14:37:37 -040043 do {
Austin Schuhc55b0172022-02-20 17:52:35 -080044 mid = first;
45 RealScalar abskey = abs(row(mid));
Brian Silverman72890c22015-09-19 14:37:37 -040046 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));
Austin Schuhc55b0172022-02-20 17:52:35 -080056
Brian Silverman72890c22015-09-19 14:37:37 -040057 if (mid > ncut) last = mid - 1;
Austin Schuhc55b0172022-02-20 17:52:35 -080058 else if (mid < ncut ) first = mid + 1;
Brian Silverman72890c22015-09-19 14:37:37 -040059 } while (mid != ncut );
Austin Schuhc55b0172022-02-20 17:52:35 -080060
61 return 0; /* mid is equal to ncut */
Brian Silverman72890c22015-09-19 14:37:37 -040062}
63
64}// end namespace internal
65
66/** \ingroup IterativeLinearSolvers_Module
67 * \class IncompleteLUT
68 * \brief Incomplete LU factorization with dual-threshold strategy
69 *
Austin Schuh189376f2018-12-20 22:11:15 +110070 * \implsparsesolverconcept
71 *
Brian Silverman72890c22015-09-19 14:37:37 -040072 * During the numerical factorization, two dropping rules are used :
73 * 1) any element whose magnitude is less than some tolerance is dropped.
Austin Schuhc55b0172022-02-20 17:52:35 -080074 * This tolerance is obtained by multiplying the input tolerance @p droptol
Brian Silverman72890c22015-09-19 14:37:37 -040075 * by the average magnitude of all the original elements in the current row.
Austin Schuhc55b0172022-02-20 17:52:35 -080076 * 2) After the elimination of the row, only the @p fill largest elements in
77 * the L part and the @p fill largest elements in the U part are kept
78 * (in addition to the diagonal element ). Note that @p fill is computed from
79 * the input parameter @p fillfactor which is used the ratio to control the fill_in
Brian Silverman72890c22015-09-19 14:37:37 -040080 * relatively to the initial number of nonzero elements.
Austin Schuhc55b0172022-02-20 17:52:35 -080081 *
Brian Silverman72890c22015-09-19 14:37:37 -040082 * The two extreme cases are when @p droptol=0 (to keep all the @p fill*2 largest elements)
Austin Schuhc55b0172022-02-20 17:52:35 -080083 * and when @p fill=n/2 with @p droptol being different to zero.
84 *
85 * References : Yousef Saad, ILUT: A dual threshold incomplete LU factorization,
Brian Silverman72890c22015-09-19 14:37:37 -040086 * Numerical Linear Algebra with Applications, 1(4), pp 387-402, 1994.
Austin Schuhc55b0172022-02-20 17:52:35 -080087 *
Brian Silverman72890c22015-09-19 14:37:37 -040088 * NOTE : The following implementation is derived from the ILUT implementation
Austin Schuhc55b0172022-02-20 17:52:35 -080089 * in the SPARSKIT package, Copyright (C) 2005, the Regents of the University of Minnesota
90 * released under the terms of the GNU LGPL:
Brian Silverman72890c22015-09-19 14:37:37 -040091 * http://www-users.cs.umn.edu/~saad/software/SPARSKIT/README
92 * However, Yousef Saad gave us permission to relicense his ILUT code to MPL2.
93 * See the Eigen mailing list archive, thread: ILUT, date: July 8, 2012:
94 * http://listengine.tuxfamily.org/lists.tuxfamily.org/eigen/2012/07/msg00064.html
95 * alternatively, on GMANE:
96 * http://comments.gmane.org/gmane.comp.lib.eigen/3302
97 */
Austin Schuh189376f2018-12-20 22:11:15 +110098template <typename _Scalar, typename _StorageIndex = int>
99class IncompleteLUT : public SparseSolverBase<IncompleteLUT<_Scalar, _StorageIndex> >
Brian Silverman72890c22015-09-19 14:37:37 -0400100{
Austin Schuh189376f2018-12-20 22:11:15 +1100101 protected:
102 typedef SparseSolverBase<IncompleteLUT> Base;
103 using Base::m_isInitialized;
104 public:
Brian Silverman72890c22015-09-19 14:37:37 -0400105 typedef _Scalar Scalar;
Austin Schuh189376f2018-12-20 22:11:15 +1100106 typedef _StorageIndex StorageIndex;
Brian Silverman72890c22015-09-19 14:37:37 -0400107 typedef typename NumTraits<Scalar>::Real RealScalar;
108 typedef Matrix<Scalar,Dynamic,1> Vector;
Austin Schuh189376f2018-12-20 22:11:15 +1100109 typedef Matrix<StorageIndex,Dynamic,1> VectorI;
110 typedef SparseMatrix<Scalar,RowMajor,StorageIndex> FactorType;
111
112 enum {
113 ColsAtCompileTime = Dynamic,
114 MaxColsAtCompileTime = Dynamic
115 };
Brian Silverman72890c22015-09-19 14:37:37 -0400116
117 public:
Austin Schuhc55b0172022-02-20 17:52:35 -0800118
Brian Silverman72890c22015-09-19 14:37:37 -0400119 IncompleteLUT()
120 : m_droptol(NumTraits<Scalar>::dummy_precision()), m_fillfactor(10),
Austin Schuh189376f2018-12-20 22:11:15 +1100121 m_analysisIsOk(false), m_factorizationIsOk(false)
Brian Silverman72890c22015-09-19 14:37:37 -0400122 {}
Austin Schuhc55b0172022-02-20 17:52:35 -0800123
Brian Silverman72890c22015-09-19 14:37:37 -0400124 template<typename MatrixType>
Austin Schuh189376f2018-12-20 22:11:15 +1100125 explicit IncompleteLUT(const MatrixType& mat, const RealScalar& droptol=NumTraits<Scalar>::dummy_precision(), int fillfactor = 10)
Brian Silverman72890c22015-09-19 14:37:37 -0400126 : m_droptol(droptol),m_fillfactor(fillfactor),
Austin Schuh189376f2018-12-20 22:11:15 +1100127 m_analysisIsOk(false),m_factorizationIsOk(false)
Brian Silverman72890c22015-09-19 14:37:37 -0400128 {
129 eigen_assert(fillfactor != 0);
Austin Schuhc55b0172022-02-20 17:52:35 -0800130 compute(mat);
Brian Silverman72890c22015-09-19 14:37:37 -0400131 }
Austin Schuhc55b0172022-02-20 17:52:35 -0800132
133 EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT { return m_lu.rows(); }
134
135 EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT { return m_lu.cols(); }
Brian Silverman72890c22015-09-19 14:37:37 -0400136
137 /** \brief Reports whether previous computation was successful.
138 *
Austin Schuhc55b0172022-02-20 17:52:35 -0800139 * \returns \c Success if computation was successful,
Brian Silverman72890c22015-09-19 14:37:37 -0400140 * \c NumericalIssue if the matrix.appears to be negative.
141 */
142 ComputationInfo info() const
143 {
144 eigen_assert(m_isInitialized && "IncompleteLUT is not initialized.");
145 return m_info;
146 }
Austin Schuhc55b0172022-02-20 17:52:35 -0800147
Brian Silverman72890c22015-09-19 14:37:37 -0400148 template<typename MatrixType>
149 void analyzePattern(const MatrixType& amat);
Austin Schuhc55b0172022-02-20 17:52:35 -0800150
Brian Silverman72890c22015-09-19 14:37:37 -0400151 template<typename MatrixType>
152 void factorize(const MatrixType& amat);
Austin Schuhc55b0172022-02-20 17:52:35 -0800153
Brian Silverman72890c22015-09-19 14:37:37 -0400154 /**
155 * Compute an incomplete LU factorization with dual threshold on the matrix mat
156 * No pivoting is done in this version
Austin Schuhc55b0172022-02-20 17:52:35 -0800157 *
Brian Silverman72890c22015-09-19 14:37:37 -0400158 **/
159 template<typename MatrixType>
Austin Schuh189376f2018-12-20 22:11:15 +1100160 IncompleteLUT& compute(const MatrixType& amat)
Brian Silverman72890c22015-09-19 14:37:37 -0400161 {
Austin Schuhc55b0172022-02-20 17:52:35 -0800162 analyzePattern(amat);
Brian Silverman72890c22015-09-19 14:37:37 -0400163 factorize(amat);
164 return *this;
165 }
166
Austin Schuhc55b0172022-02-20 17:52:35 -0800167 void setDroptol(const RealScalar& droptol);
168 void setFillfactor(int fillfactor);
169
Brian Silverman72890c22015-09-19 14:37:37 -0400170 template<typename Rhs, typename Dest>
Austin Schuh189376f2018-12-20 22:11:15 +1100171 void _solve_impl(const Rhs& b, Dest& x) const
Brian Silverman72890c22015-09-19 14:37:37 -0400172 {
Austin Schuh189376f2018-12-20 22:11:15 +1100173 x = m_Pinv * b;
Brian Silverman72890c22015-09-19 14:37:37 -0400174 x = m_lu.template triangularView<UnitLower>().solve(x);
175 x = m_lu.template triangularView<Upper>().solve(x);
Austin Schuhc55b0172022-02-20 17:52:35 -0800176 x = m_P * x;
Brian Silverman72890c22015-09-19 14:37:37 -0400177 }
178
Brian Silverman72890c22015-09-19 14:37:37 -0400179protected:
180
181 /** keeps off-diagonal entries; drops diagonal entries */
182 struct keep_diag {
183 inline bool operator() (const Index& row, const Index& col, const Scalar&) const
184 {
185 return row!=col;
186 }
187 };
188
189protected:
190
191 FactorType m_lu;
192 RealScalar m_droptol;
193 int m_fillfactor;
194 bool m_analysisIsOk;
195 bool m_factorizationIsOk;
Brian Silverman72890c22015-09-19 14:37:37 -0400196 ComputationInfo m_info;
Austin Schuh189376f2018-12-20 22:11:15 +1100197 PermutationMatrix<Dynamic,Dynamic,StorageIndex> m_P; // Fill-reducing permutation
198 PermutationMatrix<Dynamic,Dynamic,StorageIndex> m_Pinv; // Inverse permutation
Brian Silverman72890c22015-09-19 14:37:37 -0400199};
200
201/**
202 * Set control parameter droptol
Austin Schuhc55b0172022-02-20 17:52:35 -0800203 * \param droptol Drop any element whose magnitude is less than this tolerance
204 **/
Austin Schuh189376f2018-12-20 22:11:15 +1100205template<typename Scalar, typename StorageIndex>
206void IncompleteLUT<Scalar,StorageIndex>::setDroptol(const RealScalar& droptol)
Brian Silverman72890c22015-09-19 14:37:37 -0400207{
Austin Schuhc55b0172022-02-20 17:52:35 -0800208 this->m_droptol = droptol;
Brian Silverman72890c22015-09-19 14:37:37 -0400209}
210
211/**
212 * Set control parameter fillfactor
Austin Schuhc55b0172022-02-20 17:52:35 -0800213 * \param fillfactor This is used to compute the number @p fill_in of largest elements to keep on each row.
214 **/
Austin Schuh189376f2018-12-20 22:11:15 +1100215template<typename Scalar, typename StorageIndex>
216void IncompleteLUT<Scalar,StorageIndex>::setFillfactor(int fillfactor)
Brian Silverman72890c22015-09-19 14:37:37 -0400217{
Austin Schuhc55b0172022-02-20 17:52:35 -0800218 this->m_fillfactor = fillfactor;
Brian Silverman72890c22015-09-19 14:37:37 -0400219}
220
Austin Schuh189376f2018-12-20 22:11:15 +1100221template <typename Scalar, typename StorageIndex>
Brian Silverman72890c22015-09-19 14:37:37 -0400222template<typename _MatrixType>
Austin Schuh189376f2018-12-20 22:11:15 +1100223void IncompleteLUT<Scalar,StorageIndex>::analyzePattern(const _MatrixType& amat)
Brian Silverman72890c22015-09-19 14:37:37 -0400224{
225 // Compute the Fill-reducing permutation
Austin Schuh189376f2018-12-20 22:11:15 +1100226 // Since ILUT does not perform any numerical pivoting,
227 // it is highly preferable to keep the diagonal through symmetric permutations.
Austin Schuh189376f2018-12-20 22:11:15 +1100228 // To this end, let's symmetrize the pattern and perform AMD on it.
229 SparseMatrix<Scalar,ColMajor, StorageIndex> mat1 = amat;
230 SparseMatrix<Scalar,ColMajor, StorageIndex> mat2 = amat.transpose();
Brian Silverman72890c22015-09-19 14:37:37 -0400231 // FIXME for a matrix with nearly symmetric pattern, mat2+mat1 is the appropriate choice.
Austin Schuhc55b0172022-02-20 17:52:35 -0800232 // on the other hand for a really non-symmetric pattern, mat2*mat1 should be preferred...
Austin Schuh189376f2018-12-20 22:11:15 +1100233 SparseMatrix<Scalar,ColMajor, StorageIndex> AtA = mat2 + mat1;
234 AMDOrdering<StorageIndex> ordering;
235 ordering(AtA,m_P);
236 m_Pinv = m_P.inverse(); // cache the inverse permutation
Brian Silverman72890c22015-09-19 14:37:37 -0400237 m_analysisIsOk = true;
238 m_factorizationIsOk = false;
Austin Schuh189376f2018-12-20 22:11:15 +1100239 m_isInitialized = true;
Brian Silverman72890c22015-09-19 14:37:37 -0400240}
241
Austin Schuh189376f2018-12-20 22:11:15 +1100242template <typename Scalar, typename StorageIndex>
Brian Silverman72890c22015-09-19 14:37:37 -0400243template<typename _MatrixType>
Austin Schuh189376f2018-12-20 22:11:15 +1100244void IncompleteLUT<Scalar,StorageIndex>::factorize(const _MatrixType& amat)
Brian Silverman72890c22015-09-19 14:37:37 -0400245{
246 using std::sqrt;
247 using std::swap;
248 using std::abs;
Austin Schuh189376f2018-12-20 22:11:15 +1100249 using internal::convert_index;
Brian Silverman72890c22015-09-19 14:37:37 -0400250
251 eigen_assert((amat.rows() == amat.cols()) && "The factorization should be done on a square matrix");
252 Index n = amat.cols(); // Size of the matrix
253 m_lu.resize(n,n);
254 // Declare Working vectors and variables
255 Vector u(n) ; // real values of the row -- maximum size is n --
Austin Schuh189376f2018-12-20 22:11:15 +1100256 VectorI ju(n); // column position of the values in u -- maximum size is n
257 VectorI jr(n); // Indicate the position of the nonzero elements in the vector u -- A zero location is indicated by -1
Brian Silverman72890c22015-09-19 14:37:37 -0400258
259 // Apply the fill-reducing permutation
260 eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
Austin Schuh189376f2018-12-20 22:11:15 +1100261 SparseMatrix<Scalar,RowMajor, StorageIndex> mat;
Brian Silverman72890c22015-09-19 14:37:37 -0400262 mat = amat.twistedBy(m_Pinv);
263
264 // Initialization
265 jr.fill(-1);
266 ju.fill(0);
267 u.fill(0);
268
269 // number of largest elements to keep in each row:
Austin Schuh189376f2018-12-20 22:11:15 +1100270 Index fill_in = (amat.nonZeros()*m_fillfactor)/n + 1;
Brian Silverman72890c22015-09-19 14:37:37 -0400271 if (fill_in > n) fill_in = n;
272
273 // number of largest nonzero elements to keep in the L and the U part of the current row:
274 Index nnzL = fill_in/2;
275 Index nnzU = nnzL;
276 m_lu.reserve(n * (nnzL + nnzU + 1));
277
278 // global loop over the rows of the sparse matrix
279 for (Index ii = 0; ii < n; ii++)
280 {
281 // 1 - copy the lower and the upper part of the row i of mat in the working vector u
282
283 Index sizeu = 1; // number of nonzero elements in the upper part of the current row
284 Index sizel = 0; // number of nonzero elements in the lower part of the current row
Austin Schuh189376f2018-12-20 22:11:15 +1100285 ju(ii) = convert_index<StorageIndex>(ii);
Brian Silverman72890c22015-09-19 14:37:37 -0400286 u(ii) = 0;
Austin Schuh189376f2018-12-20 22:11:15 +1100287 jr(ii) = convert_index<StorageIndex>(ii);
Brian Silverman72890c22015-09-19 14:37:37 -0400288 RealScalar rownorm = 0;
289
290 typename FactorType::InnerIterator j_it(mat, ii); // Iterate through the current row ii
291 for (; j_it; ++j_it)
292 {
293 Index k = j_it.index();
294 if (k < ii)
295 {
296 // copy the lower part
Austin Schuh189376f2018-12-20 22:11:15 +1100297 ju(sizel) = convert_index<StorageIndex>(k);
Brian Silverman72890c22015-09-19 14:37:37 -0400298 u(sizel) = j_it.value();
Austin Schuh189376f2018-12-20 22:11:15 +1100299 jr(k) = convert_index<StorageIndex>(sizel);
Brian Silverman72890c22015-09-19 14:37:37 -0400300 ++sizel;
301 }
302 else if (k == ii)
303 {
304 u(ii) = j_it.value();
305 }
306 else
307 {
308 // copy the upper part
309 Index jpos = ii + sizeu;
Austin Schuh189376f2018-12-20 22:11:15 +1100310 ju(jpos) = convert_index<StorageIndex>(k);
Brian Silverman72890c22015-09-19 14:37:37 -0400311 u(jpos) = j_it.value();
Austin Schuh189376f2018-12-20 22:11:15 +1100312 jr(k) = convert_index<StorageIndex>(jpos);
Brian Silverman72890c22015-09-19 14:37:37 -0400313 ++sizeu;
314 }
315 rownorm += numext::abs2(j_it.value());
316 }
317
318 // 2 - detect possible zero row
319 if(rownorm==0)
320 {
321 m_info = NumericalIssue;
322 return;
323 }
324 // Take the 2-norm of the current row as a relative tolerance
325 rownorm = sqrt(rownorm);
326
327 // 3 - eliminate the previous nonzero rows
328 Index jj = 0;
329 Index len = 0;
330 while (jj < sizel)
331 {
332 // In order to eliminate in the correct order,
333 // we must select first the smallest column index among ju(jj:sizel)
334 Index k;
335 Index minrow = ju.segment(jj,sizel-jj).minCoeff(&k); // k is relative to the segment
336 k += jj;
337 if (minrow != ju(jj))
338 {
339 // swap the two locations
340 Index j = ju(jj);
341 swap(ju(jj), ju(k));
Austin Schuh189376f2018-12-20 22:11:15 +1100342 jr(minrow) = convert_index<StorageIndex>(jj);
343 jr(j) = convert_index<StorageIndex>(k);
Brian Silverman72890c22015-09-19 14:37:37 -0400344 swap(u(jj), u(k));
345 }
346 // Reset this location
347 jr(minrow) = -1;
348
349 // Start elimination
350 typename FactorType::InnerIterator ki_it(m_lu, minrow);
351 while (ki_it && ki_it.index() < minrow) ++ki_it;
352 eigen_internal_assert(ki_it && ki_it.col()==minrow);
353 Scalar fact = u(jj) / ki_it.value();
354
355 // drop too small elements
356 if(abs(fact) <= m_droptol)
357 {
358 jj++;
359 continue;
360 }
361
362 // linear combination of the current row ii and the row minrow
363 ++ki_it;
364 for (; ki_it; ++ki_it)
365 {
366 Scalar prod = fact * ki_it.value();
Austin Schuh189376f2018-12-20 22:11:15 +1100367 Index j = ki_it.index();
368 Index jpos = jr(j);
Brian Silverman72890c22015-09-19 14:37:37 -0400369 if (jpos == -1) // fill-in element
370 {
371 Index newpos;
372 if (j >= ii) // dealing with the upper part
373 {
374 newpos = ii + sizeu;
375 sizeu++;
376 eigen_internal_assert(sizeu<=n);
377 }
378 else // dealing with the lower part
379 {
380 newpos = sizel;
381 sizel++;
382 eigen_internal_assert(sizel<=ii);
383 }
Austin Schuh189376f2018-12-20 22:11:15 +1100384 ju(newpos) = convert_index<StorageIndex>(j);
Brian Silverman72890c22015-09-19 14:37:37 -0400385 u(newpos) = -prod;
Austin Schuh189376f2018-12-20 22:11:15 +1100386 jr(j) = convert_index<StorageIndex>(newpos);
Brian Silverman72890c22015-09-19 14:37:37 -0400387 }
388 else
389 u(jpos) -= prod;
390 }
391 // store the pivot element
Austin Schuh189376f2018-12-20 22:11:15 +1100392 u(len) = fact;
393 ju(len) = convert_index<StorageIndex>(minrow);
Brian Silverman72890c22015-09-19 14:37:37 -0400394 ++len;
395
396 jj++;
397 } // end of the elimination on the row ii
398
399 // reset the upper part of the pointer jr to zero
400 for(Index k = 0; k <sizeu; k++) jr(ju(ii+k)) = -1;
401
402 // 4 - partially sort and insert the elements in the m_lu matrix
403
404 // sort the L-part of the row
405 sizel = len;
406 len = (std::min)(sizel, nnzL);
407 typename Vector::SegmentReturnType ul(u.segment(0, sizel));
Austin Schuh189376f2018-12-20 22:11:15 +1100408 typename VectorI::SegmentReturnType jul(ju.segment(0, sizel));
Brian Silverman72890c22015-09-19 14:37:37 -0400409 internal::QuickSplit(ul, jul, len);
410
411 // store the largest m_fill elements of the L part
412 m_lu.startVec(ii);
413 for(Index k = 0; k < len; k++)
414 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))
419 u(ii) = sqrt(m_droptol) * rownorm;
420 m_lu.insertBackByOuterInnerUnordered(ii, ii) = u(ii);
421
422 // sort the U-part of the row
423 // apply the dropping rule first
424 len = 0;
425 for(Index k = 1; k < sizeu; k++)
426 {
427 if(abs(u(ii+k)) > m_droptol * rownorm )
428 {
429 ++len;
430 u(ii + len) = u(ii + k);
431 ju(ii + len) = ju(ii + k);
432 }
433 }
434 sizeu = len + 1; // +1 to take into account the diagonal element
435 len = (std::min)(sizeu, nnzU);
436 typename Vector::SegmentReturnType uu(u.segment(ii+1, sizeu-1));
Austin Schuh189376f2018-12-20 22:11:15 +1100437 typename VectorI::SegmentReturnType juu(ju.segment(ii+1, sizeu-1));
Brian Silverman72890c22015-09-19 14:37:37 -0400438 internal::QuickSplit(uu, juu, len);
439
440 // store the largest elements of the U part
441 for(Index k = ii + 1; k < ii + len; k++)
442 m_lu.insertBackByOuterInnerUnordered(ii,ju(k)) = u(k);
443 }
Brian Silverman72890c22015-09-19 14:37:37 -0400444 m_lu.finalize();
445 m_lu.makeCompressed();
446
447 m_factorizationIsOk = true;
Brian Silverman72890c22015-09-19 14:37:37 -0400448 m_info = Success;
449}
450
Brian Silverman72890c22015-09-19 14:37:37 -0400451} // end namespace Eigen
452
453#endif // EIGEN_INCOMPLETE_LUT_H