14#ifndef DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
15#define DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
20#include <dune/common/version.hh>
21#include <dune/common/parallel/mpihelper.hh>
22#include <dune/common/parametertree.hh>
27#include <dune/istl/preconditioners.hh>
28#include <dune/istl/paamg/amg.hh>
32#include <dune/istl/solvers.hh>
33#include <dune/istl/solverfactory.hh>
57template<
class LinearOperator>
58int initSolverFactoriesForMultiTypeBlockMatrix()
60#if DUNE_VERSION_LT(DUNE_ISTL,2,11)
61 using M =
typename LinearOperator::matrix_type;
62 using X =
typename LinearOperator::range_type;
63 using Y =
typename LinearOperator::domain_type;
64 using TL = Dune::TypeList<M,X,Y>;
65 auto& dsfac = Dune::DirectSolverFactory<M,X,Y>::instance();
66 Dune::addRegistryToFactory<TL>(dsfac, Dumux::MultiTypeBlockMatrixDirectSolverTag{});
67 auto& pfac = Dune::PreconditionerFactory<LinearOperator,X,Y>::instance();
68 Dune::addRegistryToFactory<TL>(pfac, Dumux::MultiTypeBlockMatrixPreconditionerTag{});
69 using TLS = Dune::TypeList<X,Y>;
70 auto& isfac = Dune::IterativeSolverFactory<X,Y>::instance();
71 return Dune::addRegistryToFactory<TLS>(isfac, Dune::IterativeSolverTag{});
73 using OpTraits = Dune::OperatorTraits<LinearOperator>;
74 auto& pfac = Dune::PreconditionerFactory<LinearOperator>::instance();
75 Dune::addRegistryToFactory<OpTraits>(pfac, Dumux::MultiTypeBlockMatrixPreconditionerTag{});
76 auto& sfac = Dune::SolverFactory<LinearOperator>::instance();
77 return Dune::addRegistryToFactory<OpTraits>(sfac, Dumux::MultiTypeBlockMatrixSolverTag{});
89template<
class Matrix,
class LinearOperator>
93 initSolverFactoriesForMultiTypeBlockMatrix<LinearOperator>();
95 Dune::initSolverFactories<LinearOperator>();
104template <
class LinearSolverTraits,
class LinearAlgebraTraits>
109 using ScalarProduct = Dune::ScalarProduct<Vector>;
111 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>,
int>;
122 , isParallel_(
Dune::MPIHelper::getCommunication().size() > 1)
125 DUNE_THROW(Dune::InvalidStateException,
"Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
128 initializeParameters_();
129 scalarProduct_ = std::make_shared<ScalarProduct>();
130 solverCategory_ = Dune::SolverCategory::sequential;
141 const typename LinearSolverTraits::DofMapper& dofMapper,
145 , isParallel_(gridView.comm().size() > 1)
149 initializeParameters_();
152 if (solverCategory_ != Dune::SolverCategory::sequential)
154 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
155 comm_ = std::make_shared<Comm>(gridView.comm(), solverCategory_);
156 scalarProduct_ = Dune::createScalarProduct<Vector>(*comm_, solverCategory_);
157 parallelHelper_->createParallelIndexSet(*comm_);
160 scalarProduct_ = std::make_shared<ScalarProduct>();
162 solverCategory_ = Dune::SolverCategory::sequential;
163 scalarProduct_ = std::make_shared<ScalarProduct>();
174 bool solve(Matrix& A, Vector& x, Vector& b)
177 solveSequentialOrParallel_(A, x, b);
179 solveSequential_(A, x, b);
182 return result_.converged;
185 const Dune::InverseOperatorResult&
result()
const
190 const std::string&
name()
const
200 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
203 using GV =
typename LinearSolverTraits::GridView;
204 using DM =
typename LinearSolverTraits::DofMapper;
205 if constexpr (
requires { LinearSolverTraits::dofCodims; })
215 return scalarProduct_->norm(y);
219 return scalarProduct_->norm(x);
224 void initializeParameters_()
227 checkMandatoryParameters_();
228 name_ = params_.get<std::string>(
"preconditioner.type") +
"-preconditioned " + params_.get<std::string>(
"type");
229 if (params_.get<
int>(
"verbose", 0) > 0)
230 std::cout <<
"Initialized linear solver of type: " << name_ << std::endl;
233 void checkMandatoryParameters_()
235 if (!params_.hasKey(
"type"))
236 DUNE_THROW(Dune::InvalidStateException,
"Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
238 if (!params_.hasKey(
"preconditioner.type"))
239 DUNE_THROW(Dune::InvalidStateException,
"Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
243 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
248 if constexpr (LinearSolverTraits::canCommunicate && !isMultiTypeBlockMatrix<Matrix>::value)
252 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
254 using PTraits =
typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
255 solveParallel_<PTraits>(A, x, b);
259 using PTraits =
typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
260 solveParallel_<PTraits>(A, x, b);
264 solveSequential_(A, x, b);
268 solveSequential_(A, x, b);
272 template<
class ParallelTraits>
273 void solveParallel_(Matrix& A, Vector& x, Vector& b)
275 using LinearOperator =
typename ParallelTraits::LinearOperator;
281 auto linearOperator = std::make_shared<LinearOperator>(A, *comm_);
284 apply_(linearOperator, x, b);
288 void solveSequential_(Matrix& A, Vector& x, Vector& b)
291 using Traits =
typename LinearSolverTraits::template Sequential<Matrix, Vector>;
292 using LinearOperator =
typename Traits::LinearOperator;
293 auto linearOperator = std::make_shared<LinearOperator>(A);
299 apply_(linearOperator, x, b);
302 template<
class LinearOperator>
303 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
306 return Dune::getSolverFromFactory(fop, params_);
307 }
catch(Dune::Exception& e) {
308 std::cerr <<
"Dune::Exception during solver construction: " << e.what() << std::endl;
309 return std::decay_t<decltype(Dune::getSolverFromFactory(fop, params_))>();
313 template<
class LinearOperator>
314 void apply_(std::shared_ptr<LinearOperator>& fop, Vector& x, Vector& b)
316 auto solver = getSolverFromFactory_(fop);
320 bool success =
static_cast<bool>(solver);
322 int successRemote = success;
324 successRemote = comm_->communicator().min(success);
327 DUNE_THROW(Dune::Exception,
"Could not create ISTL solver");
328 else if (!successRemote)
329 DUNE_THROW(Dune::Exception,
"Could not create ISTL solver on remote process");
332 DUNE_THROW(Dune::Exception,
"Could not create ISTL solver");
336 solver->apply(x, b, result_);
339 const std::string paramGroup_;
341 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
342 std::shared_ptr<Comm> comm_;
344 bool isParallel_ =
false;
347 std::shared_ptr<ScalarProduct> scalarProduct_;
348 Dune::SolverCategory::Category solverCategory_;
349 Dune::InverseOperatorResult result_;
350 Dune::ParameterTree params_;
const std::string & name() const
Definition istlsolverfactorybackend.hh:190
IstlSolverFactoryBackend(const std::string ¶mGroup="")
Construct the backend for the sequential case only.
Definition istlsolverfactorybackend.hh:120
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:140
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition istlsolverfactorybackend.hh:174
Scalar norm(const Vector &x) const
Definition istlsolverfactorybackend.hh:195
const Dune::InverseOperatorResult & result() const
Definition istlsolverfactorybackend.hh:185
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition solver.hh:78
LinearSolver(const std::string ¶mGroup="")
Construct the solver.
Definition solver.hh:43
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:580
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v, const std::bitset< numCodims > &activeCodims) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition parallelhelpers.hh:588
Definition parallelhelpers.hh:522
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition parallelhelpers.hh:530
The specialized Dumux macro and tag for the ISTL registry to choose the solver and preconditioner at ...
The specialized Dumux tag for the ISTL registry including only ISTL solvers compatible with MultiType...
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Type traits to be used with matrix types.
Dune::SolverCategory::Category solverCategory(const GridView &gridView)
Definition solvercategory.hh:20
void prepareLinearAlgebraParallel(Matrix &A, Vector &b, ParallelHelper &pHelper)
Prepare linear algebra variables for parallel solvers.
Definition parallelhelpers.hh:1477
void initSolverFactories()
Initialize the solver factories for regular matrices or MultiTypeBlockMatrices.
Definition istlsolverfactorybackend.hh:90
Definition common/pdesolver.hh:24
Provides a helper class for nonoverlapping decomposition.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
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.