3.6-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/parallel/mpihelper.hh>
32#include <dune/common/parametertree.hh>
33
35
36// preconditioners
37#include <dune/istl/preconditioners.hh>
38#include <dune/istl/paamg/amg.hh>
39#include "preconditioners.hh"
40
41// solvers
42#include <dune/istl/solvers.hh>
43#include <dune/istl/solverfactory.hh>
44
49
50namespace Dumux {
51
52// initSolverFactoriesForMultiTypeBlockMatrix differs in different compilation units,
53// so we have it in an anonymous namespace
54namespace {
55
64template<class LinearOperator>
65int initSolverFactoriesForMultiTypeBlockMatrix()
66{
67 using M = typename LinearOperator::matrix_type;
68 using X = typename LinearOperator::range_type;
69 using Y = typename LinearOperator::domain_type;
70 using TL = Dune::TypeList<M,X,Y>;
71 auto& dsfac = Dune::DirectSolverFactory<M,X,Y>::instance();
72 Dune::addRegistryToFactory<TL>(dsfac, Dumux::MultiTypeBlockMatrixDirectSolverTag{});
73 auto& pfac = Dune::PreconditionerFactory<LinearOperator,X,Y>::instance();
74 Dune::addRegistryToFactory<TL>(pfac, Dumux::MultiTypeBlockMatrixPreconditionerTag{});
75 using TLS = Dune::TypeList<X,Y>;
76 auto& isfac = Dune::IterativeSolverFactory<X,Y>::instance();
77 return Dune::addRegistryToFactory<TLS>(isfac, Dune::IterativeSolverTag{});
78}
79} // end namespace
80
88template<class Matrix, class LinearOperator>
90{
92 initSolverFactoriesForMultiTypeBlockMatrix<LinearOperator>();
93 else
94 Dune::initSolverFactories<LinearOperator>();
95}
96
103template <class LinearSolverTraits>
105{
106public:
107
113 IstlSolverFactoryBackend(const std::string& paramGroup = "")
114 : paramGroup_(paramGroup)
115 , isParallel_(Dune::MPIHelper::getCommunication().size() > 1)
116 {
117 if (isParallel_)
118 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
119
120 firstCall_ = true;
121 initializeParameters_();
122 }
123
131 IstlSolverFactoryBackend(const typename LinearSolverTraits::GridView& gridView,
132 const typename LinearSolverTraits::DofMapper& dofMapper,
133 const std::string& paramGroup = "")
134 : paramGroup_(paramGroup)
135#if HAVE_MPI
136 , isParallel_(Dune::MPIHelper::getCommunication().size() > 1)
137#endif
138 {
139 firstCall_ = true;
140 initializeParameters_();
141#if HAVE_MPI
142 if (isParallel_)
143 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
144#endif
145 }
146
153 void updateAfterGridAdaption(const typename LinearSolverTraits::GridView& gridView,
154 const typename LinearSolverTraits::DofMapper& dofMapper)
155 {
156#if HAVE_MPI
157 if (isParallel_)
158 parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
159#endif
160 }
161
169 template<class Matrix, class Vector>
170 bool solve(Matrix& A, Vector& x, Vector& b)
171 {
172#if HAVE_MPI
173 solveSequentialOrParallel_(A, x, b);
174#else
175 solveSequential_(A, x, b);
176#endif
177 firstCall_ = false;
178 return result_.converged;
179 }
180
181 const Dune::InverseOperatorResult& result() const
182 {
183 return result_;
184 }
185
186 const std::string& name() const
187 {
188 return name_;
189 }
190
191private:
192
193 void initializeParameters_()
194 {
196 checkMandatoryParameters_();
197 name_ = params_.get<std::string>("preconditioner.type") + "-preconditioned " + params_.get<std::string>("type");
198 if (params_.get<int>("verbose", 0) > 0)
199 std::cout << "Initialized linear solver of type: " << name_ << std::endl;
200 }
201
202 void checkMandatoryParameters_()
203 {
204 if (!params_.hasKey("type"))
205 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Type\" parameter to select the solver");
206
207 if (!params_.hasKey("preconditioner.type"))
208 DUNE_THROW(Dune::InvalidStateException, "Solver factory needs \"LinearSolver.Preconditioner.Type\" parameter to select the preconditioner");
209 }
210
211#if HAVE_MPI
212 template<class Matrix, class Vector>
213 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
214 {
215 // Dune::MultiTypeBlockMatrix does not provide a ColIterator which is needed by Dune::NonoverlappingSchwarzOperator.
216 // We therefore can only solve these types of systems sequentially.
217 // TODO: This can be adapted once the situation in Dune ISTL changes.
218 if constexpr (LinearSolverTraits::canCommunicate && !isMultiTypeBlockMatrix<Matrix>::value)
219 {
220 if (isParallel_)
221 {
222 if (LinearSolverTraits::isNonOverlapping(parallelHelper_->gridView()))
223 {
224 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
225 solveParallel_<PTraits>(A, x, b);
226 }
227 else
228 {
229 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
230 solveParallel_<PTraits>(A, x, b);
231 }
232 }
233 else
234 solveSequential_(A, x, b);
235 }
236 else
237 {
238 solveSequential_(A, x, b);
239 }
240 }
241
242 template<class ParallelTraits, class Matrix, class Vector>
243 void solveParallel_(Matrix& A, Vector& x, Vector& b)
244 {
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 }
263#endif // HAVE_MPI
264
265 template<class Matrix, class Vector>
266 void solveSequential_(Matrix& A, Vector& x, Vector& b)
267 {
268 // construct linear operator
269 using Traits = typename LinearSolverTraits::template Sequential<Matrix, Vector>;
270 using LinearOperator = typename Traits::LinearOperator;
271 auto linearOperator = std::make_shared<LinearOperator>(A);
272
273 if (firstCall_)
274 initSolverFactories<Matrix, LinearOperator>();
275
276 // construct solver
277 auto solver = getSolverFromFactory_(linearOperator);
278
279 // solve linear system
280 solver->apply(x, b, result_);
281 }
282
283 template<class LinearOperator>
284 auto getSolverFromFactory_(std::shared_ptr<LinearOperator>& fop)
285 {
286 try { return Dune::getSolverFromFactory(fop, params_); }
287 catch(Dune::Exception& e)
288 {
289 std::cerr << "Could not create solver with factory" << std::endl;
290 std::cerr << e.what() << std::endl;
291 throw e;
292 }
293 }
294
295 const std::string paramGroup_;
296#if HAVE_MPI
297 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> parallelHelper_;
298#endif
299 bool isParallel_ = false;
300 bool firstCall_;
301
302 Dune::InverseOperatorResult result_;
303 Dune::ParameterTree params_;
304 std::string name_;
305};
306
307} // end namespace Dumux
308
309#endif
Type traits to be used with matrix types.
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
The specialized Dumux macro and tag for the ISTL registry to choose the solver and preconditioner at ...
Provides a helper class for nonoverlapping decomposition.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
void initSolverFactories()
Initialize the solver factories for regular matrices or MultiTypeBlockMatrices.
Definition: istlsolverfactorybackend.hh:89
static constexpr bool canCommunicate
Definition: gridcapabilities.hh:63
Definition: deprecated.hh:149
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:105
const std::string & name() const
Definition: istlsolverfactorybackend.hh:186
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:131
IstlSolverFactoryBackend(const std::string &paramGroup="")
Construct the backend for the sequential case only.
Definition: istlsolverfactorybackend.hh:113
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: istlsolverfactorybackend.hh:170
const Dune::InverseOperatorResult & result() const
Definition: istlsolverfactorybackend.hh:181
void updateAfterGridAdaption(const typename LinearSolverTraits::GridView &gridView, const typename LinearSolverTraits::DofMapper &dofMapper)
Update the solver after grid adaption.
Definition: istlsolverfactorybackend.hh:153
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