parallelamgbackend.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef EWOMS_PARALLEL_AMG_BACKEND_HH
28#define EWOMS_PARALLEL_AMG_BACKEND_HH
29
30#include "linalgproperties.hh"
32#include "bicgstabsolver.hh"
33#include "combinedcriterion.hh"
35
36#include <dune/istl/bcrsmatrix.hh>
37#include <dune/istl/paamg/amg.hh>
38#include <dune/istl/paamg/pinfo.hh>
39#include <dune/istl/owneroverlapcopy.hh>
40
41#include <memory>
42#include <tuple>
43#include <utility>
44
45namespace Opm::Linear {
46
47template <class TypeTag>
48class ParallelAmgBackend;
49
50} // namespace Opm::Linear
51
52namespace Opm::Properties {
53
54// Create new type tags
55namespace TTag {
56
58{ using InheritsFrom = std::tuple<ParallelBaseLinearSolver>; };
59
60} // end namespace TTag
61
62template<class TypeTag>
63struct LinearSolverBackend<TypeTag, TTag::ParallelAmgLinearSolver>
65
66} // namespace Opm::Properties
67
68namespace Opm::Parameters {
69
72struct AmgCoarsenTarget { static constexpr int value = 5000; };
73
74}
75
76namespace Opm::Linear {
77
84template <class TypeTag>
86{
88
95
96 using ParallelOperator = typename ParentType::ParallelOperator;
97 using OverlappingVector = typename ParentType::OverlappingVector;
98 using ParallelPreconditioner = typename ParentType::ParallelPreconditioner;
99 using ParallelScalarProduct = typename ParentType::ParallelScalarProduct;
100
101 static constexpr int numEq = getPropValue<TypeTag, Properties::NumEq>();
102 using VectorBlock = Dune::FieldVector<LinearSolverScalar, numEq>;
103 using MatrixBlock = typename SparseMatrixAdapter::MatrixBlock;
104 using IstlMatrix = typename SparseMatrixAdapter::IstlMatrix;
105
106 using Vector = Dune::BlockVector<VectorBlock>;
107
108 // define the smoother used for the AMG and specify its
109 // arguments
110 using SequentialSmoother = Dune::SeqSOR<IstlMatrix, Vector, Vector>;
111// using SequentialSmoother = Dune::SeqSSOR<IstlMatrix,Vector,Vector>;
112// using SequentialSmoother = Dune::SeqJac<IstlMatrix,Vector,Vector>;
113// using SequentialSmoother = Dune::SeqILU<IstlMatrix,Vector,Vector>;
114
115#if HAVE_MPI
116 using OwnerOverlapCopyCommunication = Dune::OwnerOverlapCopyCommunication<Opm::Linear::Index>;
117 using FineOperator = Dune::OverlappingSchwarzOperator<IstlMatrix,
118 Vector,
119 Vector,
120 OwnerOverlapCopyCommunication>;
121 using FineScalarProduct = Dune::OverlappingSchwarzScalarProduct<Vector,
122 OwnerOverlapCopyCommunication>;
123 using ParallelSmoother = Dune::BlockPreconditioner<Vector,
124 Vector,
125 OwnerOverlapCopyCommunication,
126 SequentialSmoother>;
127 using AMG = Dune::Amg::AMG<FineOperator,
128 Vector,
129 ParallelSmoother,
130 OwnerOverlapCopyCommunication>;
131#else
132 using FineOperator = Dune::MatrixAdapter<IstlMatrix, Vector, Vector>;
133 using FineScalarProduct = Dune::SeqScalarProduct<Vector>;
134 using ParallelSmoother = SequentialSmoother;
135 using AMG = Dune::Amg::AMG<FineOperator, Vector, ParallelSmoother>;
136#endif
137
138 using RawLinearSolver = BiCGStabSolver<ParallelOperator,
139 OverlappingVector,
140 AMG> ;
141
142 static_assert(std::is_same<SparseMatrixAdapter, IstlSparseMatrixAdapter<MatrixBlock> >::value,
143 "The ParallelAmgBackend linear solver backend requires the IstlSparseMatrixAdapter");
144
145public:
146 ParallelAmgBackend(const Simulator& simulator)
147 : ParentType(simulator)
148 { }
149
150 static void registerParameters()
151 {
153
154 Parameters::Register<Parameters::LinearSolverMaxError<Scalar>>
155 ("The maximum residual error which the linear solver tolerates "
156 "without giving up");
157 Parameters::Register<Parameters::AmgCoarsenTarget>
158 ("The coarsening target for the agglomerations of "
159 "the AMG preconditioner");
160 }
161
162protected:
164
165 std::shared_ptr<AMG> preparePreconditioner_()
166 {
167#if HAVE_MPI
168 // create and initialize DUNE's OwnerOverlapCopyCommunication
169 // using the domestic overlap
170 istlComm_ = std::make_shared<OwnerOverlapCopyCommunication>(MPI_COMM_WORLD);
171 setupAmgIndexSet_(this->overlappingMatrix_->overlap(), istlComm_->indexSet());
172 istlComm_->remoteIndices().template rebuild<false>();
173#endif
174
175 // create the parallel scalar product and the parallel operator
176#if HAVE_MPI
177 fineOperator_ = std::make_shared<FineOperator>(*this->overlappingMatrix_, *istlComm_);
178#else
179 fineOperator_ = std::make_shared<FineOperator>(*this->overlappingMatrix_);
180#endif
181
182 setupAmg_();
183
184 return amg_;
185 }
186
188 { /* nothing to do */ }
189
190 std::shared_ptr<RawLinearSolver> prepareSolver_(ParallelOperator& parOperator,
191 ParallelScalarProduct& parScalarProduct,
192 AMG& parPreCond)
193 {
194 const auto& gridView = this->simulator_.gridView();
195 using CCC = CombinedCriterion<OverlappingVector, decltype(gridView.comm())>;
196
197 Scalar linearSolverTolerance = Parameters::Get<Parameters::LinearSolverTolerance<Scalar>>();
198 Scalar linearSolverAbsTolerance = Parameters::Get<Parameters::LinearSolverAbsTolerance<Scalar>>();
199 if (linearSolverAbsTolerance < 0.0) {
200 linearSolverAbsTolerance = this->simulator_.model().newtonMethod().tolerance()/100.0;
201 }
202
203 convCrit_.reset(new CCC(gridView.comm(),
204 /*residualReductionTolerance=*/linearSolverTolerance,
205 /*absoluteResidualTolerance=*/linearSolverAbsTolerance,
206 Parameters::Get<Parameters::LinearSolverMaxError<Scalar>>()));
207
208 auto bicgstabSolver =
209 std::make_shared<RawLinearSolver>(parPreCond, *convCrit_, parScalarProduct);
210
211 int verbosity = 0;
212 if (parOperator.overlap().myRank() == 0)
213 verbosity = Parameters::Get<Parameters::LinearSolverVerbosity>();
214 bicgstabSolver->setVerbosity(verbosity);
215 bicgstabSolver->setMaxIterations(Parameters::Get<Parameters::LinearSolverMaxIterations>());
216 bicgstabSolver->setLinearOperator(&parOperator);
217 bicgstabSolver->setRhs(this->overlappingb_);
218
219 return bicgstabSolver;
220 }
221
222 std::pair<bool,int> runSolver_(std::shared_ptr<RawLinearSolver> solver)
223 {
224 bool converged = solver->apply(*this->overlappingx_);
225 return std::make_pair(converged, int(solver->report().iterations()));
226 }
227
229 { /* nothing to do */ }
230
231#if HAVE_MPI
232 template <class ParallelIndexSet>
233 void setupAmgIndexSet_(const Overlap& overlap, ParallelIndexSet& istlIndices)
234 {
235 using GridAttributes = Dune::OwnerOverlapCopyAttributeSet;
236 using GridAttributeSet = Dune::OwnerOverlapCopyAttributeSet::AttributeSet;
237
238 // create DUNE's ParallelIndexSet from a domestic overlap
239 istlIndices.beginResize();
240 for (Index curIdx = 0; static_cast<size_t>(curIdx) < overlap.numDomestic(); ++curIdx) {
241 GridAttributeSet gridFlag =
242 overlap.iAmMasterOf(curIdx)
243 ? GridAttributes::owner
244 : GridAttributes::copy;
245
246 // an index is used by other processes if it is in the
247 // domestic or in the foreign overlap.
248 bool isShared = overlap.isInOverlap(curIdx);
249
250 assert(curIdx == overlap.globalToDomestic(overlap.domesticToGlobal(curIdx)));
251 istlIndices.add(/*globalIdx=*/overlap.domesticToGlobal(curIdx),
252 Dune::ParallelLocalIndex<GridAttributeSet>(static_cast<size_t>(curIdx),
253 gridFlag,
254 isShared));
255 }
256 istlIndices.endResize();
257 }
258#endif
259
260 // trailing return type with decltype used for detecting existence of setUseFixedOrder member function by overloading the setUseFixedOrder function
261 template <typename C>
262 auto setUseFixedOrder(C criterion, bool booleanValue) -> decltype(criterion.setUseFixedOrder(booleanValue))
263 {
264 return criterion.setUseFixedOrder(booleanValue); // Set flag to ensure that the matrices in the AMG hierarchy are constructed with deterministic indices.
265 }
266 template <typename C>
267 void setUseFixedOrder(C, ...)
268 {
269 // do nothing, since the function setUseFixedOrder does not exist yet
270 }
271
273 {
274 if (amg_)
275 amg_.reset();
276
277 int verbosity = 0;
278 if (this->simulator_.vanguard().gridView().comm().rank() == 0)
279 verbosity = Parameters::Get<Parameters::LinearSolverVerbosity>();
280
281 using SmootherArgs = typename Dune::Amg::SmootherTraits<ParallelSmoother>::Arguments;
282
283 SmootherArgs smootherArgs;
284 smootherArgs.iterations = 1;
285 smootherArgs.relaxationFactor = 1.0;
286
287 // specify the coarsen criterion:
288 //
289 // using CoarsenCriterion =
290 // Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<IstlMatrix,
291 // Dune::Amg::FirstDiagonal>>
292 using CoarsenCriterion = Dune::Amg::
293 CoarsenCriterion<Dune::Amg::SymmetricCriterion<IstlMatrix, Dune::Amg::FrobeniusNorm> >;
294 int coarsenTarget = Parameters::Get<Parameters::AmgCoarsenTarget>();
295 CoarsenCriterion coarsenCriterion(/*maxLevel=*/15, coarsenTarget);
296 coarsenCriterion.setDefaultValuesAnisotropic(GridView::dimension,
297 /*aggregateSizePerDim=*/3);
298 if (verbosity > 0)
299 coarsenCriterion.setDebugLevel(1);
300 else
301 coarsenCriterion.setDebugLevel(0); // make the AMG shut up
302
303 // reduce the minium coarsen rate (default is 1.2)
304 coarsenCriterion.setMinCoarsenRate(1.05);
305 // coarsenCriterion.setAccumulate(Dune::Amg::noAccu);
306 coarsenCriterion.setAccumulate(Dune::Amg::atOnceAccu);
307 coarsenCriterion.setSkipIsolated(false);
308 setUseFixedOrder(coarsenCriterion, true); // If possible, set flag to ensure that the matrices in the AMG hierarchy are constructed with deterministic indices.
309
310// instantiate the AMG preconditioner
311#if HAVE_MPI
312 amg_ = std::make_shared<AMG>(*fineOperator_, coarsenCriterion, smootherArgs, *istlComm_);
313#else
314 amg_ = std::make_shared<AMG>(*fineOperator_, coarsenCriterion, smootherArgs);
315#endif
316 }
317
318 std::unique_ptr<ConvergenceCriterion<OverlappingVector> > convCrit_;
319
320 std::shared_ptr<FineOperator> fineOperator_;
321 std::shared_ptr<AMG> amg_;
322
323#if HAVE_MPI
324 std::shared_ptr<OwnerOverlapCopyCommunication> istlComm_;
325#endif
326};
327
328} // namespace Opm::Linear
329
330#endif
Implements a preconditioned stabilized BiCG linear solver.
Definition: bicgstabsolver.hh:54
Convergence criterion which looks at the absolute value of the residual and fails if the linear solve...
Definition: combinedcriterion.hh:56
Provides a linear solver backend using the parallel algebraic multi-grid (AMG) linear solver from DUN...
Definition: parallelamgbackend.hh:86
std::shared_ptr< AMG > preparePreconditioner_()
Definition: parallelamgbackend.hh:165
std::shared_ptr< RawLinearSolver > prepareSolver_(ParallelOperator &parOperator, ParallelScalarProduct &parScalarProduct, AMG &parPreCond)
Definition: parallelamgbackend.hh:190
void cleanupSolver_()
Definition: parallelamgbackend.hh:228
std::pair< bool, int > runSolver_(std::shared_ptr< RawLinearSolver > solver)
Definition: parallelamgbackend.hh:222
auto setUseFixedOrder(C criterion, bool booleanValue) -> decltype(criterion.setUseFixedOrder(booleanValue))
Definition: parallelamgbackend.hh:262
std::unique_ptr< ConvergenceCriterion< OverlappingVector > > convCrit_
Definition: parallelamgbackend.hh:318
void setupAmgIndexSet_(const Overlap &overlap, ParallelIndexSet &istlIndices)
Definition: parallelamgbackend.hh:233
friend ParentType
Definition: parallelamgbackend.hh:163
std::shared_ptr< FineOperator > fineOperator_
Definition: parallelamgbackend.hh:320
static void registerParameters()
Definition: parallelamgbackend.hh:150
std::shared_ptr< AMG > amg_
Definition: parallelamgbackend.hh:321
void setUseFixedOrder(C,...)
Definition: parallelamgbackend.hh:267
void cleanupPreconditioner_()
Definition: parallelamgbackend.hh:187
ParallelAmgBackend(const Simulator &simulator)
Definition: parallelamgbackend.hh:146
void setupAmg_()
Definition: parallelamgbackend.hh:272
std::shared_ptr< OwnerOverlapCopyCommunication > istlComm_
Definition: parallelamgbackend.hh:324
Provides the common code which is required by most linear solvers.
Definition: parallelbasebackend.hh:109
GetPropType< TypeTag, Properties::OverlappingVector > OverlappingVector
Definition: parallelbasebackend.hh:122
Opm::Linear::OverlappingScalarProduct< OverlappingVector, Overlap > ParallelScalarProduct
Definition: parallelbasebackend.hh:129
OverlappingVector * overlappingx_
Definition: parallelbasebackend.hh:383
static void registerParameters()
Register all run-time parameters for the linear solver.
Definition: parallelbasebackend.hh:153
Opm::Linear::OverlappingOperator< OverlappingMatrix, OverlappingVector, OverlappingVector > ParallelOperator
Definition: parallelbasebackend.hh:132
OverlappingVector * overlappingb_
Definition: parallelbasebackend.hh:382
OverlappingMatrix * overlappingMatrix_
Definition: parallelbasebackend.hh:381
Opm::Linear::OverlappingPreconditioner< SequentialPreconditioner, Overlap > ParallelPreconditioner
Definition: parallelbasebackend.hh:128
const Simulator & simulator_
Definition: parallelbasebackend.hh:377
Declares the properties required by the black oil model.
Definition: bicgstabsolver.hh:42
int Index
The type of an index of a degree of freedom.
Definition: overlaptypes.hh:44
Definition: blackoilnewtonmethodparameters.hh:31
Definition: blackoilmodel.hh:72
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
Definition: parallelamgbackend.hh:72
static constexpr int value
Definition: parallelamgbackend.hh:72
The type of the linear solver to be used.
Definition: linalgproperties.hh:38
Definition: parallelamgbackend.hh:58
std::tuple< ParallelBaseLinearSolver > InheritsFrom
Definition: parallelamgbackend.hh:58