aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/IterativeLinearSolvers/IncompleteLUT.h
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2014-09-01 17:19:16 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2014-09-01 17:19:16 +0200
commit8a74ce922cbee1c416675e57e48393dd6b399e4e (patch)
tree96facd92fce2ed7bce9935a3e19dbc24ec143cac /Eigen/src/IterativeLinearSolvers/IncompleteLUT.h
parent863b7362bc6cbd6f5c3ba96b4d6935557f9d9486 (diff)
Make IncompleteLUT use SparseSolverBase.
Diffstat (limited to 'Eigen/src/IterativeLinearSolvers/IncompleteLUT.h')
-rw-r--r--Eigen/src/IterativeLinearSolvers/IncompleteLUT.h26
1 files changed, 8 insertions, 18 deletions
diff --git a/Eigen/src/IterativeLinearSolvers/IncompleteLUT.h b/Eigen/src/IterativeLinearSolvers/IncompleteLUT.h
index a39ed7748..f337c5fb0 100644
--- a/Eigen/src/IterativeLinearSolvers/IncompleteLUT.h
+++ b/Eigen/src/IterativeLinearSolvers/IncompleteLUT.h
@@ -94,8 +94,12 @@ Index QuickSplit(VectorV &row, VectorI &ind, Index ncut)
* http://comments.gmane.org/gmane.comp.lib.eigen/3302
*/
template <typename _Scalar>
-class IncompleteLUT : internal::noncopyable
+class IncompleteLUT : public SparseSolverBase<IncompleteLUT<_Scalar> >
{
+ protected:
+ typedef SparseSolverBase<IncompleteLUT<_Scalar> > Base;
+ using Base::m_isInitialized;
+ public:
typedef _Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar,Dynamic,1> Vector;
@@ -108,13 +112,13 @@ class IncompleteLUT : internal::noncopyable
IncompleteLUT()
: m_droptol(NumTraits<Scalar>::dummy_precision()), m_fillfactor(10),
- m_analysisIsOk(false), m_factorizationIsOk(false), m_isInitialized(false)
+ m_analysisIsOk(false), m_factorizationIsOk(false)
{}
template<typename MatrixType>
IncompleteLUT(const MatrixType& mat, const RealScalar& droptol=NumTraits<Scalar>::dummy_precision(), int fillfactor = 10)
: m_droptol(droptol),m_fillfactor(fillfactor),
- m_analysisIsOk(false),m_factorizationIsOk(false),m_isInitialized(false)
+ m_analysisIsOk(false),m_factorizationIsOk(false)
{
eigen_assert(fillfactor != 0);
compute(mat);
@@ -159,11 +163,7 @@ class IncompleteLUT : internal::noncopyable
void setFillfactor(int fillfactor);
template<typename Rhs, typename Dest>
-#ifndef EIGEN_TEST_EVALUATORS
- void _solve(const Rhs& b, Dest& x) const
-#else
void _solve_impl(const Rhs& b, Dest& x) const
-#endif
{
x = m_Pinv * b;
x = m_lu.template triangularView<UnitLower>().solve(x);
@@ -180,15 +180,6 @@ class IncompleteLUT : internal::noncopyable
&& "IncompleteLUT::solve(): invalid number of rows of the right hand side matrix b");
return internal::solve_retval<IncompleteLUT, Rhs>(*this, b.derived());
}
-#else
- template<typename Rhs> inline const Solve<IncompleteLUT, Rhs>
- solve(const MatrixBase<Rhs>& b) const
- {
- eigen_assert(m_isInitialized && "IncompleteLUT is not initialized.");
- eigen_assert(cols()==b.rows()
- && "IncompleteLUT::solve(): invalid number of rows of the right hand side matrix b");
- return Solve<IncompleteLUT, Rhs>(*this, b.derived());
- }
#endif
protected:
@@ -208,7 +199,6 @@ protected:
int m_fillfactor;
bool m_analysisIsOk;
bool m_factorizationIsOk;
- bool m_isInitialized;
ComputationInfo m_info;
PermutationMatrix<Dynamic,Dynamic,Index> m_P; // Fill-reducing permutation
PermutationMatrix<Dynamic,Dynamic,Index> m_Pinv; // Inverse permutation
@@ -473,7 +463,7 @@ struct solve_retval<IncompleteLUT<_MatrixType>, Rhs>
template<typename Dest> void evalTo(Dest& dst) const
{
- dec()._solve(rhs(),dst);
+ dec()._solve_impl(rhs(),dst);
}
};