Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::IncompleteLUT< _Scalar, _StorageIndex > Class Template Reference

Incomplete LU factorization with dual-threshold strategy. More...

#include <src/eigen/Eigen/src/IterativeLinearSolvers/IncompleteLUT.h>

+ Inheritance diagram for Eigen::IncompleteLUT< _Scalar, _StorageIndex >:
+ Collaboration diagram for Eigen::IncompleteLUT< _Scalar, _StorageIndex >:

Classes

struct  keep_diag
 

Public Types

enum  { ColsAtCompileTime = Dynamic , MaxColsAtCompileTime = Dynamic }
 
typedef _Scalar Scalar
 
typedef _StorageIndex StorageIndex
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef Matrix< Scalar, Dynamic, 1 > Vector
 
typedef Matrix< StorageIndex, Dynamic, 1 > VectorI
 
typedef SparseMatrix< Scalar, RowMajor, StorageIndexFactorType
 

Public Member Functions

 IncompleteLUT ()
 
template<typename MatrixType >
 IncompleteLUT (const MatrixType &mat, const RealScalar &droptol=NumTraits< Scalar >::dummy_precision(), int fillfactor=10)
 
Index rows () const
 
Index cols () const
 
ComputationInfo info () const
 Reports whether previous computation was successful.
 
template<typename MatrixType >
void analyzePattern (const MatrixType &amat)
 
template<typename MatrixType >
void factorize (const MatrixType &amat)
 
template<typename MatrixType >
IncompleteLUTcompute (const MatrixType &amat)
 
void setDroptol (const RealScalar &droptol)
 
void setFillfactor (int fillfactor)
 
template<typename Rhs , typename Dest >
void _solve_impl (const Rhs &b, Dest &x) const
 
template<typename _MatrixType >
void analyzePattern (const _MatrixType &amat)
 
template<typename _MatrixType >
void factorize (const _MatrixType &amat)
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Protected Types

typedef SparseSolverBase< IncompleteLUTBase
 

Protected Attributes

FactorType m_lu
 
RealScalar m_droptol
 
int m_fillfactor
 
bool m_analysisIsOk
 
bool m_factorizationIsOk
 
ComputationInfo m_info
 
PermutationMatrix< Dynamic, Dynamic, StorageIndexm_P
 
PermutationMatrix< Dynamic, Dynamic, StorageIndexm_Pinv
 
bool m_isInitialized
 

Detailed Description

template<typename _Scalar, typename _StorageIndex = int>
class Eigen::IncompleteLUT< _Scalar, _StorageIndex >

Incomplete LU factorization with dual-threshold strategy.

\implsparsesolverconcept

During the numerical factorization, two dropping rules are used : 1) any element whose magnitude is less than some tolerance is dropped. This tolerance is obtained by multiplying the input tolerance droptol by the average magnitude of all the original elements in the current row. 2) After the elimination of the row, only the fill largest elements in the L part and the fill largest elements in the U part are kept (in addition to the diagonal element ). Note that fill is computed from the input parameter fillfactor which is used the ratio to control the fill_in relatively to the initial number of nonzero elements.

The two extreme cases are when droptol=0 (to keep all the fill*2 largest elements) and when fill=n/2 with droptol being different to zero.

References : Yousef Saad, ILUT: A dual threshold incomplete LU factorization, Numerical Linear Algebra with Applications, 1(4), pp 387-402, 1994.

NOTE : The following implementation is derived from the ILUT implementation in the SPARSKIT package, Copyright (C) 2005, the Regents of the University of Minnesota released under the terms of the GNU LGPL: http://www-users.cs.umn.edu/~saad/software/SPARSKIT/README However, Yousef Saad gave us permission to relicense his ILUT code to MPL2. See the Eigen mailing list archive, thread: ILUT, date: July 8, 2012: http://listengine.tuxfamily.org/lists.tuxfamily.org/eigen/2012/07/msg00064.html alternatively, on GMANE: http://comments.gmane.org/gmane.comp.lib.eigen/3302

Member Typedef Documentation

◆ Base

template<typename _Scalar , typename _StorageIndex = int>
typedef SparseSolverBase<IncompleteLUT> Eigen::IncompleteLUT< _Scalar, _StorageIndex >::Base
protected

◆ FactorType

template<typename _Scalar , typename _StorageIndex = int>
typedef SparseMatrix<Scalar,RowMajor,StorageIndex> Eigen::IncompleteLUT< _Scalar, _StorageIndex >::FactorType

◆ RealScalar

template<typename _Scalar , typename _StorageIndex = int>
typedef NumTraits<Scalar>::Real Eigen::IncompleteLUT< _Scalar, _StorageIndex >::RealScalar

◆ Scalar

template<typename _Scalar , typename _StorageIndex = int>
typedef _Scalar Eigen::IncompleteLUT< _Scalar, _StorageIndex >::Scalar

◆ StorageIndex

template<typename _Scalar , typename _StorageIndex = int>
typedef _StorageIndex Eigen::IncompleteLUT< _Scalar, _StorageIndex >::StorageIndex

◆ Vector

template<typename _Scalar , typename _StorageIndex = int>
typedef Matrix<Scalar,Dynamic,1> Eigen::IncompleteLUT< _Scalar, _StorageIndex >::Vector

◆ VectorI

template<typename _Scalar , typename _StorageIndex = int>
typedef Matrix<StorageIndex,Dynamic,1> Eigen::IncompleteLUT< _Scalar, _StorageIndex >::VectorI

Member Enumeration Documentation

◆ anonymous enum

template<typename _Scalar , typename _StorageIndex = int>
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 
112 {
115 };
@ ColsAtCompileTime
Definition IncompleteLUT.h:113
@ MaxColsAtCompileTime
Definition IncompleteLUT.h:114
const int Dynamic
Definition Constants.h:21

Constructor & Destructor Documentation

◆ IncompleteLUT() [1/2]

template<typename _Scalar , typename _StorageIndex = int>
Eigen::IncompleteLUT< _Scalar, _StorageIndex >::IncompleteLUT ( )
inline
120 : m_droptol(NumTraits<Scalar>::dummy_precision()), m_fillfactor(10),
122 {}
bool m_analysisIsOk
Definition IncompleteLUT.h:194
bool m_factorizationIsOk
Definition IncompleteLUT.h:195
int m_fillfactor
Definition IncompleteLUT.h:193
RealScalar m_droptol
Definition IncompleteLUT.h:192

◆ IncompleteLUT() [2/2]

template<typename _Scalar , typename _StorageIndex = int>
template<typename MatrixType >
Eigen::IncompleteLUT< _Scalar, _StorageIndex >::IncompleteLUT ( const MatrixType &  mat,
const RealScalar droptol = NumTraits<Scalar>::dummy_precision(),
int  fillfactor = 10 
)
inlineexplicit
126 : m_droptol(droptol),m_fillfactor(fillfactor),
128 {
129 eigen_assert(fillfactor != 0);
130 compute(mat);
131 }
#define eigen_assert(x)
Definition Macros.h:579
IncompleteLUT & compute(const MatrixType &amat)
Definition IncompleteLUT.h:160

References Eigen::IncompleteLUT< _Scalar, _StorageIndex >::compute(), and eigen_assert.

+ Here is the call graph for this function:

Member Function Documentation

◆ _solve_impl() [1/2]

template<typename _Scalar , typename _StorageIndex = int>
template<typename Rhs , typename Dest >
void Eigen::IncompleteLUT< _Scalar, _StorageIndex >::_solve_impl ( const Rhs &  b,
Dest &  x 
) const
inline
172 {
173 x = m_Pinv * b;
174 x = m_lu.template triangularView<UnitLower>().solve(x);
175 x = m_lu.template triangularView<Upper>().solve(x);
176 x = m_P * x;
177 }
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_Pinv
Definition IncompleteLUT.h:198
FactorType m_lu
Definition IncompleteLUT.h:191
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_P
Definition IncompleteLUT.h:197
TCoord< P > x(const P &p)
Definition geometry_traits.hpp:297

References Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_lu, Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_P, and Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_Pinv.

◆ _solve_impl() [2/2]

template<typename Derived >
template<typename Rhs , typename Dest >
void Eigen::SparseSolverBase< Derived >::_solve_impl ( const SparseMatrixBase< Rhs > &  b,
SparseMatrixBase< Dest > &  dest 
) const
inlineinherited
112 {
113 internal::solve_sparse_through_dense_panels(derived(), b.derived(), dest.derived());
114 }
Derived & derived()
Definition SparseSolverBase.h:79
enable_if< Rhs::ColsAtCompileTime!=1 &&Dest::ColsAtCompileTime!=1 >::type solve_sparse_through_dense_panels(const Decomposition &dec, const Rhs &rhs, Dest &dest)
Definition SparseSolverBase.h:23

References Eigen::SparseSolverBase< Derived >::derived(), Eigen::SparseMatrixBase< Derived >::derived(), and Eigen::internal::solve_sparse_through_dense_panels().

+ Here is the call graph for this function:

◆ analyzePattern() [1/2]

template<typename _Scalar , typename _StorageIndex = int>
template<typename _MatrixType >
void Eigen::IncompleteLUT< _Scalar, _StorageIndex >::analyzePattern ( const _MatrixType &  amat)
224{
225 // Compute the Fill-reducing permutation
226 // Since ILUT does not perform any numerical pivoting,
227 // it is highly preferable to keep the diagonal through symmetric permutations.
228#ifndef EIGEN_MPL2_ONLY
229 // To this end, let's symmetrize the pattern and perform AMD on it.
230 SparseMatrix<Scalar,ColMajor, StorageIndex> mat1 = amat;
231 SparseMatrix<Scalar,ColMajor, StorageIndex> mat2 = amat.transpose();
232 // FIXME for a matrix with nearly symmetric pattern, mat2+mat1 is the appropriate choice.
233 // on the other hand for a really non-symmetric pattern, mat2*mat1 should be prefered...
234 SparseMatrix<Scalar,ColMajor, StorageIndex> AtA = mat2 + mat1;
235 AMDOrdering<StorageIndex> ordering;
236 ordering(AtA,m_P);
237 m_Pinv = m_P.inverse(); // cache the inverse permutation
238#else
239 // If AMD is not available, (MPL2-only), then let's use the slower COLAMD routine.
240 SparseMatrix<Scalar,ColMajor, StorageIndex> mat1 = amat;
241 COLAMDOrdering<StorageIndex> ordering;
242 ordering(mat1,m_Pinv);
243 m_P = m_Pinv.inverse();
244#endif
245
246 m_analysisIsOk = true;
247 m_factorizationIsOk = false;
248 m_isInitialized = true;
249}
bool m_isInitialized
Definition SparseSolverBase.h:119
InverseReturnType inverse() const
Definition PermutationMatrix.h:196
static enum @17 ordering

References ordering, and Eigen::SparseMatrixBase< Derived >::transpose().

+ Here is the call graph for this function:

◆ analyzePattern() [2/2]

template<typename _Scalar , typename _StorageIndex = int>
template<typename MatrixType >
void Eigen::IncompleteLUT< _Scalar, _StorageIndex >::analyzePattern ( const MatrixType &  amat)

Referenced by Eigen::IncompleteLUT< _Scalar, _StorageIndex >::compute().

+ Here is the caller graph for this function:

◆ cols()

template<typename _Scalar , typename _StorageIndex = int>
Index Eigen::IncompleteLUT< _Scalar, _StorageIndex >::cols ( ) const
inline
135{ return m_lu.cols(); }
Index cols() const
Definition SparseMatrix.h:138

References Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::cols(), and Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_lu.

+ Here is the call graph for this function:

◆ compute()

template<typename _Scalar , typename _StorageIndex = int>
template<typename MatrixType >
IncompleteLUT & Eigen::IncompleteLUT< _Scalar, _StorageIndex >::compute ( const MatrixType &  amat)
inline

Compute an incomplete LU factorization with dual threshold on the matrix mat No pivoting is done in this version

161 {
162 analyzePattern(amat);
163 factorize(amat);
164 return *this;
165 }
void factorize(const MatrixType &amat)
void analyzePattern(const MatrixType &amat)

References Eigen::IncompleteLUT< _Scalar, _StorageIndex >::analyzePattern(), and Eigen::IncompleteLUT< _Scalar, _StorageIndex >::factorize().

Referenced by Eigen::IncompleteLUT< _Scalar, _StorageIndex >::IncompleteLUT().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ derived() [1/2]

template<typename Derived >
Derived & Eigen::SparseSolverBase< Derived >::derived ( )
inlineinherited
79{ return *static_cast<Derived*>(this); }

Referenced by Eigen::SparseSolverBase< Derived >::_solve_impl(), and Eigen::SparseSolverBase< Derived >::solve().

+ Here is the caller graph for this function:

◆ derived() [2/2]

template<typename Derived >
const Derived & Eigen::SparseSolverBase< Derived >::derived ( ) const
inlineinherited
80{ return *static_cast<const Derived*>(this); }

◆ factorize() [1/2]

template<typename _Scalar , typename _StorageIndex = int>
template<typename _MatrixType >
void Eigen::IncompleteLUT< _Scalar, _StorageIndex >::factorize ( const _MatrixType &  amat)
254{
255 using std::sqrt;
256 using std::swap;
257 using std::abs;
259
260 eigen_assert((amat.rows() == amat.cols()) && "The factorization should be done on a square matrix");
261 Index n = amat.cols(); // Size of the matrix
262 m_lu.resize(n,n);
263 // Declare Working vectors and variables
264 Vector u(n) ; // real values of the row -- maximum size is n --
265 VectorI ju(n); // column position of the values in u -- maximum size is n
266 VectorI jr(n); // Indicate the position of the nonzero elements in the vector u -- A zero location is indicated by -1
267
268 // Apply the fill-reducing permutation
269 eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
270 SparseMatrix<Scalar,RowMajor, StorageIndex> mat;
271 mat = amat.twistedBy(m_Pinv);
272
273 // Initialization
274 jr.fill(-1);
275 ju.fill(0);
276 u.fill(0);
277
278 // number of largest elements to keep in each row:
279 Index fill_in = (amat.nonZeros()*m_fillfactor)/n + 1;
280 if (fill_in > n) fill_in = n;
281
282 // number of largest nonzero elements to keep in the L and the U part of the current row:
283 Index nnzL = fill_in/2;
284 Index nnzU = nnzL;
285 m_lu.reserve(n * (nnzL + nnzU + 1));
286
287 // global loop over the rows of the sparse matrix
288 for (Index ii = 0; ii < n; ii++)
289 {
290 // 1 - copy the lower and the upper part of the row i of mat in the working vector u
291
292 Index sizeu = 1; // number of nonzero elements in the upper part of the current row
293 Index sizel = 0; // number of nonzero elements in the lower part of the current row
294 ju(ii) = convert_index<StorageIndex>(ii);
295 u(ii) = 0;
296 jr(ii) = convert_index<StorageIndex>(ii);
297 RealScalar rownorm = 0;
298
299 typename FactorType::InnerIterator j_it(mat, ii); // Iterate through the current row ii
300 for (; j_it; ++j_it)
301 {
302 Index k = j_it.index();
303 if (k < ii)
304 {
305 // copy the lower part
306 ju(sizel) = convert_index<StorageIndex>(k);
307 u(sizel) = j_it.value();
308 jr(k) = convert_index<StorageIndex>(sizel);
309 ++sizel;
310 }
311 else if (k == ii)
312 {
313 u(ii) = j_it.value();
314 }
315 else
316 {
317 // copy the upper part
318 Index jpos = ii + sizeu;
319 ju(jpos) = convert_index<StorageIndex>(k);
320 u(jpos) = j_it.value();
321 jr(k) = convert_index<StorageIndex>(jpos);
322 ++sizeu;
323 }
324 rownorm += numext::abs2(j_it.value());
325 }
326
327 // 2 - detect possible zero row
328 if(rownorm==0)
329 {
331 return;
332 }
333 // Take the 2-norm of the current row as a relative tolerance
334 rownorm = sqrt(rownorm);
335
336 // 3 - eliminate the previous nonzero rows
337 Index jj = 0;
338 Index len = 0;
339 while (jj < sizel)
340 {
341 // In order to eliminate in the correct order,
342 // we must select first the smallest column index among ju(jj:sizel)
343 Index k;
344 Index minrow = ju.segment(jj,sizel-jj).minCoeff(&k); // k is relative to the segment
345 k += jj;
346 if (minrow != ju(jj))
347 {
348 // swap the two locations
349 Index j = ju(jj);
350 swap(ju(jj), ju(k));
351 jr(minrow) = convert_index<StorageIndex>(jj);
352 jr(j) = convert_index<StorageIndex>(k);
353 swap(u(jj), u(k));
354 }
355 // Reset this location
356 jr(minrow) = -1;
357
358 // Start elimination
359 typename FactorType::InnerIterator ki_it(m_lu, minrow);
360 while (ki_it && ki_it.index() < minrow) ++ki_it;
361 eigen_internal_assert(ki_it && ki_it.col()==minrow);
362 Scalar fact = u(jj) / ki_it.value();
363
364 // drop too small elements
365 if(abs(fact) <= m_droptol)
366 {
367 jj++;
368 continue;
369 }
370
371 // linear combination of the current row ii and the row minrow
372 ++ki_it;
373 for (; ki_it; ++ki_it)
374 {
375 Scalar prod = fact * ki_it.value();
376 Index j = ki_it.index();
377 Index jpos = jr(j);
378 if (jpos == -1) // fill-in element
379 {
380 Index newpos;
381 if (j >= ii) // dealing with the upper part
382 {
383 newpos = ii + sizeu;
384 sizeu++;
385 eigen_internal_assert(sizeu<=n);
386 }
387 else // dealing with the lower part
388 {
389 newpos = sizel;
390 sizel++;
391 eigen_internal_assert(sizel<=ii);
392 }
393 ju(newpos) = convert_index<StorageIndex>(j);
394 u(newpos) = -prod;
395 jr(j) = convert_index<StorageIndex>(newpos);
396 }
397 else
398 u(jpos) -= prod;
399 }
400 // store the pivot element
401 u(len) = fact;
402 ju(len) = convert_index<StorageIndex>(minrow);
403 ++len;
404
405 jj++;
406 } // end of the elimination on the row ii
407
408 // reset the upper part of the pointer jr to zero
409 for(Index k = 0; k <sizeu; k++) jr(ju(ii+k)) = -1;
410
411 // 4 - partially sort and insert the elements in the m_lu matrix
412
413 // sort the L-part of the row
414 sizel = len;
415 len = (std::min)(sizel, nnzL);
416 typename Vector::SegmentReturnType ul(u.segment(0, sizel));
417 typename VectorI::SegmentReturnType jul(ju.segment(0, sizel));
418 internal::QuickSplit(ul, jul, len);
419
420 // store the largest m_fill elements of the L part
421 m_lu.startVec(ii);
422 for(Index k = 0; k < len; k++)
423 m_lu.insertBackByOuterInnerUnordered(ii,ju(k)) = u(k);
424
425 // store the diagonal element
426 // apply a shifting rule to avoid zero pivots (we are doing an incomplete factorization)
427 if (u(ii) == Scalar(0))
428 u(ii) = sqrt(m_droptol) * rownorm;
430
431 // sort the U-part of the row
432 // apply the dropping rule first
433 len = 0;
434 for(Index k = 1; k < sizeu; k++)
435 {
436 if(abs(u(ii+k)) > m_droptol * rownorm )
437 {
438 ++len;
439 u(ii + len) = u(ii + k);
440 ju(ii + len) = ju(ii + k);
441 }
442 }
443 sizeu = len + 1; // +1 to take into account the diagonal element
444 len = (std::min)(sizeu, nnzU);
445 typename Vector::SegmentReturnType uu(u.segment(ii+1, sizeu-1));
446 typename VectorI::SegmentReturnType juu(ju.segment(ii+1, sizeu-1));
447 internal::QuickSplit(uu, juu, len);
448
449 // store the largest elements of the U part
450 for(Index k = ii + 1; k < ii + len; k++)
451 m_lu.insertBackByOuterInnerUnordered(ii,ju(k)) = u(k);
452 }
453 m_lu.finalize();
455
456 m_factorizationIsOk = true;
457 m_info = Success;
458}
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Definition ArrayCwiseUnaryOps.h:152
#define eigen_internal_assert(x)
Definition Macros.h:585
Matrix< StorageIndex, Dynamic, 1 > VectorI
Definition IncompleteLUT.h:109
_Scalar Scalar
Definition IncompleteLUT.h:105
ComputationInfo m_info
Definition IncompleteLUT.h:196
Matrix< Scalar, Dynamic, 1 > Vector
Definition IncompleteLUT.h:108
NumTraits< Scalar >::Real RealScalar
Definition IncompleteLUT.h:107
Base::InnerIterator InnerIterator
Definition SparseMatrix.h:112
void reserve(Index reserveSize)
Definition SparseMatrix.h:262
void startVec(Index outer)
Definition SparseMatrix.h:412
Scalar & insertBackByOuterInnerUnordered(Index outer, Index inner)
Definition SparseMatrix.h:402
void finalize()
Definition SparseMatrix.h:422
void makeCompressed()
Definition SparseMatrix.h:464
void resize(Index rows, Index cols)
Definition SparseMatrix.h:621
@ NumericalIssue
Definition Constants.h:434
@ Success
Definition Constants.h:432
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC half abs(const half &a)
Definition Half.h:445
EIGEN_DEVICE_FUNC IndexDest convert_index(const IndexSrc &idx)
Definition XprHelper.h:31
Index QuickSplit(VectorV &row, VectorI &ind, Index ncut)
Definition IncompleteLUT.h:29
void swap(scoped_array< T > &a, scoped_array< T > &b)
Definition Memory.h:602
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:33

References Eigen::SparseCompressedBase< Derived >::InnerIterator::col(), Eigen::internal::convert_index(), eigen_assert, eigen_internal_assert, Eigen::SparseCompressedBase< Derived >::InnerIterator::index(), Eigen::NumericalIssue, Eigen::internal::QuickSplit(), sqrt(), Eigen::Success, Eigen::SparseMatrixBase< Derived >::twistedBy(), and Eigen::SparseCompressedBase< Derived >::InnerIterator::value().

+ Here is the call graph for this function:

◆ factorize() [2/2]

template<typename _Scalar , typename _StorageIndex = int>
template<typename MatrixType >
void Eigen::IncompleteLUT< _Scalar, _StorageIndex >::factorize ( const MatrixType &  amat)

Referenced by Eigen::IncompleteLUT< _Scalar, _StorageIndex >::compute().

+ Here is the caller graph for this function:

◆ info()

template<typename _Scalar , typename _StorageIndex = int>
ComputationInfo Eigen::IncompleteLUT< _Scalar, _StorageIndex >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was succesful, NumericalIssue if the matrix.appears to be negative.
143 {
144 eigen_assert(m_isInitialized && "IncompleteLUT is not initialized.");
145 return m_info;
146 }

References eigen_assert, Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_info, and Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_isInitialized.

◆ rows()

template<typename _Scalar , typename _StorageIndex = int>
Index Eigen::IncompleteLUT< _Scalar, _StorageIndex >::rows ( ) const
inline
133{ return m_lu.rows(); }
Index rows() const
Definition SparseMatrix.h:136

References Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_lu, and Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::rows().

+ Here is the call graph for this function:

◆ setDroptol()

template<typename Scalar , typename StorageIndex >
void Eigen::IncompleteLUT< Scalar, StorageIndex >::setDroptol ( const RealScalar droptol)

Set control parameter droptol

Parameters
droptolDrop any element whose magnitude is less than this tolerance
207{
208 this->m_droptol = droptol;
209}

◆ setFillfactor()

template<typename Scalar , typename StorageIndex >
void Eigen::IncompleteLUT< Scalar, StorageIndex >::setFillfactor ( int  fillfactor)

Set control parameter fillfactor

Parameters
fillfactorThis is used to compute the number fill_in of largest elements to keep on each row.
217{
218 this->m_fillfactor = fillfactor;
219}

◆ solve() [1/2]

template<typename Derived >
template<typename Rhs >
const Solve< Derived, Rhs > Eigen::SparseSolverBase< Derived >::solve ( const MatrixBase< Rhs > &  b) const
inlineinherited
Returns
an expression of the solution x of $ A x = b $ using the current decomposition of A.
See also
compute()
89 {
90 eigen_assert(m_isInitialized && "Solver is not initialized.");
91 eigen_assert(derived().rows()==b.rows() && "solve(): invalid number of rows of the right hand side matrix b");
92 return Solve<Derived, Rhs>(derived(), b.derived());
93 }
bool m_isInitialized
Definition SparseSolverBase.h:119
size_t rows(const T &raster)
Definition MarchingSquares.hpp:55

References Eigen::SparseSolverBase< Derived >::derived(), eigen_assert, and Eigen::SparseSolverBase< Derived >::m_isInitialized.

Referenced by igl::embree::bone_heat(), igl::Frame_field_deformer::compute_optimal_positions(), igl::eigs(), igl::copyleft::comiso::FrameInterpolator::interpolateSymmetric(), igl::slim::solve_weighted_arap(), and igl::copyleft::comiso::NRosyField::solveNoRoundings().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ solve() [2/2]

template<typename Derived >
template<typename Rhs >
const Solve< Derived, Rhs > Eigen::SparseSolverBase< Derived >::solve ( const SparseMatrixBase< Rhs > &  b) const
inlineinherited
Returns
an expression of the solution x of $ A x = b $ using the current decomposition of A.
See also
compute()
102 {
103 eigen_assert(m_isInitialized && "Solver is not initialized.");
104 eigen_assert(derived().rows()==b.rows() && "solve(): invalid number of rows of the right hand side matrix b");
105 return Solve<Derived, Rhs>(derived(), b.derived());
106 }

Member Data Documentation

◆ m_analysisIsOk

template<typename _Scalar , typename _StorageIndex = int>
bool Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_analysisIsOk
protected

◆ m_droptol

template<typename _Scalar , typename _StorageIndex = int>
RealScalar Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_droptol
protected

◆ m_factorizationIsOk

template<typename _Scalar , typename _StorageIndex = int>
bool Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_factorizationIsOk
protected

◆ m_fillfactor

template<typename _Scalar , typename _StorageIndex = int>
int Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_fillfactor
protected

◆ m_info

template<typename _Scalar , typename _StorageIndex = int>
ComputationInfo Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_info
protected

◆ m_isInitialized

template<typename _Scalar , typename _StorageIndex = int>
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

◆ m_lu

◆ m_P

template<typename _Scalar , typename _StorageIndex = int>
PermutationMatrix<Dynamic,Dynamic,StorageIndex> Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_P
protected

◆ m_Pinv

template<typename _Scalar , typename _StorageIndex = int>
PermutationMatrix<Dynamic,Dynamic,StorageIndex> Eigen::IncompleteLUT< _Scalar, _StorageIndex >::m_Pinv
protected

The documentation for this class was generated from the following file: