aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/UmfPackSupport
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2014-09-01 15:00:19 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2014-09-01 15:00:19 +0200
commit85c765957418cd2fd24b46ca3a14e6fcbad43f05 (patch)
tree01111ad07068a8d03fab51c277e0026de2a441cb /Eigen/src/UmfPackSupport
parentbc065c75d2a8df928cb368d0352b8fcb25791fa8 (diff)
Refactoring of sparse solvers through a SparseSolverBase class and usage of the Solve<> expression. Introduce a SolveWithGuess expression on top of Solve.
Diffstat (limited to 'Eigen/src/UmfPackSupport')
-rw-r--r--Eigen/src/UmfPackSupport/UmfPackSupport.h21
1 files changed, 13 insertions, 8 deletions
diff --git a/Eigen/src/UmfPackSupport/UmfPackSupport.h b/Eigen/src/UmfPackSupport/UmfPackSupport.h
index 3a48cecf7..845c8076a 100644
--- a/Eigen/src/UmfPackSupport/UmfPackSupport.h
+++ b/Eigen/src/UmfPackSupport/UmfPackSupport.h
@@ -121,9 +121,13 @@ inline int umfpack_get_determinant(std::complex<double> *Mx, double *Ex, void *N
* \sa \ref TutorialSparseDirectSolvers
*/
template<typename _MatrixType>
-class UmfPackLU : internal::noncopyable
+class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
{
+ protected:
+ typedef SparseSolverBase<UmfPackLU<_MatrixType> > Base;
+ using Base::m_isInitialized;
public:
+ using Base::_solve_impl;
typedef _MatrixType MatrixType;
typedef typename MatrixType::Scalar Scalar;
typedef typename MatrixType::RealScalar RealScalar;
@@ -197,7 +201,8 @@ class UmfPackLU : internal::noncopyable
analyzePattern(matrix);
factorize(matrix);
}
-
+
+#ifndef EIGEN_TEST_EVALUATORS
/** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A.
*
* \sa compute()
@@ -223,6 +228,7 @@ class UmfPackLU : internal::noncopyable
&& "UmfPackLU::solve(): invalid number of rows of the right hand side matrix b");
return internal::sparse_solve_retval<UmfPackLU, Rhs>(*this, b.derived());
}
+#endif
/** Performs a symbolic decomposition on the sparcity of \a matrix.
*
@@ -274,7 +280,7 @@ class UmfPackLU : internal::noncopyable
#ifndef EIGEN_PARSED_BY_DOXYGEN
/** \internal */
template<typename BDerived,typename XDerived>
- bool _solve(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const;
+ bool _solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const;
#endif
Scalar determinant() const;
@@ -328,7 +334,6 @@ class UmfPackLU : internal::noncopyable
void* m_symbolic;
mutable ComputationInfo m_info;
- bool m_isInitialized;
int m_factorizationIsOk;
int m_analysisIsOk;
mutable bool m_extractedDataAreDirty;
@@ -376,7 +381,7 @@ typename UmfPackLU<MatrixType>::Scalar UmfPackLU<MatrixType>::determinant() cons
template<typename MatrixType>
template<typename BDerived,typename XDerived>
-bool UmfPackLU<MatrixType>::_solve(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const
+bool UmfPackLU<MatrixType>::_solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const
{
const int rhsCols = b.cols();
eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet");
@@ -396,7 +401,7 @@ bool UmfPackLU<MatrixType>::_solve(const MatrixBase<BDerived> &b, MatrixBase<XDe
return true;
}
-
+#ifndef EIGEN_TEST_EVALUATORS
namespace internal {
template<typename _MatrixType, typename Rhs>
@@ -408,7 +413,7 @@ struct solve_retval<UmfPackLU<_MatrixType>, Rhs>
template<typename Dest> void evalTo(Dest& dst) const
{
- dec()._solve(rhs(),dst);
+ dec()._solve_impl(rhs(),dst);
}
};
@@ -426,7 +431,7 @@ struct sparse_solve_retval<UmfPackLU<_MatrixType>, Rhs>
};
} // end namespace internal
-
+#endif
} // end namespace Eigen
#endif // EIGEN_UMFPACKSUPPORT_H