26#ifndef DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
27#define DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
31#include <dune/common/parallel/mpihelper.hh>
32#include <dune/common/parametertree.hh>
37#include <dune/istl/preconditioners.hh>
38#include <dune/istl/paamg/amg.hh>
42#include <dune/istl/solvers.hh>
43#include <dune/istl/solverfactory.hh>
64template<
class LinearOperator>
65int initSolverFactoriesForMultiTypeBlockMatrix()
67 using M =
typename LinearOperator::matrix_type;
68 using X =
typename LinearOperator::range_type;
69 using Y =
typename LinearOperator::domain_type;
70 using TL = Dune::TypeList<M,X,Y>;
71 auto& dsfac = Dune::DirectSolverFactory<M,X,Y>::instance();
72 Dune::addRegistryToFactory<TL>(dsfac, Dumux::MultiTypeBlockMatrixDirectSolverTag{});
73 auto& pfac = Dune::PreconditionerFactory<LinearOperator,X,Y>::instance();
74 Dune::addRegistryToFactory<TL>(pfac, Dumux::MultiTypeBlockMatrixPreconditionerTag{});
75 using TLS = Dune::TypeList<X,Y>;
76 auto& isfac = Dune::IterativeSolverFactory<X,Y>::instance();
77 return Dune::addRegistryToFactory<TLS>(isfac, Dune::IterativeSolverTag{});
88template<
class Matrix,
class LinearOperator>
92 initSolverFactoriesForMultiTypeBlockMatrix<LinearOperator>();
94 Dune::initSolverFactories<LinearOperator>();
103template <
class LinearSolverTraits>
115 , isParallel_(
Dune::MPIHelper::getCommunication().size() > 1)
118 DUNE_THROW(Dune::InvalidStateException,
"Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
121 initializeParameters_();
132 const typename LinearSolverTraits::DofMapper& dofMapper,
136 , isParallel_(
Dune::MPIHelper::getCommunication().size() > 1)
140 initializeParameters_();
143 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
154 const typename LinearSolverTraits::DofMapper& dofMapper)
158 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
169 template<
class Matrix,
class Vector>
170 bool solve(Matrix& A, Vector& x, Vector& b)
173 solveSequentialOrParallel_(A, x, b);
175 solveSequential_(A, x, b);
178 return result_.converged;
181 const Dune::InverseOperatorResult&
result()
const
186 const std::string&
name()
const
193 void initializeParameters_()
196 checkMandatoryParameters_();
197 name_ = params_.get<std::string>(
"preconditioner.type") +
"-preconditioned " + params_.get<std::string>(
"type");
198 if (params_.get<
int>(
"verbose", 0) > 0)
199 std::cout <<
"Initialized linear solver of type: " << name_ << std::endl;
202 void checkMandatoryParameters_()
204 if (!params_.hasKey(
"type"))
205 DUNE_THROW(Dune::InvalidStateException,
"Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
207 if (!params_.hasKey(
"preconditioner.type"))
208 DUNE_THROW(Dune::InvalidStateException,
"Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
212 template<
class Matrix,
class Vector>
213 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
222 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
224 using PTraits =
typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
225 solveParallel_<PTraits>(A, x, b);
229 using PTraits =
typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
230 solveParallel_<PTraits>(A, x, b);
234 solveSequential_(A, x, b);
238 solveSequential_(A, x, b);
242 template<
class ParallelTraits,
class Matrix,
class Vector>
243 void solveParallel_(Matrix& A, Vector& x, Vector& b)
245 using Comm =
typename ParallelTraits::Comm;
246 using LinearOperator =
typename ParallelTraits::LinearOperator;
247 using ScalarProduct =
typename ParallelTraits::ScalarProduct;
250 initSolverFactories<Matrix, LinearOperator>();
252 std::shared_ptr<Comm> comm;
253 std::shared_ptr<LinearOperator> linearOperator;
254 std::shared_ptr<ScalarProduct> scalarProduct;
255 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, comm, linearOperator, scalarProduct, *parallelHelper_);
258 auto solver = getSolverFromFactory_(linearOperator);
261 solver->apply(x, b, result_);
265 template<
class Matrix,
class Vector>
266 void solveSequential_(Matrix& A, Vector& x, Vector& b)
269 using Traits =
typename LinearSolverTraits::template Sequential<Matrix, Vector>;
270 using LinearOperator =
typename Traits::LinearOperator;
271 auto linearOperator = std::make_shared<LinearOperator>(A);
274 initSolverFactories<Matrix, LinearOperator>();
277 auto solver = getSolverFromFactory_(linearOperator);
280 solver->apply(x, b, result_);
283 template<
class LinearOperator>
284 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
286 try {
return Dune::getSolverFromFactory(fop, params_); }
287 catch(Dune::Exception& e)
289 std::cerr <<
"Could not create solver with factory" << std::endl;
290 std::cerr << e.what() << std::endl;
295 const std::string paramGroup_;
297 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
299 bool isParallel_ =
false;
302 Dune::InverseOperatorResult result_;
303 Dune::ParameterTree params_;
Type traits to be used with matrix types.
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
The specialized Dumux macro and tag for the ISTL registry to choose the solver and preconditioner at ...
Provides a helper class for nonoverlapping decomposition.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
void initSolverFactories()
Initialize the solver factories for regular matrices or MultiTypeBlockMatrices.
Definition: istlsolverfactorybackend.hh:89
static constexpr bool canCommunicate
Definition: gridcapabilities.hh:63
Definition: deprecated.hh:149
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:49
A linear solver using the dune-istl solver factory to choose the solver and preconditioner at runtime...
Definition: istlsolverfactorybackend.hh:105
const std::string & name() const
Definition: istlsolverfactorybackend.hh:186
IstlSolverFactoryBackend(const typename LinearSolverTraits::GridView &gridView, const typename LinearSolverTraits::DofMapper &dofMapper, const std::string ¶mGroup="")
Construct the backend for parallel or sequential runs.
Definition: istlsolverfactorybackend.hh:131
IstlSolverFactoryBackend(const std::string ¶mGroup="")
Construct the backend for the sequential case only.
Definition: istlsolverfactorybackend.hh:113
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: istlsolverfactorybackend.hh:170
const Dune::InverseOperatorResult & result() const
Definition: istlsolverfactorybackend.hh:181
void updateAfterGridAdaption(const typename LinearSolverTraits::GridView &gridView, const typename LinearSolverTraits::DofMapper &dofMapper)
Update the solver after grid adaption.
Definition: istlsolverfactorybackend.hh:153
static Dune::ParameterTree createParameterTree(const std::string ¶mGroup="")
Create a tree containing parameters required for the linear solvers and precondioners of the Dune IST...
Definition: linearsolverparameters.hh:49
Base class for linear solvers.
Definition: solver.hh:37
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:79