3.3.0
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#include <dune/common/parallel/mpicommunication.hh>
33#include <dune/grid/common/capabilities.hh>
34#include <dune/istl/paamg/amg.hh>
35#include <dune/istl/paamg/pinfo.hh>
36#include <dune/istl/solvers.hh>
37
40
41namespace Dumux {
42
48template <class LinearSolverTraits>
50{
51public:
57 AMGBiCGSTABBackend(const std::string& paramGroup = "")
59 , isParallel_(Dune::MPIHelper::getCollectiveCommunication().size() > 1)
60 {
61 if (isParallel_)
62 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
63 }
64
72 AMGBiCGSTABBackend(const typename LinearSolverTraits::GridView& gridView,
73 const typename LinearSolverTraits::DofMapper& dofMapper,
74 const std::string& paramGroup = "")
76#if HAVE_MPI
77 , isParallel_(Dune::MPIHelper::getCollectiveCommunication().size() > 1)
78#endif
79 {
80#if HAVE_MPI
81 if (isParallel_)
82 phelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
83#endif
84 }
85
93 template<class Matrix, class Vector>
94 bool solve(Matrix& A, Vector& x, Vector& b)
95 {
96#if HAVE_MPI
97 solveSequentialOrParallel_(A, x, b);
98#else
99 solveSequential_(A, x, b);
100#endif
101 return result_.converged;
102 }
103
107 std::string name() const
108 {
109 return "AMG-preconditioned BiCGSTAB solver";
110 }
111
115 const Dune::InverseOperatorResult& result() const
116 {
117 return result_;
118 }
119
120private:
121
122#if HAVE_MPI
123 template<class Matrix, class Vector>
124 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
125 {
126 if constexpr (LinearSolverTraits::canCommunicate)
127 {
128 if (isParallel_)
129 {
130 if (LinearSolverTraits::isNonOverlapping(phelper_->gridView()))
131 {
132 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
133 solveParallel_<PTraits>(A, x, b);
134 }
135 else
136 {
137 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
138 solveParallel_<PTraits>(A, x, b);
139 }
140 }
141 else
142 solveSequential_(A, x, b);
143 }
144 else
145 {
146 solveSequential_(A, x, b);
147 }
148 }
149
150 template<class ParallelTraits, class Matrix, class Vector>
151 void solveParallel_(Matrix& A, Vector& x, Vector& b)
152 {
153 using Comm = typename ParallelTraits::Comm;
154 using LinearOperator = typename ParallelTraits::LinearOperator;
155 using ScalarProduct = typename ParallelTraits::ScalarProduct;
156
157 std::shared_ptr<Comm> comm;
158 std::shared_ptr<LinearOperator> linearOperator;
159 std::shared_ptr<ScalarProduct> scalarProduct;
160 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, comm, linearOperator, scalarProduct, *phelper_);
161
162 using SeqSmoother = Dune::SeqSSOR<Matrix, Vector, Vector>;
163 using Smoother = typename ParallelTraits::template Preconditioner<SeqSmoother>;
164 solveWithAmg_<Smoother>(A, x, b, linearOperator, comm, scalarProduct);
165 }
166#endif // HAVE_MPI
167
168 template<class Matrix, class Vector>
169 void solveSequential_(Matrix& A, Vector& x, Vector& b)
170 {
171 using Comm = Dune::Amg::SequentialInformation;
172 using Traits = typename LinearSolverTraits::template Sequential<Matrix, Vector>;
173 using LinearOperator = typename Traits::LinearOperator;
174 using ScalarProduct = typename Traits::ScalarProduct;
175
176 auto comm = std::make_shared<Comm>();
177 auto linearOperator = std::make_shared<LinearOperator>(A);
178 auto scalarProduct = std::make_shared<ScalarProduct>();
179
180 using Smoother = Dune::SeqSSOR<Matrix, Vector, Vector>;
181 solveWithAmg_<Smoother>(A, x, b, linearOperator, comm, scalarProduct);
182 }
183
184 template<class Smoother, class Matrix, class Vector, class LinearOperator, class Comm, class ScalarProduct>
185 void solveWithAmg_(Matrix& A, Vector& x, Vector& b,
186 std::shared_ptr<LinearOperator>& linearOperator,
187 std::shared_ptr<Comm>& comm,
188 std::shared_ptr<ScalarProduct>& scalarProduct)
189 {
190 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
191 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<Matrix, Dune::Amg::FirstDiagonal>>;
192
195 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
196 params.setDefaultValuesIsotropic(LinearSolverTraits::GridView::dimension);
197 params.setDebugLevel(this->verbosity());
198 Criterion criterion(params);
199 SmootherArgs smootherArgs;
200 smootherArgs.iterations = 1;
201 smootherArgs.relaxationFactor = 1;
202
203 using Amg = Dune::Amg::AMG<LinearOperator, Vector, Smoother, Comm>;
204 auto amg = std::make_shared<Amg>(*linearOperator, criterion, smootherArgs, *comm);
205
206 Dune::BiCGSTABSolver<Vector> solver(*linearOperator, *scalarProduct, *amg, this->residReduction(), this->maxIter(),
207 comm->communicator().rank() == 0 ? this->verbosity() : 0);
208
209 solver.apply(x, b, result_);
210 }
211
212#if HAVE_MPI
213 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> phelper_;
214#endif
215 Dune::InverseOperatorResult result_;
216 bool isParallel_ = false;
217};
218
219} // end namespace Dumux
220
221#endif
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:50
const Dune::InverseOperatorResult & result() const
The result containing the convergence history.
Definition: amgbackend.hh:115
std::string name() const
The name of the solver.
Definition: amgbackend.hh:107
AMGBiCGSTABBackend(const std::string &paramGroup="")
Construct the backend for the sequential case only.
Definition: amgbackend.hh:57
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:72
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: amgbackend.hh:94
Base class for linear solvers.
Definition: solver.hh:37
double residReduction() const
the linear solver residual reduction
Definition: solver.hh:99
int maxIter() const
the maximum number of linear solver iterations
Definition: solver.hh:91
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:79
int verbosity() const
the verbosity level
Definition: solver.hh:83