11#ifndef EIGEN_INCOMPLETE_LUT_H
12#define EIGEN_INCOMPLETE_LUT_H
28template <
typename VectorV,
typename VectorI>
31 typedef typename VectorV::RealScalar
RealScalar;
98template <
typename _Scalar,
typename _StorageIndex =
int>
124 template<
typename MatrixType>
148 template<
typename MatrixType>
151 template<
typename MatrixType>
159 template<
typename MatrixType>
170 template<
typename Rhs,
typename Dest>
174 x =
m_lu.template triangularView<UnitLower>().solve(
x);
175 x =
m_lu.template triangularView<Upper>().solve(
x);
205template<
typename Scalar,
typename StorageIndex>
208 this->m_droptol = droptol;
215template<
typename Scalar,
typename StorageIndex>
218 this->m_fillfactor = fillfactor;
221template <
typename Scalar,
typename StorageIndex>
222template<
typename _MatrixType>
236 m_Pinv = m_P.inverse();
237 m_analysisIsOk =
true;
238 m_factorizationIsOk =
false;
239 m_isInitialized =
true;
242template <
typename Scalar,
typename StorageIndex>
243template<
typename _MatrixType>
251 eigen_assert((amat.rows() == amat.cols()) &&
"The factorization should be done on a square matrix");
260 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
270 Index fill_in = (amat.nonZeros()*m_fillfactor)/
n + 1;
271 if (fill_in >
n) fill_in =
n;
274 Index nnzL = fill_in/2;
276 m_lu.reserve(
n * (nnzL + nnzU + 1));
279 for (
Index ii = 0; ii <
n; ii++)
285 ju(ii) = convert_index<StorageIndex>(ii);
287 jr(ii) = convert_index<StorageIndex>(ii);
297 ju(sizel) = convert_index<StorageIndex>(k);
298 u(sizel) = j_it.
value();
299 jr(k) = convert_index<StorageIndex>(sizel);
304 u(ii) = j_it.
value();
309 Index jpos = ii + sizeu;
310 ju(jpos) = convert_index<StorageIndex>(k);
311 u(jpos) = j_it.
value();
312 jr(k) = convert_index<StorageIndex>(jpos);
325 rownorm = sqrt(rownorm);
335 Index minrow = ju.segment(jj,sizel-jj).minCoeff(&k);
337 if (minrow != ju(jj))
342 jr(minrow) = convert_index<StorageIndex>(jj);
343 jr(
j) = convert_index<StorageIndex>(k);
351 while (ki_it && ki_it.
index() < minrow) ++ki_it;
356 if(
abs(fact) <= m_droptol)
364 for (; ki_it; ++ki_it)
384 ju(newpos) = convert_index<StorageIndex>(
j);
386 jr(
j) = convert_index<StorageIndex>(newpos);
393 ju(len) = convert_index<StorageIndex>(minrow);
400 for(
Index k = 0; k <sizeu; k++) jr(ju(ii+k)) = -1;
406 len = (std::min)(sizel, nnzL);
407 typename Vector::SegmentReturnType ul(u.segment(0, sizel));
408 typename VectorI::SegmentReturnType jul(ju.segment(0, sizel));
413 for(
Index k = 0; k < len; k++)
414 m_lu.insertBackByOuterInnerUnordered(ii,ju(k)) = u(k);
419 u(ii) = sqrt(m_droptol) * rownorm;
420 m_lu.insertBackByOuterInnerUnordered(ii, ii) = u(ii);
425 for(
Index k = 1; k < sizeu; k++)
427 if(
abs(u(ii+k)) > m_droptol * rownorm )
430 u(ii + len) = u(ii + k);
431 ju(ii + len) = ju(ii + k);
435 len = (std::min)(sizeu, nnzU);
436 typename Vector::SegmentReturnType uu(u.segment(ii+1, sizeu-1));
437 typename VectorI::SegmentReturnType juu(ju.segment(ii+1, sizeu-1));
441 for(
Index k = ii + 1; k < ii + len; k++)
442 m_lu.insertBackByOuterInnerUnordered(ii,ju(k)) = u(k);
445 m_lu.makeCompressed();
447 m_factorizationIsOk =
true;
int n
Definition BiCGSTAB_simple.cpp:1
#define eigen_internal_assert(x)
Definition Macros.h:1043
#define EIGEN_NOEXCEPT
Definition Macros.h:1418
#define EIGEN_CONSTEXPR
Definition Macros.h:787
#define eigen_assert(x)
Definition Macros.h:1037
std::vector< int > ind
Definition Slicing_stdvector_cxx11.cpp:1
MatrixXf mat
Definition Tutorial_AdvancedInitialization_CommaTemporary.cpp:1
MatrixXd mat1(size, size)
Scalar * b
Definition benchVecAdd.cpp:17
SCALAR Scalar
Definition bench_gemm.cpp:46
NumTraits< Scalar >::Real RealScalar
Definition bench_gemm.cpp:47
MatrixXf MatrixType
Definition benchmark-blocking-sizes.cpp:52
Incomplete LU factorization with dual-threshold strategy.
Definition IncompleteLUT.h:100
_StorageIndex StorageIndex
Definition IncompleteLUT.h:106
Matrix< StorageIndex, Dynamic, 1 > VectorI
Definition IncompleteLUT.h:109
_Scalar Scalar
Definition IncompleteLUT.h:105
void factorize(const MatrixType &amat)
void setFillfactor(int fillfactor)
Definition IncompleteLUT.h:216
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition IncompleteLUT.h:135
void _solve_impl(const Rhs &b, Dest &x) const
Definition IncompleteLUT.h:171
IncompleteLUT(const MatrixType &mat, const RealScalar &droptol=NumTraits< Scalar >::dummy_precision(), int fillfactor=10)
Definition IncompleteLUT.h:125
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
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_Pinv
Definition IncompleteLUT.h:198
IncompleteLUT()
Definition IncompleteLUT.h:119
ComputationInfo m_info
Definition IncompleteLUT.h:196
SparseSolverBase< IncompleteLUT > Base
Definition IncompleteLUT.h:102
SparseMatrix< Scalar, RowMajor, StorageIndex > FactorType
Definition IncompleteLUT.h:110
ComputationInfo info() const
Reports whether previous computation was successful.
Definition IncompleteLUT.h:142
void setDroptol(const RealScalar &droptol)
Definition IncompleteLUT.h:206
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition IncompleteLUT.h:133
IncompleteLUT & compute(const MatrixType &amat)
Definition IncompleteLUT.h:160
bool m_isInitialized
Definition SparseSolverBase.h:119
Matrix< Scalar, Dynamic, 1 > Vector
Definition IncompleteLUT.h:108
@ ColsAtCompileTime
Definition IncompleteLUT.h:113
@ MaxColsAtCompileTime
Definition IncompleteLUT.h:114
FactorType m_lu
Definition IncompleteLUT.h:191
void analyzePattern(const MatrixType &amat)
NumTraits< Scalar >::Real RealScalar
Definition IncompleteLUT.h:107
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_P
Definition IncompleteLUT.h:197
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
Permutation matrix.
Definition PermutationMatrix.h:298
Definition SparseCompressedBase.h:159
Index col() const
Definition SparseCompressedBase.h:225
StorageIndex index() const
Definition SparseCompressedBase.h:222
const Scalar & value() const
Definition SparseCompressedBase.h:219
TransposeReturnType transpose()
Definition SparseMatrixBase.h:354
SparseSymmetricPermutationProduct< Derived, Upper|Lower > twistedBy(const PermutationMatrix< Dynamic, Dynamic, StorageIndex > &perm) const
Definition SparseMatrixBase.h:329
A versatible sparse matrix representation.
Definition SparseMatrix.h:98
Index rows() const
Definition SparseMatrix.h:138
Index cols() const
Definition SparseMatrix.h:140
A base class for sparse solvers.
Definition SparseSolverBase.h:68
bool m_isInitialized
Definition SparseSolverBase.h:119
#define abs(x)
Definition datatypes.h:17
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition gnuplot_common_settings.hh:12
ComputationInfo
Definition Constants.h:440
@ NumericalIssue
Definition Constants.h:444
@ Success
Definition Constants.h:442
int EIGEN_BLAS_FUNC() swap(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
Definition level1_impl.h:130
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
EIGEN_CONSTEXPR Index first(const T &x) EIGEN_NOEXCEPT
Definition IndexedViewHelper.h:81
void swap(scoped_array< T > &a, scoped_array< T > &b)
Definition Memory.h:709
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition MathFunctions.h:1292
Namespace containing all symbols from the Eigen library.
Definition bench_norm.cpp:85
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
const Product< Lhs, Rhs > prod(const Lhs &lhs, const Rhs &rhs)
Definition evaluators.cpp:8
const int Dynamic
Definition Constants.h:22
Definition BandTriangularSolver.h:13
Definition IncompleteLUT.h:182
bool operator()(const Index &row, const Index &col, const Scalar &) const
Definition IncompleteLUT.h:183
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition NumTraits.h:233
Definition ForwardDeclarations.h:17
std::ptrdiff_t j
Definition tut_arithmetic_redux_minmax.cpp:2