Dune Core Modules (unstable)

kamg.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5#ifndef DUNE_AMG_KAMG_HH
6#define DUNE_AMG_KAMG_HH
7
9#include "amg.hh"
10
11namespace Dune
12{
13 namespace Amg
14 {
15
30 template<class AMG>
32 : public Preconditioner<typename AMG::Domain,typename AMG::Range>
33 {
35 typedef typename AMG::Domain Domain;
37 typedef typename AMG::Range Range;
38 public:
39
42 {
43 return amg_.category();
44 };
45
54 : amg_(amg), coarseSolver_(coarseSolver)
55 {}
56
58 void pre([[maybe_unused]] typename AMG::Domain& x, [[maybe_unused]] typename AMG::Range& b)
59 {}
60
62 void post([[maybe_unused]] typename AMG::Domain& x)
63 {}
64
66 void apply(typename AMG::Domain& v, const typename AMG::Range& d)
67 {
68 // Copy data
69 *levelContext_->update=0;
70 *levelContext_->rhs = d;
71 *levelContext_->lhs = v;
72
73 presmooth(*levelContext_, amg_.preSteps_);
74 bool processFineLevel =
75 amg_.moveToCoarseLevel(*levelContext_);
76
77 if(processFineLevel) {
78 typename AMG::Range b=*levelContext_->rhs;
79 typename AMG::Domain x=*levelContext_->update;
81 coarseSolver_->apply(x, b, res);
82 *levelContext_->update=x;
83 }
84
85 amg_.moveToFineLevel(*levelContext_, processFineLevel);
86
87 postsmooth(*levelContext_, amg_.postSteps_);
88 v=*levelContext_->update;
89 }
90
96 {
97 return coarseSolver_;
98 }
99
104 void setLevelContext(std::shared_ptr<typename AMG::LevelContext> p)
105 {
106 levelContext_=p;
107 }
108
111 {}
112
113 private:
115 AMG& amg_;
117 std::shared_ptr<InverseOperator<Domain,Range> > coarseSolver_;
119 std::shared_ptr<typename AMG::LevelContext> levelContext_;
120 };
121
122
123
137 template<class M, class X, class S, class PI=SequentialInformation,
138 class K=GeneralizedPCGSolver<X>, class A=std::allocator<X> >
139 class KAMG : public Preconditioner<X,X>
140 {
141 public:
145 typedef K KrylovSolver;
155 typedef typename Amg::Operator Operator;
157 typedef typename Amg::Domain Domain;
159 typedef typename Amg::Range Range;
164
167 {
168 return amg.category();
169 };
170
182 KAMG(OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
183 const SmootherArgs& smootherArgs, const Parameters& parms,
184 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1);
185
199 template<class C>
200 KAMG(const Operator& fineOperator, const C& criterion,
201 const SmootherArgs& smootherArgs=SmootherArgs(),
202 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1,
204
206 void pre(Domain& x, Range& b);
208 void post(Domain& x);
210 void apply(Domain& v, const Range& d);
211
212 std::size_t maxlevels();
213
214 private:
216 Amg amg;
217
219 std::size_t maxLevelKrylovSteps;
220
222 double levelDefectReduction;
223
225 std::vector<std::shared_ptr<typename Amg::ScalarProduct> > scalarproducts;
226
228 std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > > ksolvers;
229 };
230
231
232 template<class M, class X, class S, class P, class K, class A>
234 const SmootherArgs& smootherArgs, const Parameters& params,
235 std::size_t ksteps, double reduction)
236 : amg(matrices, coarseSolver, smootherArgs, params),
237 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
238 {}
239
240
241 template<class M, class X, class S, class P, class K, class A>
242 template<class C>
243 KAMG<M,X,S,P,K,A>::KAMG(const Operator& fineOperator, const C& criterion,
244 const SmootherArgs& smootherArgs,
245 std::size_t ksteps, double reduction,
246 const ParallelInformation& pinfo)
247 : amg(fineOperator, criterion, smootherArgs, pinfo),
248 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
249 {}
250
251
252 template<class M, class X, class S, class P, class K, class A>
254 {
255 amg.pre(x,b);
256 scalarproducts.reserve(amg.levels());
257 ksolvers.reserve(amg.levels());
258
259 typename OperatorHierarchy::ParallelMatrixHierarchy::Iterator
260 matrix = amg.matrices_->matrices().coarsest();
262 pinfo = amg.matrices_->parallelInformation().coarsest();
263 bool hasCoarsest=(amg.levels()==amg.maxlevels());
264
265 if(hasCoarsest) {
266 if(matrix==amg.matrices_->matrices().finest())
267 return;
268 --matrix;
269 --pinfo;
270 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, amg.solver_)));
271 }else
272 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, std::shared_ptr<InverseOperator<Domain,Range> >())));
273
274 std::ostringstream s;
275
276 if(matrix!=amg.matrices_->matrices().finest())
277 while(true) {
278 scalarproducts.push_back(createScalarProduct<X>(*pinfo,category()));
279 std::shared_ptr<InverseOperator<Domain,Range> > ks =
280 std::shared_ptr<InverseOperator<Domain,Range> >(new KrylovSolver(*matrix, *(scalarproducts.back()),
281 *(ksolvers.back()), levelDefectReduction,
282 maxLevelKrylovSteps, 0));
283 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, ks)));
284 --matrix;
285 --pinfo;
286 if(matrix==amg.matrices_->matrices().finest())
287 break;
288 }
289 }
290
291
292 template<class M, class X, class S, class P, class K, class A>
294 {
295 amg.post(x);
296
297 }
298
299 template<class M, class X, class S, class P, class K, class A>
301 {
302 if(ksolvers.size()==0)
303 {
304 Range td=d;
306 amg.solver_->apply(v,td,res);
307 }else
308 {
309 typedef typename Amg::LevelContext LevelContext;
310 std::shared_ptr<LevelContext> levelContext(new LevelContext);
311 amg.initIteratorsWithFineLevel(*levelContext);
312 typedef typename std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > >::iterator Iter;
313 for(Iter solver=ksolvers.begin(); solver!=ksolvers.end(); ++solver)
314 (*solver)->setLevelContext(levelContext);
315 ksolvers.back()->apply(v,d);
316 }
317 }
318
319 template<class M, class X, class S, class P, class K, class A>
320 std::size_t KAMG<M,X,S,P,K,A>::maxlevels()
321 {
322 return amg.maxlevels();
323 }
324
326 } // Amg
327} // Dune
328
329#endif
The AMG preconditioner.
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:66
LevelIterator< Hierarchy< ParallelInformation, Allocator >, ParallelInformation > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:220
an algebraic multigrid method using a Krylov-cycle.
Definition: kamg.hh:140
Amg::Domain Domain
the type of the domain.
Definition: kamg.hh:157
Amg::SmootherArgs SmootherArgs
The type of the arguments for construction of the smoothers.
Definition: kamg.hh:153
Amg::ParallelInformation ParallelInformation
the type of the parallelinformation to use.
Definition: kamg.hh:151
Amg::CoarseSolver CoarseSolver
The type of the coarse solver.
Definition: kamg.hh:149
Amg::OperatorHierarchy OperatorHierarchy
The type of the hierarchy of operators.
Definition: kamg.hh:147
Amg::Range Range
The type of the range.
Definition: kamg.hh:159
Amg::ScalarProduct ScalarProduct
The type of the scalar product.
Definition: kamg.hh:163
AMG< M, X, S, PI, A > Amg
The type of the underlying AMG.
Definition: kamg.hh:143
Amg::Operator Operator
the type of the lineatr operator.
Definition: kamg.hh:155
Amg::ParallelInformationHierarchy ParallelInformationHierarchy
The type of the hierarchy of parallel information.
Definition: kamg.hh:161
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: kamg.hh:166
K KrylovSolver
The type of the Krylov solver for the cycle.
Definition: kamg.hh:145
Two grid operator for AMG with Krylov cycle.
Definition: kamg.hh:33
void pre(typename AMG::Domain &x, typename AMG::Range &b)
Prepare the preconditioner.
Definition: kamg.hh:58
KAmgTwoGrid(AMG &amg, std::shared_ptr< InverseOperator< Domain, Range > > coarseSolver)
Constructor.
Definition: kamg.hh:53
~KAmgTwoGrid()
Destructor.
Definition: kamg.hh:110
void post(typename AMG::Domain &x)
Clean up.
Definition: kamg.hh:62
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: kamg.hh:41
void setLevelContext(std::shared_ptr< typename AMG::LevelContext > p)
Set the level context pointer.
Definition: kamg.hh:104
InverseOperator< Domain, Range > * coarseSolver()
Get a pointer to the coarse grid solver.
Definition: kamg.hh:95
void apply(typename AMG::Domain &v, const typename AMG::Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: kamg.hh:66
The hierarchies build by the coarsening process.
Definition: matrixhierarchy.hh:61
All parameters for AMG.
Definition: parameters.hh:416
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1307
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:33
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:52
void apply(Domain &v, const Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: kamg.hh:300
X Domain
The domain type.
Definition: amg.hh:88
KAMG(OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const SmootherArgs &smootherArgs, const Parameters &parms, std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1)
Construct a new amg with a specific coarse solver.
Definition: kamg.hh:233
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:101
M Operator
The matrix operator type.
Definition: amg.hh:74
void post(Domain &x)
Clean up.
Definition: kamg.hh:293
X Range
The range type.
Definition: amg.hh:90
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:406
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:428
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: kamg.hh:253
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: amg.hh:195
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:81
Dune namespace.
Definition: alignedallocator.hh:13
Define general preconditioner interface.
Statistics about the application of an inverse operator.
Definition: solver.hh:50
Category
Definition: solvercategory.hh:23
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)