Dune Core Modules (2.6.0)

kamg.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_AMG_KAMG_HH
4#define DUNE_AMG_KAMG_HH
5
7#include "amg.hh"
8
9namespace Dune
10{
11 namespace Amg
12 {
13
28 template<class AMG>
30 : public Preconditioner<typename AMG::Domain,typename AMG::Range>
31 {
33 typedef typename AMG::Domain Domain;
35 typedef typename AMG::Range Range;
36 public:
37
40 {
41 return amg_.category();
42 };
43
52 : amg_(amg), coarseSolver_(coarseSolver)
53 {}
54
56 void pre(typename AMG::Domain& x, typename AMG::Range& b)
57 {
59 }
60
62 void post(typename AMG::Domain& x)
63 {
65 }
66
68 void apply(typename AMG::Domain& v, const typename AMG::Range& d)
69 {
70 // Copy data
71 *levelContext_->update=0;
72 *levelContext_->rhs = d;
73 *levelContext_->lhs = v;
74
75 presmooth(*levelContext_, amg_.preSteps_);
76 bool processFineLevel =
77 amg_.moveToCoarseLevel(*levelContext_);
78
79 if(processFineLevel) {
80 typename AMG::Range b=*levelContext_->rhs;
81 typename AMG::Domain x=*levelContext_->update;
83 coarseSolver_->apply(x, b, res);
84 *levelContext_->update=x;
85 }
86
87 amg_.moveToFineLevel(*levelContext_, processFineLevel);
88
89 postsmooth(*levelContext_, amg_.postSteps_);
90 v=*levelContext_->update;
91 }
92
98 {
99 return coarseSolver_;
100 }
101
106 void setLevelContext(std::shared_ptr<typename AMG::LevelContext> p)
107 {
108 levelContext_=p;
109 }
110
113 {}
114
115 private:
117 AMG& amg_;
119 std::shared_ptr<InverseOperator<Domain,Range> > coarseSolver_;
121 std::shared_ptr<typename AMG::LevelContext> levelContext_;
122 };
123
124
125
139 template<class M, class X, class S, class PI=SequentialInformation,
140 class K=GeneralizedPCGSolver<X>, class A=std::allocator<X> >
141 class KAMG : public Preconditioner<X,X>
142 {
143 public:
147 typedef K KrylovSolver;
157 typedef typename Amg::Operator Operator;
159 typedef typename Amg::Domain Domain;
161 typedef typename Amg::Range Range;
166
169 {
170 return amg.category();
171 };
172
186 KAMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
187 const SmootherArgs& smootherArgs, std::size_t gamma,
188 std::size_t preSmoothingSteps=1, std::size_t postSmoothingSteps=1,
189 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1) DUNE_DEPRECATED;
190
202 KAMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
203 const SmootherArgs& smootherArgs, const Parameters& parms,
204 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1);
205
226 template<class C>
227 KAMG(const Operator& fineOperator, const C& criterion,
228 const SmootherArgs& smootherArgs, std::size_t gamma,
229 std::size_t preSmoothingSteps=1, std::size_t postSmoothingSteps=1,
230 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1,
232
246 template<class C>
247 KAMG(const Operator& fineOperator, const C& criterion,
248 const SmootherArgs& smootherArgs=SmootherArgs(),
249 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1,
251
253 void pre(Domain& x, Range& b);
255 void post(Domain& x);
257 void apply(Domain& v, const Range& d);
258
259 std::size_t maxlevels();
260
261 private:
263 Amg amg;
264
266 std::size_t maxLevelKrylovSteps;
267
269 double levelDefectReduction;
270
272 std::vector<std::shared_ptr<typename Amg::ScalarProduct> > scalarproducts;
273
275 std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > > ksolvers;
276 };
277
278 template<class M, class X, class S, class P, class K, class A>
279 KAMG<M,X,S,P,K,A>::KAMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
280 const SmootherArgs& smootherArgs,
281 std::size_t gamma, std::size_t preSmoothingSteps,
282 std::size_t postSmoothingSteps,
283 std::size_t ksteps, double reduction)
284 : amg(matrices, coarseSolver, smootherArgs, gamma, preSmoothingSteps,
285 postSmoothingSteps), maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
286 {}
287
288 template<class M, class X, class S, class P, class K, class A>
290 const SmootherArgs& smootherArgs, const Parameters& params,
291 std::size_t ksteps, double reduction)
292 : amg(matrices, coarseSolver, smootherArgs, params),
293 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
294 {}
295
296 template<class M, class X, class S, class P, class K, class A>
297 template<class C>
298 KAMG<M,X,S,P,K,A>::KAMG(const Operator& fineOperator, const C& criterion,
299 const SmootherArgs& smootherArgs, std::size_t gamma,
300 std::size_t preSmoothingSteps, std::size_t postSmoothingSteps,
301 std::size_t ksteps, double reduction,
302 const ParallelInformation& pinfo)
303 : amg(fineOperator, criterion, smootherArgs, gamma, preSmoothingSteps,
304 postSmoothingSteps, false, pinfo), maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
305 {}
306
307 template<class M, class X, class S, class P, class K, class A>
308 template<class C>
309 KAMG<M,X,S,P,K,A>::KAMG(const Operator& fineOperator, const C& criterion,
310 const SmootherArgs& smootherArgs,
311 std::size_t ksteps, double reduction,
312 const ParallelInformation& pinfo)
313 : amg(fineOperator, criterion, smootherArgs, pinfo),
314 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
315 {}
316
317
318 template<class M, class X, class S, class P, class K, class A>
320 {
321 amg.pre(x,b);
322 scalarproducts.reserve(amg.levels());
323 ksolvers.reserve(amg.levels());
324
325 typename OperatorHierarchy::ParallelMatrixHierarchy::Iterator
326 matrix = amg.matrices_->matrices().coarsest();
328 pinfo = amg.matrices_->parallelInformation().coarsest();
329 bool hasCoarsest=(amg.levels()==amg.maxlevels());
330
331 if(hasCoarsest) {
332 if(matrix==amg.matrices_->matrices().finest())
333 return;
334 --matrix;
335 --pinfo;
336 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, amg.solver_)));
337 }else
338 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, std::shared_ptr<InverseOperator<Domain,Range> >())));
339
340 std::ostringstream s;
341
342 if(matrix!=amg.matrices_->matrices().finest())
343 while(true) {
344 scalarproducts.push_back(createScalarProduct<X>(*pinfo,category()));
345 std::shared_ptr<InverseOperator<Domain,Range> > ks =
346 std::shared_ptr<InverseOperator<Domain,Range> >(new KrylovSolver(*matrix, *(scalarproducts.back()),
347 *(ksolvers.back()), levelDefectReduction,
348 maxLevelKrylovSteps, 0));
349 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, ks)));
350 --matrix;
351 --pinfo;
352 if(matrix==amg.matrices_->matrices().finest())
353 break;
354 }
355 }
356
357
358 template<class M, class X, class S, class P, class K, class A>
360 {
361 amg.post(x);
362
363 }
364
365 template<class M, class X, class S, class P, class K, class A>
367 {
368 if(ksolvers.size()==0)
369 {
370 Range td=d;
372 amg.solver_->apply(v,td,res);
373 }else
374 {
375 typedef typename Amg::LevelContext LevelContext;
376 std::shared_ptr<LevelContext> levelContext(new LevelContext);
377 amg.initIteratorsWithFineLevel(*levelContext);
378 typedef typename std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > >::iterator Iter;
379 for(Iter solver=ksolvers.begin(); solver!=ksolvers.end(); ++solver)
380 (*solver)->setLevelContext(levelContext);
381 ksolvers.back()->apply(v,d);
382 }
383 }
384
385 template<class M, class X, class S, class P, class K, class A>
386 std::size_t KAMG<M,X,S,P,K,A>::maxlevels()
387 {
388 return amg.maxlevels();
389 }
390
392 } // Amg
393} // Dune
394
395#endif
The AMG preconditioner.
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:59
LevelIterator< Hierarchy< ParallelInformation, Allocator >, ParallelInformation > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:250
an algebraic multigrid method using a Krylov-cycle.
Definition: kamg.hh:142
Amg::Domain Domain
the type of the domain.
Definition: kamg.hh:159
Amg::SmootherArgs SmootherArgs
The type of the arguments for construction of the smoothers.
Definition: kamg.hh:155
Amg::ParallelInformation ParallelInformation
the type of the parallelinformation to use.
Definition: kamg.hh:153
Amg::CoarseSolver CoarseSolver
The type of the coarse solver.
Definition: kamg.hh:151
Amg::OperatorHierarchy OperatorHierarchy
The type of the hierarchy of operators.
Definition: kamg.hh:149
Amg::Range Range
The type of the range.
Definition: kamg.hh:161
Amg::ScalarProduct ScalarProduct
The type of the scalar product.
Definition: kamg.hh:165
AMG< M, X, S, PI, A > Amg
The type of the underlying AMG.
Definition: kamg.hh:145
Amg::Operator Operator
the type of the lineatr operator.
Definition: kamg.hh:157
Amg::ParallelInformationHierarchy ParallelInformationHierarchy
The type of the hierarchy of parallel information.
Definition: kamg.hh:163
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: kamg.hh:168
K KrylovSolver
The type of the Krylov solver for the cycle.
Definition: kamg.hh:147
Two grid operator for AMG with Krylov cycle.
Definition: kamg.hh:31
void pre(typename AMG::Domain &x, typename AMG::Range &b)
Prepare the preconditioner.
Definition: kamg.hh:56
KAmgTwoGrid(AMG &amg, std::shared_ptr< InverseOperator< Domain, Range > > coarseSolver)
Constructor.
Definition: kamg.hh:51
~KAmgTwoGrid()
Destructor.
Definition: kamg.hh:112
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:39
void setLevelContext(std::shared_ptr< typename AMG::LevelContext > p)
Set the level context pointer.
Definition: kamg.hh:106
InverseOperator< Domain, Range > * coarseSolver()
Get a pointer to the coarse grid solver.
Definition: kamg.hh:97
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:68
The hierarchies build by the coarsening process.
Definition: hierarchy.hh:310
All parameters for AMG.
Definition: parameters.hh:391
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1363
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:30
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:46
#define DUNE_DEPRECATED
Mark some entity as deprecated.
Definition: deprecated.hh:84
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
void apply(Domain &v, const Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: kamg.hh:366
X Domain
The domain type.
Definition: amg.hh:81
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:94
M Operator
The matrix operator type.
Definition: amg.hh:67
void post(Domain &x)
Clean up.
Definition: kamg.hh:359
X Range
The range type.
Definition: amg.hh:83
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:454
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:476
KAMG(const OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const SmootherArgs &smootherArgs, std::size_t gamma, std::size_t preSmoothingSteps=1, std::size_t postSmoothingSteps=1, std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1)
Construct a new amg with a specific coarse solver.
Definition: kamg.hh:279
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: kamg.hh:319
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: amg.hh:138
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:74
Dune namespace.
Definition: alignedallocator.hh:10
STL namespace.
Define general preconditioner interface.
Statistics about the application of an inverse operator.
Definition: solver.hh:41
Category
Definition: solvercategory.hh:21
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)