3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
27#define DUMUX_LINEAR_ISTL_SOLVERFACTORYBACKEND_HH
28
29#include <memory>
30
31#include <dune/common/version.hh>
32#include <dune/common/parallel/mpihelper.hh>
33#include <dune/common/parametertree.hh>
34
36
37// preconditioners
38#include <dune/istl/preconditioners.hh>
39#include <dune/istl/paamg/amg.hh>
40#include "preconditioners.hh"
41
42// solvers
43#include <dune/istl/solvers.hh>
44#include <dune/istl/solverfactory.hh>
45
50
51namespace Dumux {
52
53// initSolverFactoriesForMultiTypeBlockMatrix differs in different compilation units,
54// so we have it in an anonymous namespace
55namespace {
56
65template<class LinearOperator>
66int initSolverFactoriesForMultiTypeBlockMatrix()
67{
68 using M = typename LinearOperator::matrix_type;
69 using X = typename LinearOperator::range_type;
70 using Y = typename LinearOperator::domain_type;
71 using TL = Dune::TypeList<M,X,Y>;
72 auto& dsfac = Dune::DirectSolverFactory<M,X,Y>::instance();
73 Dune::addRegistryToFactory<TL>(dsfac, Dumux::MultiTypeBlockMatrixDirectSolverTag{});
74#if DUNE_VERSION_GTE(DUNE_ISTL,2,8)
75 auto& pfac = Dune::PreconditionerFactory<LinearOperator,X,Y>::instance();
76#else
77 auto& pfac = Dune::PreconditionerFactory<M,X,Y>::instance();
78#endif
79 Dune::addRegistryToFactory<TL>(pfac, Dumux::MultiTypeBlockMatrixPreconditionerTag{});
80 using TLS = Dune::TypeList<X,Y>;
81 auto& isfac = Dune::IterativeSolverFactory<X,Y>::instance();
82 return Dune::addRegistryToFactory<TLS>(isfac, Dune::IterativeSolverTag{});
83}
84} // end namespace
85
93template<class Matrix, class LinearOperator>
95{
97 initSolverFactoriesForMultiTypeBlockMatrix<LinearOperator>();
98 else
99#if DUNE_VERSION_GTE(DUNE_ISTL,2,8)
100 Dune::initSolverFactories<LinearOperator>();
101#else
102 {
103 using X = typename LinearOperator::range_type;
104 using Y = typename LinearOperator::domain_type;
105 Dune::initSolverFactories<Matrix, X, Y>();
106 }
107#endif
108}
109
117template <class LinearSolverTraits>
119{
120public:
121
127 IstlSolverFactoryBackend(const std::string& paramGroup = "")
128 : paramGroup_(paramGroup)
129 , isParallel_(Dune::MPIHelper::getCollectiveCommunication().size() > 1)
130 {
131 if (isParallel_)
132 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
133
134 firstCall_ = true;
135 initializeParameters_();
136 }
137
145 IstlSolverFactoryBackend(const typename LinearSolverTraits::GridView& gridView,
146 const typename LinearSolverTraits::DofMapper& dofMapper,
147 const std::string& paramGroup = "")
148 : paramGroup_(paramGroup)
149#if HAVE_MPI
150 , isParallel_(Dune::MPIHelper::getCollectiveCommunication().size() > 1)
151#endif
152 {
153 firstCall_ = true;
154 initializeParameters_();
155#if HAVE_MPI
156 if (isParallel_)
157 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
158#endif
159 }
160
168 template<class Matrix, class Vector>
169 bool solve(Matrix& A, Vector& x, Vector& b)
170 {
171#if HAVE_MPI
172 solveSequentialOrParallel_(A, x, b);
173#else
174 solveSequential_(A, x, b);
175#endif
176 firstCall_ = false;
177 return result_.converged;
178 }
179
180 const Dune::InverseOperatorResult& result() const
181 {
182 return result_;
183 }
184
185 const std::string& name() const
186 {
187 return name_;
188 }
189
190private:
191
192 void initializeParameters_()
193 {
195 checkMandatoryParameters_();
196 name_ = params_.get<std::string>("preconditioner.type") + "-preconditioned " + params_.get<std::string>("type");
197 if (params_.get<int>("verbose", 0) > 0)
198 std::cout << "Initialized linear solver of type: " << name_ << std::endl;
199 }
200
201 void checkMandatoryParameters_()
202 {
203 if (!params_.hasKey("type"))
204 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
205
206 if (!params_.hasKey("preconditioner.type"))
207 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
208 }
209
210#if HAVE_MPI
211 template<class Matrix, class Vector>
212 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
213 {
214 // Dune::MultiTypeBlockMatrix does not provide a ColIterator which is needed by Dune::NonoverlappingSchwarzOperator.
215 // We therefore can only solve these types of systems sequentially.
216 // TODO: This can be adapted once the situation in Dune ISTL changes.
217 if constexpr (LinearSolverTraits::canCommunicate && !isMultiTypeBlockMatrix<Matrix>::value)
218 {
219 if (isParallel_)
220 {
221 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
222 {
223 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
224 solveParallel_<PTraits>(A, x, b);
225 }
226 else
227 {
228 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
229 solveParallel_<PTraits>(A, x, b);
230 }
231 }
232 else
233 solveSequential_(A, x, b);
234 }
235 else
236 {
237 solveSequential_(A, x, b);
238 }
239 }
240
241 template<class ParallelTraits, class Matrix, class Vector>
242 void solveParallel_(Matrix& A, Vector& x, Vector& b)
243 {
244#if DUNE_VERSION_GTE(DUNE_ISTL,2,8)
245 using Comm = typename ParallelTraits::Comm;
246 using LinearOperator = typename ParallelTraits::LinearOperator;
247 using ScalarProduct = typename ParallelTraits::ScalarProduct;
248
249 if (firstCall_)
250 initSolverFactories<Matrix, LinearOperator>();
251
252 std::shared_ptr<Comm> comm;
253 std::shared_ptr<LinearOperator> linearOperator;
254 std::shared_ptr<ScalarProduct> scalarProduct;
255 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, comm, linearOperator, scalarProduct, *parallelHelper_);
256
257 // construct solver
258 auto solver = getSolverFromFactory_(linearOperator);
259
260 // solve linear system
261 solver->apply(x, b, result_);
262#else
263 DUNE_THROW(Dune::NotImplemented, "Parallel solvers only available for dune-istl > 2.7.0");
264#endif
265 }
266#endif // HAVE_MPI
267
268 template<class Matrix, class Vector>
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 // construct solver
280 auto solver = getSolverFromFactory_(linearOperator);
281
282 // solve linear system
283 solver->apply(x, b, result_);
284 }
285
286 template<class LinearOperator>
287 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
288 {
289 try { return Dune::getSolverFromFactory(fop, params_); }
290 catch(Dune::Exception& e)
291 {
292 std::cerr << "Could not create solver with factory" << std::endl;
293 std::cerr << e.what() << std::endl;
294 throw e;
295 }
296 }
297
298 const std::string paramGroup_;
299#if HAVE_MPI
300 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
301#endif
302 bool isParallel_ = false;
303 bool firstCall_;
304
305 Dune::InverseOperatorResult result_;
306 Dune::ParameterTree params_;
307 std::string name_;
308};
309
310} // end namespace Dumux
311
312#endif
Type traits to be used with matrix types.
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.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
Provides a helper class for nonoverlapping decomposition.
Definition: adapt.hh:29
void initSolverFactories()
Initialize the solver factories for regular matrices or MultiTypeBlockMatrices.
Definition: istlsolverfactorybackend.hh:94
static constexpr bool canCommunicate
Definition: gridcapabilities.hh:63
Definition: common/pdesolver.hh:36
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:49
A linear solver using the dune-istl solver factory to choose the solver and preconditioner at runtime...
Definition: istlsolverfactorybackend.hh:119
const std::string & name() const
Definition: istlsolverfactorybackend.hh:185
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:145
IstlSolverFactoryBackend(const std::string &paramGroup="")
Construct the backend for the sequential case only.
Definition: istlsolverfactorybackend.hh:127
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: istlsolverfactorybackend.hh:169
const Dune::InverseOperatorResult & result() const
Definition: istlsolverfactorybackend.hh:180
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:46
Base class for linear solvers.
Definition: solver.hh:37
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:79