version 3.9
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
93template <class LinearSolverTraits, class LinearAlgebraTraits>
95{
96 using Matrix = typename LinearAlgebraTraits::Matrix;
97 using Vector = typename LinearAlgebraTraits::Vector;
98 using ScalarProduct = Dune::ScalarProduct<Vector>;
99#if HAVE_MPI
100 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>, int>;
101#endif
102public:
103
109 IstlSolverFactoryBackend(const std::string& paramGroup = "")
110 : paramGroup_(paramGroup)
111 , isParallel_(Dune::MPIHelper::getCommunication().size() > 1)
112 {
113 if (isParallel_)
114 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
115
116 firstCall_ = true;
117 initializeParameters_();
118 scalarProduct_ = std::make_shared<ScalarProduct>();
119 solverCategory_ = Dune::SolverCategory::sequential;
120 }
121
129 IstlSolverFactoryBackend(const typename LinearSolverTraits::GridView& gridView,
130 const typename LinearSolverTraits::DofMapper& dofMapper,
131 const std::string& paramGroup = "")
132 : paramGroup_(paramGroup)
133#if HAVE_MPI
134 , isParallel_(gridView.comm().size() > 1)
135#endif
136 {
137 firstCall_ = true;
138 initializeParameters_();
139#if HAVE_MPI
140 solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView);
141 if (solverCategory_ != Dune::SolverCategory::sequential)
142 {
143 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
144 comm_ = std::make_shared<Comm>(gridView.comm(), solverCategory_);
145 scalarProduct_ = Dune::createScalarProduct<Vector>(*comm_, solverCategory_);
146 parallelHelper_->createParallelIndexSet(*comm_);
147 }
148 else
149 scalarProduct_ = std::make_shared<ScalarProduct>();
150#else
151 solverCategory_ = Dune::SolverCategory::sequential;
152 scalarProduct_ = std::make_shared<ScalarProduct>();
153#endif
154 }
155
163 bool solve(Matrix& A, Vector& x, Vector& b)
164 {
165#if HAVE_MPI
166 solveSequentialOrParallel_(A, x, b);
167#else
168 solveSequential_(A, x, b);
169#endif
170 firstCall_ = false;
171 return result_.converged;
172 }
173
174 const Dune::InverseOperatorResult& result() const
175 {
176 return result_;
177 }
178
179 const std::string& name() const
180 {
181 return name_;
182 }
183
184 Scalar norm(const Vector& x) const
185 {
186#if HAVE_MPI
188 {
189 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
190 {
191 auto y(x); // make a copy because the vector needs to be made consistent
192 using GV = typename LinearSolverTraits::GridView;
193 using DM = typename LinearSolverTraits::DofMapper;
194 ParallelVectorHelper<GV, DM, LinearSolverTraits::dofCodim> vectorHelper(parallelHelper_->gridView(), parallelHelper_->dofMapper());
195 vectorHelper.makeNonOverlappingConsistent(y);
196 return scalarProduct_->norm(y);
197 }
198 }
199#endif
200 return scalarProduct_->norm(x);
201 }
202
203private:
204
205 void initializeParameters_()
206 {
208 checkMandatoryParameters_();
209 name_ = params_.get<std::string>("preconditioner.type") + "-preconditioned " + params_.get<std::string>("type");
210 if (params_.get<int>("verbose", 0) > 0)
211 std::cout << "Initialized linear solver of type: " << name_ << std::endl;
212 }
213
214 void checkMandatoryParameters_()
215 {
216 if (!params_.hasKey("type"))
217 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
218
219 if (!params_.hasKey("preconditioner.type"))
220 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
221 }
222
223#if HAVE_MPI
224 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
225 {
226 // Dune::MultiTypeBlockMatrix does not provide a ColIterator which is needed by Dune::NonoverlappingSchwarzOperator.
227 // We therefore can only solve these types of systems sequentially.
228 // TODO: This can be adapted once the situation in Dune ISTL changes.
229 if constexpr (LinearSolverTraits::canCommunicate && !isMultiTypeBlockMatrix<Matrix>::value)
230 {
231 if (isParallel_)
232 {
233 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
234 {
235 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
236 solveParallel_<PTraits>(A, x, b);
237 }
238 else
239 {
240 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
241 solveParallel_<PTraits>(A, x, b);
242 }
243 }
244 else
245 solveSequential_(A, x, b);
246 }
247 else
248 {
249 solveSequential_(A, x, b);
250 }
251 }
252
253 template<class ParallelTraits>
254 void solveParallel_(Matrix& A, Vector& x, Vector& b)
255 {
256 using LinearOperator = typename ParallelTraits::LinearOperator;
257
258 if (firstCall_)
259 initSolverFactories<Matrix, LinearOperator>();
260
261 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, *parallelHelper_);
262 auto linearOperator = std::make_shared<LinearOperator>(A, *comm_);
263
264 // solve linear system
265 apply_(linearOperator, x, b);
266 }
267#endif // HAVE_MPI
268
269 void solveSequential_(Matrix& A, Vector& x, Vector& b)
270 {
271 // construct linear operator
272 using Traits = typename LinearSolverTraits::template Sequential<Matrix, Vector>;
273 using LinearOperator = typename Traits::LinearOperator;
274 auto linearOperator = std::make_shared<LinearOperator>(A);
275
276 if (firstCall_)
277 initSolverFactories<Matrix, LinearOperator>();
278
279 // solve linear system
280 apply_(linearOperator, x, b);
281 }
282
283 template<class LinearOperator>
284 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
285 {
286 try {
287 return Dune::getSolverFromFactory(fop, params_);
288 } catch(Dune::Exception& e) {
289 std::cerr << "Dune::Exception during solver construction: " << e.what() << std::endl;
290 return std::decay_t<decltype(Dune::getSolverFromFactory(fop, params_))>();
291 }
292 }
293
294 template<class LinearOperator>
295 void apply_(std::shared_ptr<LinearOperator>& fop, Vector& x, Vector& b)
296 {
297 auto solver = getSolverFromFactory_(fop);
298
299 // make solver constructor failure recoverable by throwing an exception
300 // on all processes if one or more processes fail
301 bool success = static_cast<bool>(solver);
302#if HAVE_MPI
303 int successRemote = success;
304 if (isParallel_)
305 successRemote = comm_->communicator().min(success);
306
307 if (!success)
308 DUNE_THROW(Dune::Exception, "Could not create ISTL solver");
309 else if (!successRemote)
310 DUNE_THROW(Dune::Exception, "Could not create ISTL solver on remote process");
311#else
312 if (!success)
313 DUNE_THROW(Dune::Exception, "Could not create ISTL solver");
314#endif
315
316 // solve linear system (here we assume that either all processes are successful or all fail)
317 solver->apply(x, b, result_);
318 }
319
320 const std::string paramGroup_;
321#if HAVE_MPI
322 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
323 std::shared_ptr<Comm> comm_;
324#endif
325 bool isParallel_ = false;
326 bool firstCall_;
327
328 std::shared_ptr<ScalarProduct> scalarProduct_;
329 Dune::SolverCategory::Category solverCategory_;
330 Dune::InverseOperatorResult result_;
331 Dune::ParameterTree params_;
332 std::string name_;
333};
334
335} // end namespace Dumux
336
337#endif
A linear solver using the dune-istl solver factory to choose the solver and preconditioner at runtime...
Definition: istlsolverfactorybackend.hh:95
const std::string & name() const
Definition: istlsolverfactorybackend.hh:179
IstlSolverFactoryBackend(const std::string &paramGroup="")
Construct the backend for the sequential case only.
Definition: istlsolverfactorybackend.hh:109
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:129
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: istlsolverfactorybackend.hh:163
Scalar norm(const Vector &x) const
Definition: istlsolverfactorybackend.hh:184
const Dune::InverseOperatorResult & result() const
Definition: istlsolverfactorybackend.hh:174
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 ...
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.
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.