3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
amgbackend.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 *****************************************************************************/
25#ifndef DUMUX_PARALLEL_AMGBACKEND_HH
26#define DUMUX_PARALLEL_AMGBACKEND_HH
27
28#include <memory>
29
30#include <dune/common/exceptions.hh>
31#include <dune/common/parallel/indexset.hh>
32
33#include <dune/common/version.hh>
34#if DUNE_VERSION_LT(DUNE_COMMON,2,7)
35#include <dune/common/parallel/mpicollectivecommunication.hh>
36#else
37#include <dune/common/parallel/mpicommunication.hh>
38#endif
39
40#include <dune/grid/common/capabilities.hh>
41#include <dune/istl/paamg/amg.hh>
42#include <dune/istl/paamg/pinfo.hh>
43#include <dune/istl/solvers.hh>
44
47
48namespace Dumux {
49
55template <class LinearSolverTraits>
57{
58public:
64 AMGBiCGSTABBackend(const std::string& paramGroup = "")
66 , isParallel_(Dune::MPIHelper::getCollectiveCommunication().size() > 1)
67 {
68 if (isParallel_)
69 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
70
71 reset();
72 }
73
81 AMGBiCGSTABBackend(const typename LinearSolverTraits::GridView& gridView,
82 const typename LinearSolverTraits::DofMapper& dofMapper,
83 const std::string& paramGroup = "")
85#if HAVE_MPI
86 , phelper_(std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper))
87 , isParallel_(Dune::MPIHelper::getCollectiveCommunication().size() > 1)
88#endif
89 {
90 reset();
91 }
92
100 template<class Matrix, class Vector>
101 bool solve(Matrix& A, Vector& x, Vector& b)
102 {
103#if HAVE_MPI
104 solveSequentialOrParallel_(A, x, b);
105#else
106 solveSequential_(A, x, b);
107#endif
108 firstCall_ = false;
109 return result_.converged;
110 }
111
113 void reset()
114 {
115 firstCall_ = true;
116 }
117
121 std::string name() const
122 {
123 return "AMG-preconditioned BiCGSTAB solver";
124 }
125
129 const Dune::InverseOperatorResult& result() const
130 {
131 return result_;
132 }
133
134private:
135
136#if HAVE_MPI
137 template<class Matrix, class Vector>
138 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
139 {
140 if constexpr (LinearSolverTraits::canCommunicate)
141 {
142 if (isParallel_)
143 {
144 if (LinearSolverTraits::isNonOverlapping(phelper_->gridView()))
145 {
146 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
147 solveParallel_<PTraits>(A, x, b);
148 }
149 else
150 {
151 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
152 solveParallel_<PTraits>(A, x, b);
153 }
154 }
155 else
156 solveSequential_(A, x, b);
157 }
158 else
159 {
160 solveSequential_(A, x, b);
161 }
162 }
163
164 template<class ParallelTraits, class Matrix, class Vector>
165 void solveParallel_(Matrix& A, Vector& x, Vector& b)
166 {
167 using Comm = typename ParallelTraits::Comm;
168 using LinearOperator = typename ParallelTraits::LinearOperator;
169 using ScalarProduct = typename ParallelTraits::ScalarProduct;
170
171 if (firstCall_)
172 phelper_->initGhostsAndOwners();
173
174 std::shared_ptr<Comm> comm;
175 std::shared_ptr<LinearOperator> linearOperator;
176 std::shared_ptr<ScalarProduct> scalarProduct;
177 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, comm, linearOperator, scalarProduct, *phelper_);
178
179 using SeqSmoother = Dune::SeqSSOR<Matrix, Vector, Vector>;
180 using Smoother = typename ParallelTraits::template Preconditioner<SeqSmoother>;
181 solveWithAmg_<Smoother>(A, x, b, linearOperator, comm, scalarProduct);
182 }
183#endif // HAVE_MPI
184
185 template<class Matrix, class Vector>
186 void solveSequential_(Matrix& A, Vector& x, Vector& b)
187 {
188 using Comm = Dune::Amg::SequentialInformation;
189 using Traits = typename LinearSolverTraits::template Sequential<Matrix, Vector>;
190 using LinearOperator = typename Traits::LinearOperator;
191 using ScalarProduct = typename Traits::ScalarProduct;
192
193 auto comm = std::make_shared<Comm>();
194 auto linearOperator = std::make_shared<LinearOperator>(A);
195 auto scalarProduct = std::make_shared<ScalarProduct>();
196
197 using Smoother = Dune::SeqSSOR<Matrix, Vector, Vector>;
198 solveWithAmg_<Smoother>(A, x, b, linearOperator, comm, scalarProduct);
199 }
200
201 template<class Smoother, class Matrix, class Vector, class LinearOperator, class Comm, class ScalarProduct>
202 void solveWithAmg_(Matrix& A, Vector& x, Vector& b,
203 std::shared_ptr<LinearOperator>& linearOperator,
204 std::shared_ptr<Comm>& comm,
205 std::shared_ptr<ScalarProduct>& scalarProduct)
206 {
207 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
208 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<Matrix, Dune::Amg::FirstDiagonal>>;
209
212 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
213 params.setDefaultValuesIsotropic(LinearSolverTraits::GridView::dimension);
214 params.setDebugLevel(this->verbosity());
215 Criterion criterion(params);
216 SmootherArgs smootherArgs;
217 smootherArgs.iterations = 1;
218 smootherArgs.relaxationFactor = 1;
219
220 using Amg = Dune::Amg::AMG<LinearOperator, Vector, Smoother, Comm>;
221 auto amg = std::make_shared<Amg>(*linearOperator, criterion, smootherArgs, *comm);
222
223 Dune::BiCGSTABSolver<Vector> solver(*linearOperator, *scalarProduct, *amg, this->residReduction(), this->maxIter(),
224 comm->communicator().rank() == 0 ? this->verbosity() : 0);
225
226 solver.apply(x, b, result_);
227 }
228
229#if HAVE_MPI
230 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> phelper_;
231#endif
232 Dune::InverseOperatorResult result_;
233 bool isParallel_ = false;
234 bool firstCall_;
235};
236
237// deprecation helper -> will be removed after 3.2
238template<class GridView, class Traits>
241
242} // namespace Dumux
243
246
247namespace Dumux {
248
255template<class TypeTag>
256using AMGBackend [[deprecated("Use AMGBiCGSTABBackend instead. Will be removed after 3.2!")]]
259
260} // namespace Dumux
261
262#endif // DUMUX_AMGBACKEND_HH
Define traits for linear solvers.
Provides a helper class for nonoverlapping decomposition.
Base class for linear solvers.
Definition: adapt.hh:29
Definition: common/pdesolver.hh:35
A linear solver based on the ISTL AMG preconditioner and the ISTL BiCGSTAB solver.
Definition: amgbackend.hh:57
const Dune::InverseOperatorResult & result() const
The result containing the convergence history.
Definition: amgbackend.hh:129
void reset()
reset the linear solver
Definition: amgbackend.hh:113
std::string name() const
The name of the solver.
Definition: amgbackend.hh:121
AMGBiCGSTABBackend(const std::string &paramGroup="")
Construct the backend for the sequential case only.
Definition: amgbackend.hh:64
AMGBiCGSTABBackend(const typename LinearSolverTraits::GridView &gridView, const typename LinearSolverTraits::DofMapper &dofMapper, const std::string &paramGroup="")
Construct the backend for parallel or sequential runs.
Definition: amgbackend.hh:81
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: amgbackend.hh:101
The implementation is specialized for the different discretizations.
Definition: linearsolvertraits.hh:68
A parallel helper class providing a nonoverlapping decomposition of all degrees of freedom.
Definition: parallelhelpers.hh:47
Base class for linear solvers.
Definition: solver.hh:37
double residReduction() const
the linear solver residual reduction
Definition: solver.hh:131
int maxIter() const
the maximum number of linear solver iterations
Definition: solver.hh:123
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:111
int verbosity() const
the verbosity level
Definition: solver.hh:115
Declares all properties used in Dumux.