3.5-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#include <dune/common/parallel/mpicommunication.hh>
33#include <dune/common/parallel/mpihelper.hh>
34#include <dune/grid/common/capabilities.hh>
35#include <dune/istl/paamg/amg.hh>
36#include <dune/istl/paamg/pinfo.hh>
37#include <dune/istl/solvers.hh>
38
41
42namespace Dumux {
43
49template <class LinearSolverTraits>
51{
52public:
58 AMGBiCGSTABBackend(const std::string& paramGroup = "")
60 , isParallel_(Dune::MPIHelper::getCommunication().size() > 1)
61 {
62 if (isParallel_)
63 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
64
65 checkAvailabilityOfDirectSolver_();
66 }
67
75 AMGBiCGSTABBackend(const typename LinearSolverTraits::GridView& gridView,
76 const typename LinearSolverTraits::DofMapper& dofMapper,
77 const std::string& paramGroup = "")
79#if HAVE_MPI
80 , isParallel_(Dune::MPIHelper::getCommunication().size() > 1)
81#endif
82 {
83#if HAVE_MPI
84 if (isParallel_)
85 phelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
86#endif
87 checkAvailabilityOfDirectSolver_();
88 }
89
96 void updateAfterGridAdaption(const typename LinearSolverTraits::GridView& gridView,
97 const typename LinearSolverTraits::DofMapper& dofMapper)
98 {
99#if HAVE_MPI
100 if (isParallel_)
101 phelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
102#endif
103 }
104
112 template<class Matrix, class Vector>
113 bool solve(Matrix& A, Vector& x, Vector& b)
114 {
115#if HAVE_MPI
116 solveSequentialOrParallel_(A, x, b);
117#else
118 solveSequential_(A, x, b);
119#endif
120 return result_.converged;
121 }
122
126 std::string name() const
127 {
128 return "AMG-preconditioned BiCGSTAB solver";
129 }
130
134 const Dune::InverseOperatorResult& result() const
135 {
136 return result_;
137 }
138
139private:
141 void checkAvailabilityOfDirectSolver_()
142 {
143#if !HAVE_SUPERLU && !HAVE_UMFPACK
144 std::cout << "\nAMGBiCGSTABBackend: No direct solver backend found. Using iterative solver as coarse grid solver.\n"
145 << "Note that dune-istl currently hard-codes a tolerance of 1e-2 for the iterative coarse grid solver.\n"
146 << "This may result in reduced accuracy or performance depending on your setup.\nConsider installing "
147 << "UMFPack (SuiteSparse) or SuperLU or apply the istl patch, see dumux/patches/README.md." << std::endl;
148#endif
149 }
150
151#if HAVE_MPI
152 template<class Matrix, class Vector>
153 void solveSequentialOrParallel_(Matrix& A, Vector& x, Vector& b)
154 {
156 {
157 if (isParallel_)
158 {
159 if (LinearSolverTraits::isNonOverlapping(phelper_->gridView()))
160 {
161 using PTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, Vector>;
162 solveParallel_<PTraits>(A, x, b);
163 }
164 else
165 {
166 using PTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, Vector>;
167 solveParallel_<PTraits>(A, x, b);
168 }
169 }
170 else
171 solveSequential_(A, x, b);
172 }
173 else
174 {
175 solveSequential_(A, x, b);
176 }
177 }
178
179 template<class ParallelTraits, class Matrix, class Vector>
180 void solveParallel_(Matrix& A, Vector& x, Vector& b)
181 {
182 using Comm = typename ParallelTraits::Comm;
183 using LinearOperator = typename ParallelTraits::LinearOperator;
184 using ScalarProduct = typename ParallelTraits::ScalarProduct;
185
186 std::shared_ptr<Comm> comm;
187 std::shared_ptr<LinearOperator> linearOperator;
188 std::shared_ptr<ScalarProduct> scalarProduct;
189 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, comm, linearOperator, scalarProduct, *phelper_);
190
191 using SeqSmoother = Dune::SeqSSOR<Matrix, Vector, Vector>;
192 using Smoother = typename ParallelTraits::template Preconditioner<SeqSmoother>;
193 solveWithAmg_<Smoother>(A, x, b, linearOperator, comm, scalarProduct);
194 }
195#endif // HAVE_MPI
196
197 template<class Matrix, class Vector>
198 void solveSequential_(Matrix& A, Vector& x, Vector& b)
199 {
200 using Comm = Dune::Amg::SequentialInformation;
201 using Traits = typename LinearSolverTraits::template Sequential<Matrix, Vector>;
202 using LinearOperator = typename Traits::LinearOperator;
203 using ScalarProduct = typename Traits::ScalarProduct;
204
205 auto comm = std::make_shared<Comm>();
206 auto linearOperator = std::make_shared<LinearOperator>(A);
207 auto scalarProduct = std::make_shared<ScalarProduct>();
208
209 using Smoother = Dune::SeqSSOR<Matrix, Vector, Vector>;
210 solveWithAmg_<Smoother>(A, x, b, linearOperator, comm, scalarProduct);
211 }
212
213 template<class Smoother, class Matrix, class Vector, class LinearOperator, class Comm, class ScalarProduct>
214 void solveWithAmg_(Matrix& A, Vector& x, Vector& b,
215 std::shared_ptr<LinearOperator>& linearOperator,
216 std::shared_ptr<Comm>& comm,
217 std::shared_ptr<ScalarProduct>& scalarProduct)
218 {
219 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
220 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<Matrix, Dune::Amg::FirstDiagonal>>;
221
224 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
225 params.setDefaultValuesIsotropic(LinearSolverTraits::GridView::dimension);
226 params.setDebugLevel(this->verbosity());
227 Criterion criterion(params);
228 SmootherArgs smootherArgs;
229 smootherArgs.iterations = 1;
230 smootherArgs.relaxationFactor = 1;
231
232 using Amg = Dune::Amg::AMG<LinearOperator, Vector, Smoother, Comm>;
233 auto amg = std::make_shared<Amg>(*linearOperator, criterion, smootherArgs, *comm);
234
235 Dune::BiCGSTABSolver<Vector> solver(*linearOperator, *scalarProduct, *amg, this->residReduction(), this->maxIter(),
236 comm->communicator().rank() == 0 ? this->verbosity() : 0);
237
238 solver.apply(x, b, result_);
239 }
240
241#if HAVE_MPI
242 std::unique_ptr<ParallelISTLHelper<LinearSolverTraits>> phelper_;
243#endif
244 Dune::InverseOperatorResult result_;
245 bool isParallel_ = false;
246};
247
248} // end namespace Dumux
249
250#endif
Provides a helper class for nonoverlapping decomposition.
Base class for linear solvers.
Definition: adapt.hh:29
static constexpr bool canCommunicate
Definition: gridcapabilities.hh:63
Definition: common/pdesolver.hh:36
A linear solver based on the ISTL AMG preconditioner and the ISTL BiCGSTAB solver.
Definition: amgbackend.hh:51
const Dune::InverseOperatorResult & result() const
The result containing the convergence history.
Definition: amgbackend.hh:134
std::string name() const
The name of the solver.
Definition: amgbackend.hh:126
AMGBiCGSTABBackend(const std::string &paramGroup="")
Construct the backend for the sequential case only.
Definition: amgbackend.hh:58
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:75
bool solve(Matrix &A, Vector &x, Vector &b)
Solve a linear system.
Definition: amgbackend.hh:113
void updateAfterGridAdaption(const typename LinearSolverTraits::GridView &gridView, const typename LinearSolverTraits::DofMapper &dofMapper)
Update the solver after grid adaption.
Definition: amgbackend.hh:96
Base class for linear solvers.
Definition: solver.hh:37
Scalar 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