version 3.8
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/parallel/mpihelper.hh>
20#include <dune/common/parametertree.hh>
21
23
24// preconditioners
25#include <dune/istl/preconditioners.hh>
26#include <dune/istl/paamg/amg.hh>
27#include "preconditioners.hh"
28
29// solvers
30#include <dune/istl/solvers.hh>
31#include <dune/istl/solverfactory.hh>
32
39
40namespace Dumux {
41
42// initSolverFactoriesForMultiTypeBlockMatrix differs in different compilation units,
43// so we have it in an anonymous namespace
44namespace {
45
54template<class LinearOperator>
55int initSolverFactoriesForMultiTypeBlockMatrix()
56{
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{});
68}
69} // end namespace
70
78template<class Matrix, class LinearOperator>
80{
82 initSolverFactoriesForMultiTypeBlockMatrix<LinearOperator>();
83 else
84 Dune::initSolverFactories<LinearOperator>();
85}
86
93namespace Detail {
94
95template <class LinearSolverTraits>
97{
98public:
99
106 void updateAfterGridAdaption(const typename LinearSolverTraits::GridView& gridView,
107 const typename LinearSolverTraits::DofMapper& dofMapper)
108 {
109#if HAVE_MPI
110 if (isParallel_)
111 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
112#endif
113 }
114
122 template<class Matrix, class Vector>
123 bool solve(Matrix& A, Vector& x, Vector& b)
124 {
125#if HAVE_MPI
126 solveSequentialOrParallel_(A, x, b);
127#else
128 solveSequential_(A, x, b);
129#endif
130 firstCall_ = false;
131 return result_.converged;
132 }
133
134 const Dune::InverseOperatorResult& result() const
135 {
136 return result_;
137 }
138
139 const std::string& name() const
140 {
141 return name_;
142 }
143
144private:
145
146 void initializeParameters_()
147 {
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;
153 }
154
155 void checkMandatoryParameters_()
156 {
157 if (!params_.hasKey("type"))
158 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
159
160 if (!params_.hasKey("preconditioner.type"))
161 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
162 }
163
164#if HAVE_MPI
165 template<class Matrix, class Vector>
166 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
167 {
168 // Dune::MultiTypeBlockMatrix does not provide a ColIterator which is needed by Dune::NonoverlappingSchwarzOperator.
169 // We therefore can only solve these types of systems sequentially.
170 // TODO: This can be adapted once the situation in Dune ISTL changes.
171 if constexpr (LinearSolverTraits::canCommunicate && !isMultiTypeBlockMatrix<Matrix>::value)
172 {
173 if (isParallel_)
174 {
175 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
176 {
177 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
178 solveParallel_<PTraits>(A, x, b);
179 }
180 else
181 {
182 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
183 solveParallel_<PTraits>(A, x, b);
184 }
185 }
186 else
187 solveSequential_(A, x, b);
188 }
189 else
190 {
191 solveSequential_(A, x, b);
192 }
193 }
194
195 template<class ParallelTraits, class Matrix, class Vector>
196 void solveParallel_(Matrix& A, Vector& x, Vector& b)
197 {
198 using Comm = typename ParallelTraits::Comm;
199 using LinearOperator = typename ParallelTraits::LinearOperator;
200 using ScalarProduct = typename ParallelTraits::ScalarProduct;
201
202 if (firstCall_)
203 initSolverFactories<Matrix, LinearOperator>();
204
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_);
209
210 // construct solver
211 auto solver = getSolverFromFactory_(linearOperator);
212
213 // solve linear system
214 solver->apply(x, b, result_);
215 }
216#endif // HAVE_MPI
217
218 template<class Matrix, class Vector>
219 void solveSequential_(Matrix& A, Vector& x, Vector& b)
220 {
221 // construct linear operator
222 using Traits = typename LinearSolverTraits::template Sequential<Matrix, Vector>;
223 using LinearOperator = typename Traits::LinearOperator;
224 auto linearOperator = std::make_shared<LinearOperator>(A);
225
226 if (firstCall_)
227 initSolverFactories<Matrix, LinearOperator>();
228
229 // construct solver
230 auto solver = getSolverFromFactory_(linearOperator);
231
232 // solve linear system
233 solver->apply(x, b, result_);
234 }
235
236 template<class LinearOperator>
237 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
238 {
239 try { return Dune::getSolverFromFactory(fop, params_); }
240 catch(Dune::Exception& e)
241 {
242 std::cerr << "Could not create solver with factory" << std::endl;
243 std::cerr << e.what() << std::endl;
244 throw e;
245 }
246 }
247
248 const std::string paramGroup_;
249#if HAVE_MPI
250 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
251#endif
252 bool isParallel_ = false;
253 bool firstCall_;
254
255 Dune::InverseOperatorResult result_;
256 Dune::ParameterTree params_;
257 std::string name_;
258};
259
260template <class LinearSolverTraits, class LinearAlgebraTraits>
262{
263 using Matrix = typename LinearAlgebraTraits::Matrix;
264 using Vector = typename LinearAlgebraTraits::Vector;
265 using ScalarProduct = Dune::ScalarProduct<Vector>;
266#if HAVE_MPI
267 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>, int>;
268#endif
269public:
270
276 IstlSolverFactoryBackend(const std::string& paramGroup = "")
277 : paramGroup_(paramGroup)
278 , isParallel_(Dune::MPIHelper::getCommunication().size() > 1)
279 {
280 if (isParallel_)
281 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
282
283 firstCall_ = true;
284 initializeParameters_();
285 scalarProduct_ = std::make_shared<ScalarProduct>();
286 solverCategory_ = Dune::SolverCategory::sequential;
287 }
288
296 IstlSolverFactoryBackend(const typename LinearSolverTraits::GridView& gridView,
297 const typename LinearSolverTraits::DofMapper& dofMapper,
298 const std::string& paramGroup = "")
299 : paramGroup_(paramGroup)
300#if HAVE_MPI
301 , isParallel_(gridView.comm().size() > 1)
302#endif
303 {
304 firstCall_ = true;
305 initializeParameters_();
306#if HAVE_MPI
307 solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView);
308 if (solverCategory_ != Dune::SolverCategory::sequential)
309 {
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_);
314 }
315 else
316 scalarProduct_ = std::make_shared<ScalarProduct>();
317#else
318 solverCategory_ = Dune::SolverCategory::sequential;
319 scalarProduct_ = std::make_shared<ScalarProduct>();
320#endif
321 }
322
330 bool solve(Matrix& A, Vector& x, Vector& b)
331 {
332#if HAVE_MPI
333 solveSequentialOrParallel_(A, x, b);
334#else
335 solveSequential_(A, x, b);
336#endif
337 firstCall_ = false;
338 return result_.converged;
339 }
340
341 const Dune::InverseOperatorResult& result() const
342 {
343 return result_;
344 }
345
346 const std::string& name() const
347 {
348 return name_;
349 }
350
351 Scalar norm(const Vector& x) const
352 {
353#if HAVE_MPI
355 {
356 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
357 {
358 auto y(x); // make a copy because the vector needs to be made consistent
359 using GV = typename LinearSolverTraits::GridView;
360 using DM = typename LinearSolverTraits::DofMapper;
361 ParallelVectorHelper<GV, DM, LinearSolverTraits::dofCodim> vectorHelper(parallelHelper_->gridView(), parallelHelper_->dofMapper());
362 vectorHelper.makeNonOverlappingConsistent(y);
363 return scalarProduct_->norm(y);
364 }
365 }
366#endif
367 return scalarProduct_->norm(x);
368 }
369
370private:
371
372 void initializeParameters_()
373 {
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;
379 }
380
381 void checkMandatoryParameters_()
382 {
383 if (!params_.hasKey("type"))
384 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
385
386 if (!params_.hasKey("preconditioner.type"))
387 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
388 }
389
390#if HAVE_MPI
391 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
392 {
393 // Dune::MultiTypeBlockMatrix does not provide a ColIterator which is needed by Dune::NonoverlappingSchwarzOperator.
394 // We therefore can only solve these types of systems sequentially.
395 // TODO: This can be adapted once the situation in Dune ISTL changes.
396 if constexpr (LinearSolverTraits::canCommunicate && !isMultiTypeBlockMatrix<Matrix>::value)
397 {
398 if (isParallel_)
399 {
400 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
401 {
402 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
403 solveParallel_<PTraits>(A, x, b);
404 }
405 else
406 {
407 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
408 solveParallel_<PTraits>(A, x, b);
409 }
410 }
411 else
412 solveSequential_(A, x, b);
413 }
414 else
415 {
416 solveSequential_(A, x, b);
417 }
418 }
419
420 template<class ParallelTraits>
421 void solveParallel_(Matrix& A, Vector& x, Vector& b)
422 {
423 using LinearOperator = typename ParallelTraits::LinearOperator;
424
425 if (firstCall_)
426 initSolverFactories<Matrix, LinearOperator>();
427
428 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, *parallelHelper_);
429 auto linearOperator = std::make_shared<LinearOperator>(A, *comm_);
430
431 // construct solver
432 auto solver = getSolverFromFactory_(linearOperator);
433
434 // solve linear system
435 solver->apply(x, b, result_);
436 }
437#endif // HAVE_MPI
438
439 void solveSequential_(Matrix& A, Vector& x, Vector& b)
440 {
441 // construct linear operator
442 using Traits = typename LinearSolverTraits::template Sequential<Matrix, Vector>;
443 using LinearOperator = typename Traits::LinearOperator;
444 auto linearOperator = std::make_shared<LinearOperator>(A);
445
446 if (firstCall_)
447 initSolverFactories<Matrix, LinearOperator>();
448
449 // construct solver
450 auto solver = getSolverFromFactory_(linearOperator);
451
452 // solve linear system
453 solver->apply(x, b, result_);
454 }
455
456 template<class LinearOperator>
457 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
458 {
459 try { return Dune::getSolverFromFactory(fop, params_); }
460 catch(Dune::Exception& e)
461 {
462 std::cerr << "Could not create solver with factory" << std::endl;
463 std::cerr << e.what() << std::endl;
464 throw e;
465 }
466 }
467
468 const std::string paramGroup_;
469#if HAVE_MPI
470 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
471 std::shared_ptr<Comm> comm_;
472#endif
473 bool isParallel_ = false;
474 bool firstCall_;
475
476 std::shared_ptr<ScalarProduct> scalarProduct_;
477 Dune::SolverCategory::Category solverCategory_;
478 Dune::InverseOperatorResult result_;
479 Dune::ParameterTree params_;
480 std::string name_;
481};
482
483template<int i>
485{
486 template<class ...T>
488};
489
490template<>
492{
493 template<class T>
495};
496
497} // end namespace Detail
498
506template<class ...T>
507using IstlSolverFactoryBackend = typename Detail::IstlSolverFactoryBackendImplHelper<sizeof...(T)>::template type<T...>;
508
509} // end namespace Dumux
510
511#endif
Definition: istlsolverfactorybackend.hh:262
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: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 &paramGroup="")
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 &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
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
Definition: adapt.hh:17
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.
Solver category.
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.