25#ifndef DUMUX_PARALLEL_AMGBACKEND_HH
26#define DUMUX_PARALLEL_AMGBACKEND_HH
30#include <dune/common/exceptions.hh>
31#include <dune/common/parallel/indexset.hh>
32#include <dune/common/parallel/mpicommunication.hh>
33#include <dune/common/parallel/mpihelper.hh>
34#include <dune/grid/common/capabilities.hh>
35#include <dune/istl/paamg/amg.hh>
36#include <dune/istl/paamg/pinfo.hh>
37#include <dune/istl/solvers.hh>
49template <
class LinearSolverTraits>
60 , isParallel_(
Dune::MPIHelper::getCommunication().size() > 1)
63 DUNE_THROW(Dune::InvalidStateException,
"Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
65 checkAvailabilityOfDirectSolver_();
76 const typename LinearSolverTraits::DofMapper& dofMapper,
80 , isParallel_(
Dune::MPIHelper::getCommunication().size() > 1)
85 phelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
87 checkAvailabilityOfDirectSolver_();
97 const typename LinearSolverTraits::DofMapper& dofMapper)
101 phelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
112 template<
class Matrix,
class Vector>
113 bool solve(Matrix& A, Vector& x, Vector& b)
116 solveSequentialOrParallel_(A, x, b);
118 solveSequential_(A, x, b);
120 return result_.converged;
128 return "AMG-preconditioned BiCGSTAB solver";
134 const Dune::InverseOperatorResult&
result()
const
141 void checkAvailabilityOfDirectSolver_()
143#if !HAVE_SUPERLU && !HAVE_UMFPACK
144 std::cout <<
"\nAMGBiCGSTABBackend: No direct solver backend found. Using iterative solver as coarse grid solver.\n"
145 <<
"Note that dune-istl currently hard-codes a tolerance of 1e-2 for the iterative coarse grid solver.\n"
146 <<
"This may result in reduced accuracy or performance depending on your setup.\nConsider installing "
147 <<
"UMFPack (SuiteSparse) or SuperLU or apply the istl patch, see dumux/patches/README.md." << std::endl;
152 template<
class Matrix,
class Vector>
153 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
159 if (LinearSolverTraits::isNonOverlapping(phelper_->gridView()))
161 using PTraits =
typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
162 solveParallel_<PTraits>(A, x, b);
166 using PTraits =
typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
167 solveParallel_<PTraits>(A, x, b);
171 solveSequential_(A, x, b);
175 solveSequential_(A, x, b);
179 template<
class ParallelTraits,
class Matrix,
class Vector>
180 void solveParallel_(Matrix& A, Vector& x, Vector& b)
182 using Comm =
typename ParallelTraits::Comm;
183 using LinearOperator =
typename ParallelTraits::LinearOperator;
184 using ScalarProduct =
typename ParallelTraits::ScalarProduct;
186 std::shared_ptr<Comm> comm;
187 std::shared_ptr<LinearOperator> linearOperator;
188 std::shared_ptr<ScalarProduct> scalarProduct;
189 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, comm, linearOperator, scalarProduct, *phelper_);
191 using SeqSmoother = Dune::SeqSSOR<Matrix, Vector, Vector>;
192 using Smoother =
typename ParallelTraits::template Preconditioner<SeqSmoother>;
193 solveWithAmg_<Smoother>(A, x, b, linearOperator, comm, scalarProduct);
197 template<
class Matrix,
class Vector>
198 void solveSequential_(Matrix& A, Vector& x, Vector& b)
200 using Comm = Dune::Amg::SequentialInformation;
201 using Traits =
typename LinearSolverTraits::template Sequential<Matrix, Vector>;
202 using LinearOperator =
typename Traits::LinearOperator;
203 using ScalarProduct =
typename Traits::ScalarProduct;
205 auto comm = std::make_shared<Comm>();
206 auto linearOperator = std::make_shared<LinearOperator>(A);
207 auto scalarProduct = std::make_shared<ScalarProduct>();
209 using Smoother = Dune::SeqSSOR<Matrix, Vector, Vector>;
210 solveWithAmg_<Smoother>(A, x, b, linearOperator, comm, scalarProduct);
213 template<
class Smoother,
class Matrix,
class Vector,
class LinearOperator,
class Comm,
class ScalarProduct>
214 void solveWithAmg_(Matrix& A, Vector& x, Vector& b,
215 std::shared_ptr<LinearOperator>& linearOperator,
216 std::shared_ptr<Comm>& comm,
217 std::shared_ptr<ScalarProduct>& scalarProduct)
219 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
220 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<Matrix, Dune::Amg::FirstDiagonal>>;
224 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
225 params.setDefaultValuesIsotropic(LinearSolverTraits::GridView::dimension);
227 Criterion criterion(params);
228 SmootherArgs smootherArgs;
229 smootherArgs.iterations = 1;
230 smootherArgs.relaxationFactor = 1;
232 using Amg = Dune::Amg::AMG<LinearOperator, Vector, Smoother, Comm>;
233 auto amg = std::make_shared<Amg>(*linearOperator, criterion, smootherArgs, *comm);
235 Dune::BiCGSTABSolver<Vector> solver(*linearOperator, *scalarProduct, *amg, this->
residReduction(), this->
maxIter(),
236 comm->communicator().rank() == 0 ? this->verbosity() : 0);
238 solver.apply(x, b, result_);
242 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> phelper_;
244 Dune::InverseOperatorResult result_;
245 bool isParallel_ =
false;
Provides a helper class for nonoverlapping decomposition.
Base class for linear solvers.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
static constexpr bool canCommunicate
Definition: gridcapabilities.hh:63
Definition: deprecated.hh:149
A linear solver based on the ISTL AMG preconditioner and the ISTL BiCGSTAB solver.
Definition: amgbackend.hh:51
const Dune::InverseOperatorResult & result() const
The result containing the convergence history.
Definition: amgbackend.hh:134
std::string name() const
The name of the solver.
Definition: amgbackend.hh:126
AMGBiCGSTABBackend(const std::string ¶mGroup="")
Construct the backend for the sequential case only.
Definition: amgbackend.hh:58
AMGBiCGSTABBackend(const typename LinearSolverTraits::GridView &gridView, const typename LinearSolverTraits::DofMapper &dofMapper, const std::string ¶mGroup="")
Construct the backend for parallel or sequential runs.
Definition: amgbackend.hh:75
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: amgbackend.hh:113
void updateAfterGridAdaption(const typename LinearSolverTraits::GridView &gridView, const typename LinearSolverTraits::DofMapper &dofMapper)
Update the solver after grid adaption.
Definition: amgbackend.hh:96
Base class for linear solvers.
Definition: solver.hh:37
Scalar residReduction() const
the linear solver residual reduction
Definition: solver.hh:99
int maxIter() const
the maximum number of linear solver iterations
Definition: solver.hh:91
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:79
int verbosity() const
the verbosity level
Definition: solver.hh:83