From f0fb95135dfd2d109a793fe8793b13c401f36bf4 Mon Sep 17 00:00:00 2001 From: Gael Guennebaud Date: Mon, 15 Oct 2018 23:47:46 +0200 Subject: Iterative solvers: unify and fix handling of multiple rhs. m_info was not properly computed and the logic was repeated in several places. --- unsupported/Eigen/src/IterativeSolvers/DGMRES.h | 31 +++++++------------------ 1 file changed, 9 insertions(+), 22 deletions(-) (limited to 'unsupported/Eigen/src/IterativeSolvers/DGMRES.h') diff --git a/unsupported/Eigen/src/IterativeSolvers/DGMRES.h b/unsupported/Eigen/src/IterativeSolvers/DGMRES.h index 9ae766995..2ab56b5e7 100644 --- a/unsupported/Eigen/src/IterativeSolvers/DGMRES.h +++ b/unsupported/Eigen/src/IterativeSolvers/DGMRES.h @@ -109,6 +109,7 @@ class DGMRES : public IterativeSolverBase > using Base::m_tolerance; public: using Base::_solve_impl; + using Base::_solve_with_guess_impl; typedef _MatrixType MatrixType; typedef typename MatrixType::Scalar Scalar; typedef typename MatrixType::StorageIndex StorageIndex; @@ -141,30 +142,16 @@ class DGMRES : public IterativeSolverBase > /** \internal */ template - void _solve_with_guess_impl(const Rhs& b, Dest& x) const - { - bool failed = false; - for(Index j=0; j - void _solve_impl(const Rhs& b, MatrixBase& x) const + void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const { - x = b; - _solve_with_guess_impl(b,x.derived()); + EIGEN_STATIC_ASSERT(Rhs::ColsAtCompileTime==1 || Dest::ColsAtCompileTime==1, YOU_TRIED_CALLING_A_VECTOR_METHOD_ON_A_MATRIX); + + m_iterations = Base::maxIterations(); + m_error = Base::m_tolerance; + + dgmres(matrix(), b, x, Base::m_preconditioner); } + /** * Get the restart value */ -- cgit v1.2.3