version 3.10-dev
istlsolverfactorybackend.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
15#define DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
16
17#include <memory>
18
19#include <dune/common/version.hh>
20#include <dune/common/parallel/mpihelper.hh>
21#include <dune/common/parametertree.hh>
22
24
25// preconditioners
26#include <dune/istl/preconditioners.hh>
27#include <dune/istl/paamg/amg.hh>
28#include "preconditioners.hh"
29
30// solvers
31#include <dune/istl/solvers.hh>
32#include <dune/istl/solverfactory.hh>
33
41
42namespace Dumux {
43
44// initSolverFactoriesForMultiTypeBlockMatrix differs in different compilation units,
45// so we have it in an anonymous namespace
46namespace {
47
56template<class LinearOperator>
57int initSolverFactoriesForMultiTypeBlockMatrix()
58{
59#if DUNE_VERSION_LT(DUNE_ISTL,2,11)
60 using M = typename LinearOperator::matrix_type;
61 using X = typename LinearOperator::range_type;
62 using Y = typename LinearOperator::domain_type;
63 using TL = Dune::TypeList<M,X,Y>;
64 auto& dsfac = Dune::DirectSolverFactory<M,X,Y>::instance();
65 Dune::addRegistryToFactory<TL>(dsfac, Dumux::MultiTypeBlockMatrixDirectSolverTag{});
66 auto& pfac = Dune::PreconditionerFactory<LinearOperator,X,Y>::instance();
67 Dune::addRegistryToFactory<TL>(pfac, Dumux::MultiTypeBlockMatrixPreconditionerTag{});
68 using TLS = Dune::TypeList<X,Y>;
69 auto& isfac = Dune::IterativeSolverFactory<X,Y>::instance();
70 return Dune::addRegistryToFactory<TLS>(isfac, Dune::IterativeSolverTag{});
71#else
72 using OpTraits = Dune::OperatorTraits<LinearOperator>;
73 auto& pfac = Dune::PreconditionerFactory<LinearOperator>::instance();
74 Dune::addRegistryToFactory<OpTraits>(pfac, Dumux::MultiTypeBlockMatrixPreconditionerTag{});
75 auto& sfac = Dune::SolverFactory<LinearOperator>::instance();
76 return Dune::addRegistryToFactory<OpTraits>(sfac, Dumux::MultiTypeBlockMatrixSolverTag{});
77#endif
78}
79} // end namespace
80
88template<class Matrix, class LinearOperator>
90{
92 initSolverFactoriesForMultiTypeBlockMatrix<LinearOperator>();
93 else
94 Dune::initSolverFactories<LinearOperator>();
95}
96
103template <class LinearSolverTraits, class LinearAlgebraTraits>
105{
106 using Matrix = typename LinearAlgebraTraits::Matrix;
107 using Vector = typename LinearAlgebraTraits::Vector;
108 using ScalarProduct = Dune::ScalarProduct<Vector>;
109#if HAVE_MPI
110 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>, int>;
111#endif
112public:
113
119 IstlSolverFactoryBackend(const std::string& paramGroup = "")
120 : paramGroup_(paramGroup)
121 , isParallel_(Dune::MPIHelper::getCommunication().size() > 1)
122 {
123 if (isParallel_)
124 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
125
126 firstCall_ = true;
127 initializeParameters_();
128 scalarProduct_ = std::make_shared<ScalarProduct>();
129 solverCategory_ = Dune::SolverCategory::sequential;
130 }
131
139 IstlSolverFactoryBackend(const typename LinearSolverTraits::GridView& gridView,
140 const typename LinearSolverTraits::DofMapper& dofMapper,
141 const std::string& paramGroup = "")
142 : paramGroup_(paramGroup)
143#if HAVE_MPI
144 , isParallel_(gridView.comm().size() > 1)
145#endif
146 {
147 firstCall_ = true;
148 initializeParameters_();
149#if HAVE_MPI
150 solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView);
151 if (solverCategory_ != Dune::SolverCategory::sequential)
152 {
153 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
154 comm_ = std::make_shared<Comm>(gridView.comm(), solverCategory_);
155 scalarProduct_ = Dune::createScalarProduct<Vector>(*comm_, solverCategory_);
156 parallelHelper_->createParallelIndexSet(*comm_);
157 }
158 else
159 scalarProduct_ = std::make_shared<ScalarProduct>();
160#else
161 solverCategory_ = Dune::SolverCategory::sequential;
162 scalarProduct_ = std::make_shared<ScalarProduct>();
163#endif
164 }
165
173 bool solve(Matrix& A, Vector& x, Vector& b)
174 {
175#if HAVE_MPI
176 solveSequentialOrParallel_(A, x, b);
177#else
178 solveSequential_(A, x, b);
179#endif
180 firstCall_ = false;
181 return result_.converged;
182 }
183
184 const Dune::InverseOperatorResult& result() const
185 {
186 return result_;
187 }
188
189 const std::string& name() const
190 {
191 return name_;
192 }
193
194 Scalar norm(const Vector& x) const
195 {
196#if HAVE_MPI
198 {
199 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
200 {
201 auto y(x); // make a copy because the vector needs to be made consistent
202 using GV = typename LinearSolverTraits::GridView;
203 using DM = typename LinearSolverTraits::DofMapper;
204 ParallelVectorHelper<GV, DM, LinearSolverTraits::dofCodim> vectorHelper(parallelHelper_->gridView(), parallelHelper_->dofMapper());
205 vectorHelper.makeNonOverlappingConsistent(y);
206 return scalarProduct_->norm(y);
207 }
208 }
209#endif
210 return scalarProduct_->norm(x);
211 }
212
213private:
214
215 void initializeParameters_()
216 {
218 checkMandatoryParameters_();
219 name_ = params_.get<std::string>("preconditioner.type") + "-preconditioned " + params_.get<std::string>("type");
220 if (params_.get<int>("verbose", 0) > 0)
221 std::cout << "Initialized linear solver of type: " << name_ << std::endl;
222 }
223
224 void checkMandatoryParameters_()
225 {
226 if (!params_.hasKey("type"))
227 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
228
229 if (!params_.hasKey("preconditioner.type"))
230 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
231 }
232
233#if HAVE_MPI
234 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
235 {
236 // Dune::MultiTypeBlockMatrix does not provide a ColIterator which is needed by Dune::NonoverlappingSchwarzOperator.
237 // We therefore can only solve these types of systems sequentially.
238 // TODO: This can be adapted once the situation in Dune ISTL changes.
239 if constexpr (LinearSolverTraits::canCommunicate && !isMultiTypeBlockMatrix<Matrix>::value)
240 {
241 if (isParallel_)
242 {
243 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
244 {
245 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
246 solveParallel_<PTraits>(A, x, b);
247 }
248 else
249 {
250 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
251 solveParallel_<PTraits>(A, x, b);
252 }
253 }
254 else
255 solveSequential_(A, x, b);
256 }
257 else
258 {
259 solveSequential_(A, x, b);
260 }
261 }
262
263 template<class ParallelTraits>
264 void solveParallel_(Matrix& A, Vector& x, Vector& b)
265 {
266 using LinearOperator = typename ParallelTraits::LinearOperator;
267
268 if (firstCall_)
269 initSolverFactories<Matrix, LinearOperator>();
270
271 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, *parallelHelper_);
272 auto linearOperator = std::make_shared<LinearOperator>(A, *comm_);
273
274 // solve linear system
275 apply_(linearOperator, x, b);
276 }
277#endif // HAVE_MPI
278
279 void solveSequential_(Matrix& A, Vector& x, Vector& b)
280 {
281 // construct linear operator
282 using Traits = typename LinearSolverTraits::template Sequential<Matrix, Vector>;
283 using LinearOperator = typename Traits::LinearOperator;
284 auto linearOperator = std::make_shared<LinearOperator>(A);
285
286 if (firstCall_)
287 initSolverFactories<Matrix, LinearOperator>();
288
289 // solve linear system
290 apply_(linearOperator, x, b);
291 }
292
293 template<class LinearOperator>
294 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
295 {
296 try {
297 return Dune::getSolverFromFactory(fop, params_);
298 } catch(Dune::Exception& e) {
299 std::cerr << "Dune::Exception during solver construction: " << e.what() << std::endl;
300 return std::decay_t<decltype(Dune::getSolverFromFactory(fop, params_))>();
301 }
302 }
303
304 template<class LinearOperator>
305 void apply_(std::shared_ptr<LinearOperator>& fop, Vector& x, Vector& b)
306 {
307 auto solver = getSolverFromFactory_(fop);
308
309 // make solver constructor failure recoverable by throwing an exception
310 // on all processes if one or more processes fail
311 bool success = static_cast<bool>(solver);
312#if HAVE_MPI
313 int successRemote = success;
314 if (isParallel_)
315 successRemote = comm_->communicator().min(success);
316
317 if (!success)
318 DUNE_THROW(Dune::Exception, "Could not create ISTL solver");
319 else if (!successRemote)
320 DUNE_THROW(Dune::Exception, "Could not create ISTL solver on remote process");
321#else
322 if (!success)
323 DUNE_THROW(Dune::Exception, "Could not create ISTL solver");
324#endif
325
326 // solve linear system (here we assume that either all processes are successful or all fail)
327 solver->apply(x, b, result_);
328 }
329
330 const std::string paramGroup_;
331#if HAVE_MPI
332 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
333 std::shared_ptr<Comm> comm_;
334#endif
335 bool isParallel_ = false;
336 bool firstCall_;
337
338 std::shared_ptr<ScalarProduct> scalarProduct_;
339 Dune::SolverCategory::Category solverCategory_;
340 Dune::InverseOperatorResult result_;
341 Dune::ParameterTree params_;
342 std::string name_;
343};
344
345} // end namespace Dumux
346
347#endif
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:189
IstlSolverFactoryBackend(const std::string &paramGroup="")
Construct the backend for the sequential case only.
Definition: istlsolverfactorybackend.hh:119
IstlSolverFactoryBackend(const typename LinearSolverTraits::GridView &gridView, const typename LinearSolverTraits::DofMapper &dofMapper, const std::string &paramGroup="")
Construct the backend for parallel or sequential runs.
Definition: istlsolverfactorybackend.hh:139
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: istlsolverfactorybackend.hh:173
Scalar norm(const Vector &x) const
Definition: istlsolverfactorybackend.hh:194
const Dune::InverseOperatorResult & result() const
Definition: istlsolverfactorybackend.hh:184
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 &paramGroup="")
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
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.
static constexpr bool canCommunicate
Definition: gridcapabilities.hh:51
Definition: adapt.hh:17
void initSolverFactories()
Initialize the solver factories for regular matrices or MultiTypeBlockMatrices.
Definition: istlsolverfactorybackend.hh:89
Definition: common/pdesolver.hh:24
Provides a helper class for nonoverlapping decomposition.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
Solver category.
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.