version 3.11-dev
Loading...
Searching...
No Matches
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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13
14#ifndef DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
15#define DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
16
17#include <memory>
18#include <type_traits>
19
20#include <dune/common/version.hh>
21#include <dune/common/parallel/mpihelper.hh>
22#include <dune/common/parametertree.hh>
23
25
26// preconditioners
27#include <dune/istl/preconditioners.hh>
28#include <dune/istl/paamg/amg.hh>
29#include "preconditioners.hh"
30
31// solvers
32#include <dune/istl/solvers.hh>
33#include <dune/istl/solverfactory.hh>
34
42
43namespace Dumux {
44
45// initSolverFactoriesForMultiTypeBlockMatrix differs in different compilation units,
46// so we have it in an anonymous namespace
47namespace {
48
57template<class LinearOperator>
58int initSolverFactoriesForMultiTypeBlockMatrix()
59{
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{});
72#else
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{});
78#endif
79}
80} // end namespace
81
89template<class Matrix, class LinearOperator>
91{
93 initSolverFactoriesForMultiTypeBlockMatrix<LinearOperator>();
94 else
95 Dune::initSolverFactories<LinearOperator>();
96}
97
104template <class LinearSolverTraits, class LinearAlgebraTraits>
106{
107 using Matrix = typename LinearAlgebraTraits::Matrix;
108 using Vector = typename LinearAlgebraTraits::Vector;
109 using ScalarProduct = Dune::ScalarProduct<Vector>;
110#if HAVE_MPI
111 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>, int>;
112#endif
113public:
114
120 IstlSolverFactoryBackend(const std::string& paramGroup = "")
121 : paramGroup_(paramGroup)
122 , isParallel_(Dune::MPIHelper::getCommunication().size() > 1)
123 {
124 if (isParallel_)
125 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
126
127 firstCall_ = true;
128 initializeParameters_();
129 scalarProduct_ = std::make_shared<ScalarProduct>();
130 solverCategory_ = Dune::SolverCategory::sequential;
131 }
132
140 IstlSolverFactoryBackend(const typename LinearSolverTraits::GridView& gridView,
141 const typename LinearSolverTraits::DofMapper& dofMapper,
142 const std::string& paramGroup = "")
143 : paramGroup_(paramGroup)
144#if HAVE_MPI
145 , isParallel_(gridView.comm().size() > 1)
146#endif
147 {
148 firstCall_ = true;
149 initializeParameters_();
150#if HAVE_MPI
151 solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView);
152 if (solverCategory_ != Dune::SolverCategory::sequential)
153 {
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_);
158 }
159 else
160 scalarProduct_ = std::make_shared<ScalarProduct>();
161#else
162 solverCategory_ = Dune::SolverCategory::sequential;
163 scalarProduct_ = std::make_shared<ScalarProduct>();
164#endif
165 }
166
174 bool solve(Matrix& A, Vector& x, Vector& b)
175 {
176#if HAVE_MPI
177 solveSequentialOrParallel_(A, x, b);
178#else
179 solveSequential_(A, x, b);
180#endif
181 firstCall_ = false;
182 return result_.converged;
183 }
184
185 const Dune::InverseOperatorResult& result() const
186 {
187 return result_;
188 }
189
190 const std::string& name() const
191 {
192 return name_;
193 }
194
195 Scalar norm(const Vector& x) const
196 {
197#if HAVE_MPI
198 if constexpr (LinearSolverTraits::canCommunicate && !isMultiTypeBlockVector<Vector>::value)
199 {
200 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
201 {
202 auto y(x); // make a copy because the vector needs to be made consistent
203 using GV = typename LinearSolverTraits::GridView;
204 using DM = typename LinearSolverTraits::DofMapper;
205 if constexpr (requires { LinearSolverTraits::dofCodims; })
206 {
207 MultiCodimParallelVectorHelper<GV, DM> vectorHelper(parallelHelper_->gridView(), parallelHelper_->dofMapper());
208 vectorHelper.makeNonOverlappingConsistent(y, LinearSolverTraits::dofCodims);
209 }
210 else
211 {
212 ParallelVectorHelper<GV, DM, LinearSolverTraits::dofCodim> vectorHelper(parallelHelper_->gridView(), parallelHelper_->dofMapper());
213 vectorHelper.makeNonOverlappingConsistent(y);
214 }
215 return scalarProduct_->norm(y);
216 }
217 }
218#endif
219 return scalarProduct_->norm(x);
220 }
221
222private:
223
224 void initializeParameters_()
225 {
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;
231 }
232
233 void checkMandatoryParameters_()
234 {
235 if (!params_.hasKey("type"))
236 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
237
238 if (!params_.hasKey("preconditioner.type"))
239 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
240 }
241
242#if HAVE_MPI
243 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
244 {
245 // Dune::MultiTypeBlockMatrix does not provide a ColIterator which is needed by Dune::NonoverlappingSchwarzOperator.
246 // We therefore can only solve these types of systems sequentially.
247 // TODO: This can be adapted once the situation in Dune ISTL changes.
248 if constexpr (LinearSolverTraits::canCommunicate && !isMultiTypeBlockMatrix<Matrix>::value)
249 {
250 if (isParallel_)
251 {
252 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
253 {
254 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
255 solveParallel_<PTraits>(A, x, b);
256 }
257 else
258 {
259 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
260 solveParallel_<PTraits>(A, x, b);
261 }
262 }
263 else
264 solveSequential_(A, x, b);
265 }
266 else
267 {
268 solveSequential_(A, x, b);
269 }
270 }
271
272 template<class ParallelTraits>
273 void solveParallel_(Matrix& A, Vector& x, Vector& b)
274 {
275 using LinearOperator = typename ParallelTraits::LinearOperator;
276
277 if (firstCall_)
279
281 auto linearOperator = std::make_shared<LinearOperator>(A, *comm_);
282
283 // solve linear system
284 apply_(linearOperator, x, b);
285 }
286#endif // HAVE_MPI
287
288 void solveSequential_(Matrix& A, Vector& x, Vector& b)
289 {
290 // construct linear operator
291 using Traits = typename LinearSolverTraits::template Sequential<Matrix, Vector>;
292 using LinearOperator = typename Traits::LinearOperator;
293 auto linearOperator = std::make_shared<LinearOperator>(A);
294
295 if (firstCall_)
297
298 // solve linear system
299 apply_(linearOperator, x, b);
300 }
301
302 template<class LinearOperator>
303 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
304 {
305 try {
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_))>();
310 }
311 }
312
313 template<class LinearOperator>
314 void apply_(std::shared_ptr<LinearOperator>& fop, Vector& x, Vector& b)
315 {
316 auto solver = getSolverFromFactory_(fop);
317
318 // make solver constructor failure recoverable by throwing an exception
319 // on all processes if one or more processes fail
320 bool success = static_cast<bool>(solver);
321#if HAVE_MPI
322 int successRemote = success;
323 if (isParallel_)
324 successRemote = comm_->communicator().min(success);
325
326 if (!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");
330#else
331 if (!success)
332 DUNE_THROW(Dune::Exception, "Could not create ISTL solver");
333#endif
334
335 // solve linear system (here we assume that either all processes are successful or all fail)
336 solver->apply(x, b, result_);
337 }
338
339 const std::string paramGroup_;
340#if HAVE_MPI
341 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
342 std::shared_ptr<Comm> comm_;
343#endif
344 bool isParallel_ = false;
345 bool firstCall_;
346
347 std::shared_ptr<ScalarProduct> scalarProduct_;
348 Dune::SolverCategory::Category solverCategory_;
349 Dune::InverseOperatorResult result_;
350 Dune::ParameterTree params_;
351 std::string name_;
352};
353
354} // end namespace Dumux
355
356#endif
const std::string & name() const
Definition istlsolverfactorybackend.hh:190
IstlSolverFactoryBackend(const std::string &paramGroup="")
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 &paramGroup="")
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 &paramGroup="")
Construct the solver.
Definition solver.hh:43
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: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
Definition adapt.hh:17
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.
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.