14#ifndef DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
15#define DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
19#include <dune/common/parallel/mpihelper.hh>
20#include <dune/common/parametertree.hh>
25#include <dune/istl/preconditioners.hh>
26#include <dune/istl/paamg/amg.hh>
30#include <dune/istl/solvers.hh>
31#include <dune/istl/solverfactory.hh>
54template<
class LinearOperator>
55int initSolverFactoriesForMultiTypeBlockMatrix()
57 using M =
typename LinearOperator::matrix_type;
58 using X =
typename LinearOperator::range_type;
59 using Y =
typename LinearOperator::domain_type;
60 using TL = Dune::TypeList<M,X,Y>;
61 auto& dsfac = Dune::DirectSolverFactory<M,X,Y>::instance();
62 Dune::addRegistryToFactory<TL>(dsfac, Dumux::MultiTypeBlockMatrixDirectSolverTag{});
63 auto& pfac = Dune::PreconditionerFactory<LinearOperator,X,Y>::instance();
64 Dune::addRegistryToFactory<TL>(pfac, Dumux::MultiTypeBlockMatrixPreconditionerTag{});
65 using TLS = Dune::TypeList<X,Y>;
66 auto& isfac = Dune::IterativeSolverFactory<X,Y>::instance();
67 return Dune::addRegistryToFactory<TLS>(isfac, Dune::IterativeSolverTag{});
78template<
class Matrix,
class LinearOperator>
82 initSolverFactoriesForMultiTypeBlockMatrix<LinearOperator>();
84 Dune::initSolverFactories<LinearOperator>();
95template <
class LinearSolverTraits>
107 const typename LinearSolverTraits::DofMapper& dofMapper)
111 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
122 template<
class Matrix,
class Vector>
123 bool solve(Matrix& A, Vector& x, Vector& b)
126 solveSequentialOrParallel_(A, x, b);
128 solveSequential_(A, x, b);
131 return result_.converged;
134 const Dune::InverseOperatorResult&
result()
const
139 const std::string&
name()
const
146 void initializeParameters_()
149 checkMandatoryParameters_();
150 name_ = params_.get<std::string>(
"preconditioner.type") +
"-preconditioned " + params_.get<std::string>(
"type");
151 if (params_.get<
int>(
"verbose", 0) > 0)
152 std::cout <<
"Initialized linear solver of type: " << name_ << std::endl;
155 void checkMandatoryParameters_()
157 if (!params_.hasKey(
"type"))
158 DUNE_THROW(Dune::InvalidStateException,
"Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
160 if (!params_.hasKey(
"preconditioner.type"))
161 DUNE_THROW(Dune::InvalidStateException,
"Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
165 template<
class Matrix,
class Vector>
166 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
175 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
177 using PTraits =
typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
178 solveParallel_<PTraits>(A, x, b);
182 using PTraits =
typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
183 solveParallel_<PTraits>(A, x, b);
187 solveSequential_(A, x, b);
191 solveSequential_(A, x, b);
195 template<
class ParallelTraits,
class Matrix,
class Vector>
196 void solveParallel_(Matrix& A, Vector& x, Vector& b)
198 using Comm =
typename ParallelTraits::Comm;
199 using LinearOperator =
typename ParallelTraits::LinearOperator;
200 using ScalarProduct =
typename ParallelTraits::ScalarProduct;
203 initSolverFactories<Matrix, LinearOperator>();
205 std::shared_ptr<Comm> comm;
206 std::shared_ptr<LinearOperator> linearOperator;
207 std::shared_ptr<ScalarProduct> scalarProduct;
208 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, comm, linearOperator, scalarProduct, *parallelHelper_);
211 auto solver = getSolverFromFactory_(linearOperator);
214 solver->apply(x, b, result_);
218 template<
class Matrix,
class Vector>
219 void solveSequential_(Matrix& A, Vector& x, Vector& b)
222 using Traits =
typename LinearSolverTraits::template Sequential<Matrix, Vector>;
223 using LinearOperator =
typename Traits::LinearOperator;
224 auto linearOperator = std::make_shared<LinearOperator>(A);
227 initSolverFactories<Matrix, LinearOperator>();
230 auto solver = getSolverFromFactory_(linearOperator);
233 solver->apply(x, b, result_);
236 template<
class LinearOperator>
237 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
239 try {
return Dune::getSolverFromFactory(fop, params_); }
240 catch(Dune::Exception& e)
242 std::cerr <<
"Could not create solver with factory" << std::endl;
243 std::cerr << e.what() << std::endl;
248 const std::string paramGroup_;
250 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
252 bool isParallel_ =
false;
255 Dune::InverseOperatorResult result_;
256 Dune::ParameterTree params_;
260template <
class LinearSolverTraits,
class LinearAlgebraTraits>
265 using ScalarProduct = Dune::ScalarProduct<Vector>;
267 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>,
int>;
278 , isParallel_(
Dune::MPIHelper::getCommunication().size() > 1)
281 DUNE_THROW(Dune::InvalidStateException,
"Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
284 initializeParameters_();
285 scalarProduct_ = std::make_shared<ScalarProduct>();
286 solverCategory_ = Dune::SolverCategory::sequential;
297 const typename LinearSolverTraits::DofMapper& dofMapper,
301 , isParallel_(gridView.comm().size() > 1)
305 initializeParameters_();
307 solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView);
308 if (solverCategory_ != Dune::SolverCategory::sequential)
310 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
311 comm_ = std::make_shared<Comm>(gridView.comm(), solverCategory_);
312 scalarProduct_ = Dune::createScalarProduct<Vector>(*comm_, solverCategory_);
313 parallelHelper_->createParallelIndexSet(*comm_);
316 scalarProduct_ = std::make_shared<ScalarProduct>();
318 solverCategory_ = Dune::SolverCategory::sequential;
319 scalarProduct_ = std::make_shared<ScalarProduct>();
330 bool solve(Matrix& A, Vector& x, Vector& b)
333 solveSequentialOrParallel_(A, x, b);
335 solveSequential_(A, x, b);
338 return result_.converged;
341 const Dune::InverseOperatorResult&
result()
const
346 const std::string&
name()
const
356 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
359 using GV =
typename LinearSolverTraits::GridView;
360 using DM =
typename LinearSolverTraits::DofMapper;
363 return scalarProduct_->norm(y);
367 return scalarProduct_->norm(x);
372 void initializeParameters_()
375 checkMandatoryParameters_();
376 name_ = params_.get<std::string>(
"preconditioner.type") +
"-preconditioned " + params_.get<std::string>(
"type");
377 if (params_.get<
int>(
"verbose", 0) > 0)
378 std::cout <<
"Initialized linear solver of type: " << name_ << std::endl;
381 void checkMandatoryParameters_()
383 if (!params_.hasKey(
"type"))
384 DUNE_THROW(Dune::InvalidStateException,
"Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
386 if (!params_.hasKey(
"preconditioner.type"))
387 DUNE_THROW(Dune::InvalidStateException,
"Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
391 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
400 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
402 using PTraits =
typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
403 solveParallel_<PTraits>(A, x, b);
407 using PTraits =
typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
408 solveParallel_<PTraits>(A, x, b);
412 solveSequential_(A, x, b);
416 solveSequential_(A, x, b);
420 template<
class ParallelTraits>
421 void solveParallel_(Matrix& A, Vector& x, Vector& b)
423 using LinearOperator =
typename ParallelTraits::LinearOperator;
426 initSolverFactories<Matrix, LinearOperator>();
428 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, *parallelHelper_);
429 auto linearOperator = std::make_shared<LinearOperator>(A, *comm_);
432 auto solver = getSolverFromFactory_(linearOperator);
435 solver->apply(x, b, result_);
439 void solveSequential_(Matrix& A, Vector& x, Vector& b)
442 using Traits =
typename LinearSolverTraits::template Sequential<Matrix, Vector>;
443 using LinearOperator =
typename Traits::LinearOperator;
444 auto linearOperator = std::make_shared<LinearOperator>(A);
447 initSolverFactories<Matrix, LinearOperator>();
450 auto solver = getSolverFromFactory_(linearOperator);
453 solver->apply(x, b, result_);
456 template<
class LinearOperator>
457 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
459 try {
return Dune::getSolverFromFactory(fop, params_); }
460 catch(Dune::Exception& e)
462 std::cerr <<
"Could not create solver with factory" << std::endl;
463 std::cerr << e.what() << std::endl;
468 const std::string paramGroup_;
470 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
471 std::shared_ptr<Comm> comm_;
473 bool isParallel_ =
false;
476 std::shared_ptr<ScalarProduct> scalarProduct_;
477 Dune::SolverCategory::Category solverCategory_;
478 Dune::InverseOperatorResult result_;
479 Dune::ParameterTree params_;
Definition: istlsolverfactorybackend.hh:262
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:296
const std::string & name() const
Definition: istlsolverfactorybackend.hh:346
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: istlsolverfactorybackend.hh:330
const Dune::InverseOperatorResult & result() const
Definition: istlsolverfactorybackend.hh:341
IstlSolverFactoryBackend(const std::string ¶mGroup="")
Construct the backend for the sequential case only.
Definition: istlsolverfactorybackend.hh:276
Scalar norm(const Vector &x) const
Definition: istlsolverfactorybackend.hh:351
Definition: istlsolverfactorybackend.hh:97
const std::string & name() const
Definition: istlsolverfactorybackend.hh:139
const Dune::InverseOperatorResult & result() const
Definition: istlsolverfactorybackend.hh:134
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: istlsolverfactorybackend.hh:123
void updateAfterGridAdaption(const typename LinearSolverTraits::GridView &gridView, const typename LinearSolverTraits::DofMapper &dofMapper)
Update the solver after grid adaption.
Definition: istlsolverfactorybackend.hh:106
Base class for linear solvers.
Definition: solver.hh:27
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:78
double Scalar
Definition: solver.hh:31
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:47
Definition: parallelhelpers.hh:476
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition: parallelhelpers.hh:484
typename Detail::IstlSolverFactoryBackendImplHelper< sizeof...(T)>::template type< T... > IstlSolverFactoryBackend
A linear solver using the dune-istl solver factory to choose the solver and preconditioner at runtime...
Definition: istlsolverfactorybackend.hh:507
The specialized Dumux macro and tag for the ISTL registry to choose the solver and preconditioner at ...
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Type traits to be used with matrix types.
static constexpr bool canCommunicate
Definition: gridcapabilities.hh:51
void initSolverFactories()
Initialize the solver factories for regular matrices or MultiTypeBlockMatrices.
Definition: istlsolverfactorybackend.hh:79
Definition: common/pdesolver.hh:24
Provides a helper class for nonoverlapping decomposition.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
Definition: istlsolverfactorybackend.hh:485
V Vector
Definition: linearalgebratraits.hh:53
M Matrix
Definition: linearalgebratraits.hh:52
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:37
Helper type to determine whether a given type is a Dune::MultiTypeBlockVector.
Definition: vector.hh:22
Type traits to be used with vector types.