From 85c765957418cd2fd24b46ca3a14e6fcbad43f05 Mon Sep 17 00:00:00 2001 From: Gael Guennebaud Date: Mon, 1 Sep 2014 15:00:19 +0200 Subject: Refactoring of sparse solvers through a SparseSolverBase class and usage of the Solve<> expression. Introduce a SolveWithGuess expression on top of Solve. --- Eigen/src/UmfPackSupport/UmfPackSupport.h | 21 +++++++++++++-------- 1 file changed, 13 insertions(+), 8 deletions(-) (limited to 'Eigen/src/UmfPackSupport') 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 *Mx, double *Ex, void *N * \sa \ref TutorialSparseDirectSolvers */ template -class UmfPackLU : internal::noncopyable +class UmfPackLU : public SparseSolverBase > { + protected: + typedef SparseSolverBase > 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(*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 - bool _solve(const MatrixBase &b, MatrixBase &x) const; + bool _solve_impl(const MatrixBase &b, MatrixBase &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::Scalar UmfPackLU::determinant() cons template template -bool UmfPackLU::_solve(const MatrixBase &b, MatrixBase &x) const +bool UmfPackLU::_solve_impl(const MatrixBase &b, MatrixBase &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::_solve(const MatrixBase &b, MatrixBase @@ -408,7 +413,7 @@ struct solve_retval, Rhs> template void evalTo(Dest& dst) const { - dec()._solve(rhs(),dst); + dec()._solve_impl(rhs(),dst); } }; @@ -426,7 +431,7 @@ struct sparse_solve_retval, Rhs> }; } // end namespace internal - +#endif } // end namespace Eigen #endif // EIGEN_UMFPACKSUPPORT_H -- cgit v1.2.3