From c409d892e5fb31e5b4e89b40355bbb448633fd74 Mon Sep 17 00:00:00 2001 From: Kilian <kilian.weishaupt@iws.uni-stuttgart.de> Date: Mon, 30 Mar 2020 15:39:00 +0200 Subject: [PATCH] [linear][seqsolverbackend] Add UzawaBiCGSTABBackend and solveWithParamTree function --- dumux/linear/seqsolverbackend.hh | 54 ++++++++++++++++++++++++++++++++ 1 file changed, 54 insertions(+) diff --git a/dumux/linear/seqsolverbackend.hh b/dumux/linear/seqsolverbackend.hh index d06a0a737f..00db276e88 100644 --- a/dumux/linear/seqsolverbackend.hh +++ b/dumux/linear/seqsolverbackend.hh @@ -32,6 +32,8 @@ #include <dune/istl/solvers.hh> #include <dune/istl/superlu.hh> #include <dune/istl/umfpack.hh> +#include <dune/istl/io.hh> +#include <dune/common/indices.hh> #include <dune/common/version.hh> #include <dune/common/hybridutilities.hh> @@ -40,6 +42,8 @@ #include <dumux/common/typetraits/utility.hh> #include <dumux/linear/solver.hh> #include <dumux/linear/amgbackend.hh> +#include <dumux/linear/preconditioners.hh> +#include <dumux/linear/linearsolverparameters.hh> namespace Dumux { @@ -147,6 +151,29 @@ public: return result.converged; } + +#if DUNE_VERSION_GTE(DUNE_ISTL,2,7) + // solve with generic parameter tree + template<class Preconditioner, class Solver, class Matrix, class Vector> + static bool solveWithParamTree(const Matrix& A, Vector& x, const Vector& b, + const Dune::ParameterTree& params) + { + + // make a linear operator from a matrix + using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>; + const auto linearOperator = std::make_shared<MatrixAdapter>(A); + + auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub("preconditioner")); + Solver solver(linearOperator, precond, params); + + Vector bTmp(b); + + Dune::InverseOperatorResult result; + solver.apply(x, bTmp, result); + + return result.converged; + } +#endif }; /*! @@ -876,6 +903,33 @@ private: */ // \{ +#if DUNE_VERSION_GTE(DUNE_ISTL,2,7) +/*! + * \ingroup Linear + * \brief A Uzawa preconditioned BiCGSTAB solver for saddle-point problems + */ +template <class LinearSolverTraits> +class UzawaBiCGSTABBackend : public LinearSolver +{ +public: + using LinearSolver::LinearSolver; + + template<class Matrix, class Vector> + bool solve(const Matrix& A, Vector& x, const Vector& b) + { + using Preconditioner = SeqUzawa<Matrix, Vector, Vector>; + using Solver = Dune::BiCGSTABSolver<Vector>; + static const auto solverParams = LinearSolverParameters<LinearSolverTraits>::createParameterTree(this->paramGroup()); + return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams); + } + + std::string name() const + { + return "Uzawa preconditioned BiCGSTAB solver"; + } +}; +#endif + /*! * \ingroup Linear * \brief A simple ilu0 block diagonal preconditioner -- GitLab