3.5-git
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::getCommunication().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::getCommunication().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
167 void updateAfterGridAdaption(const typename LinearSolverTraits::GridView& gridView,
168 const typename LinearSolverTraits::DofMapper& dofMapper)
169 {
170#if HAVE_MPI
171 if (isParallel_)
172 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
173#endif
174 }
175
183 template<class Matrix, class Vector>
184 bool solve(Matrix& A, Vector& x, Vector& b)
185 {
186#if HAVE_MPI
187 solveSequentialOrParallel_(A, x, b);
188#else
189 solveSequential_(A, x, b);
190#endif
191 firstCall_ = false;
192 return result_.converged;
193 }
194
195 const Dune::InverseOperatorResult& result() const
196 {
197 return result_;
198 }
199
200 const std::string& name() const
201 {
202 return name_;
203 }
204
205private:
206
207 void initializeParameters_()
208 {
210 checkMandatoryParameters_();
211 name_ = params_.get<std::string>("preconditioner.type") + "-preconditioned " + params_.get<std::string>("type");
212 if (params_.get<int>("verbose", 0) > 0)
213 std::cout << "Initialized linear solver of type: " << name_ << std::endl;
214 }
215
216 void checkMandatoryParameters_()
217 {
218 if (!params_.hasKey("type"))
219 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
220
221 if (!params_.hasKey("preconditioner.type"))
222 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
223 }
224
225#if HAVE_MPI
226 template<class Matrix, class Vector>
227 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
228 {
229 // Dune::MultiTypeBlockMatrix does not provide a ColIterator which is needed by Dune::NonoverlappingSchwarzOperator.
230 // We therefore can only solve these types of systems sequentially.
231 // TODO: This can be adapted once the situation in Dune ISTL changes.
232 if constexpr (LinearSolverTraits::canCommunicate && !isMultiTypeBlockMatrix<Matrix>::value)
233 {
234 if (isParallel_)
235 {
236 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
237 {
238 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
239 solveParallel_<PTraits>(A, x, b);
240 }
241 else
242 {
243 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
244 solveParallel_<PTraits>(A, x, b);
245 }
246 }
247 else
248 solveSequential_(A, x, b);
249 }
250 else
251 {
252 solveSequential_(A, x, b);
253 }
254 }
255
256 template<class ParallelTraits, class Matrix, class Vector>
257 void solveParallel_(Matrix& A, Vector& x, Vector& b)
258 {
259#if DUNE_VERSION_GTE(DUNE_ISTL,2,8)
260 using Comm = typename ParallelTraits::Comm;
261 using LinearOperator = typename ParallelTraits::LinearOperator;
262 using ScalarProduct = typename ParallelTraits::ScalarProduct;
263
264 if (firstCall_)
265 initSolverFactories<Matrix, LinearOperator>();
266
267 std::shared_ptr<Comm> comm;
268 std::shared_ptr<LinearOperator> linearOperator;
269 std::shared_ptr<ScalarProduct> scalarProduct;
270 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, comm, linearOperator, scalarProduct, *parallelHelper_);
271
272 // construct solver
273 auto solver = getSolverFromFactory_(linearOperator);
274
275 // solve linear system
276 solver->apply(x, b, result_);
277#else
278 DUNE_THROW(Dune::NotImplemented, "Parallel solvers only available for dune-istl > 2.7.0");
279#endif
280 }
281#endif // HAVE_MPI
282
283 template<class Matrix, class Vector>
284 void solveSequential_(Matrix& A, Vector& x, Vector& b)
285 {
286 // construct linear operator
287 using Traits = typename LinearSolverTraits::template Sequential<Matrix, Vector>;
288 using LinearOperator = typename Traits::LinearOperator;
289 auto linearOperator = std::make_shared<LinearOperator>(A);
290
291 if (firstCall_)
292 initSolverFactories<Matrix, LinearOperator>();
293
294 // construct solver
295 auto solver = getSolverFromFactory_(linearOperator);
296
297 // solve linear system
298 solver->apply(x, b, result_);
299 }
300
301 template<class LinearOperator>
302 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
303 {
304 try { return Dune::getSolverFromFactory(fop, params_); }
305 catch(Dune::Exception& e)
306 {
307 std::cerr << "Could not create solver with factory" << std::endl;
308 std::cerr << e.what() << std::endl;
309 throw e;
310 }
311 }
312
313 const std::string paramGroup_;
314#if HAVE_MPI
315 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
316#endif
317 bool isParallel_ = false;
318 bool firstCall_;
319
320 Dune::InverseOperatorResult result_;
321 Dune::ParameterTree params_;
322 std::string name_;
323};
324
325} // end namespace Dumux
326
327#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 ...
Dumux preconditioners for iterative solvers.
Provides a helper class for nonoverlapping decomposition.
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Base class for linear solvers.
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:200
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:184
const Dune::InverseOperatorResult & result() const
Definition: istlsolverfactorybackend.hh:195
void updateAfterGridAdaption(const typename LinearSolverTraits::GridView &gridView, const typename LinearSolverTraits::DofMapper &dofMapper)
Update the solver after grid adaption.
Definition: istlsolverfactorybackend.hh:167
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:49
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