Dune Core Modules (2.5.2)

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
38 enum {
41 };
42
51 : amg_(amg), coarseSolver_(coarseSolver)
52 {}
53
55 void pre(typename AMG::Domain& x, typename AMG::Range& b)
56 {
58 }
59
61 void post(typename AMG::Domain& x)
62 {
64 }
65
67 void apply(typename AMG::Domain& v, const typename AMG::Range& d)
68 {
69 // Copy data
70 *levelContext_->update=0;
71 *levelContext_->rhs = d;
72 *levelContext_->lhs = v;
73
74 presmooth(*levelContext_, amg_.preSteps_);
75 bool processFineLevel =
76 amg_.moveToCoarseLevel(*levelContext_);
77
78 if(processFineLevel) {
79 typename AMG::Range b=*levelContext_->rhs;
80 typename AMG::Domain x=*levelContext_->update;
82 coarseSolver_->apply(x, b, res);
83 *levelContext_->update=x;
84 }
85
86 amg_.moveToFineLevel(*levelContext_, processFineLevel);
87
88 postsmooth(*levelContext_, amg_.postSteps_);
89 v=*levelContext_->update;
90 }
91
97 {
98 return coarseSolver_;
99 }
100
105 void setLevelContext(std::shared_ptr<typename AMG::LevelContext> p)
106 {
107 levelContext_=p;
108 }
109
112 {}
113
114 private:
116 AMG& amg_;
118 std::shared_ptr<InverseOperator<Domain,Range> > coarseSolver_;
120 std::shared_ptr<typename AMG::LevelContext> levelContext_;
121 };
122
123
124
138 template<class M, class X, class S, class PI=SequentialInformation,
139 class K=GeneralizedPCGSolver<X>, class A=std::allocator<X> >
140 class KAMG : public Preconditioner<X,X>
141 {
142 public:
146 typedef K KrylovSolver;
156 typedef typename Amg::Operator Operator;
158 typedef typename Amg::Domain Domain;
160 typedef typename Amg::Range Range;
164 typedef typename Amg::ScalarProduct ScalarProduct;
165
166 enum {
169 };
183 KAMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
184 const SmootherArgs& smootherArgs, std::size_t gamma,
185 std::size_t preSmoothingSteps=1, std::size_t postSmoothingSteps=1,
186 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1) DUNE_DEPRECATED;
187
199 KAMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
200 const SmootherArgs& smootherArgs, const Parameters& parms,
201 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1);
202
223 template<class C>
224 KAMG(const Operator& fineOperator, const C& criterion,
225 const SmootherArgs& smootherArgs, std::size_t gamma,
226 std::size_t preSmoothingSteps=1, std::size_t postSmoothingSteps=1,
227 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1,
229
243 template<class C>
244 KAMG(const Operator& fineOperator, const C& criterion,
245 const SmootherArgs& smootherArgs=SmootherArgs(),
246 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1,
248
250 void pre(Domain& x, Range& b);
252 void post(Domain& x);
254 void apply(Domain& v, const Range& d);
255
256 std::size_t maxlevels();
257
258 private:
260 Amg amg;
261
263 std::size_t maxLevelKrylovSteps;
264
266 double levelDefectReduction;
267
269 std::vector<std::shared_ptr<typename Amg::ScalarProduct> > scalarproducts;
270
272 std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > > ksolvers;
273 };
274
275 template<class M, class X, class S, class P, class K, class A>
276 KAMG<M,X,S,P,K,A>::KAMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
277 const SmootherArgs& smootherArgs,
278 std::size_t gamma, std::size_t preSmoothingSteps,
279 std::size_t postSmoothingSteps,
280 std::size_t ksteps, double reduction)
281 : amg(matrices, coarseSolver, smootherArgs, gamma, preSmoothingSteps,
282 postSmoothingSteps), maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
283 {}
284
285 template<class M, class X, class S, class P, class K, class A>
287 const SmootherArgs& smootherArgs, const Parameters& params,
288 std::size_t ksteps, double reduction)
289 : amg(matrices, coarseSolver, smootherArgs, params),
290 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
291 {}
292
293 template<class M, class X, class S, class P, class K, class A>
294 template<class C>
295 KAMG<M,X,S,P,K,A>::KAMG(const Operator& fineOperator, const C& criterion,
296 const SmootherArgs& smootherArgs, std::size_t gamma,
297 std::size_t preSmoothingSteps, std::size_t postSmoothingSteps,
298 std::size_t ksteps, double reduction,
299 const ParallelInformation& pinfo)
300 : amg(fineOperator, criterion, smootherArgs, gamma, preSmoothingSteps,
301 postSmoothingSteps, false, pinfo), maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
302 {}
303
304 template<class M, class X, class S, class P, class K, class A>
305 template<class C>
306 KAMG<M,X,S,P,K,A>::KAMG(const Operator& fineOperator, const C& criterion,
307 const SmootherArgs& smootherArgs,
308 std::size_t ksteps, double reduction,
309 const ParallelInformation& pinfo)
310 : amg(fineOperator, criterion, smootherArgs, pinfo),
311 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
312 {}
313
314
315 template<class M, class X, class S, class P, class K, class A>
317 {
318 amg.pre(x,b);
319 scalarproducts.reserve(amg.levels());
320 ksolvers.reserve(amg.levels());
321
322 typename OperatorHierarchy::ParallelMatrixHierarchy::Iterator
323 matrix = amg.matrices_->matrices().coarsest();
325 pinfo = amg.matrices_->parallelInformation().coarsest();
326 bool hasCoarsest=(amg.levels()==amg.maxlevels());
327
328 if(hasCoarsest) {
329 if(matrix==amg.matrices_->matrices().finest())
330 return;
331 --matrix;
332 --pinfo;
333 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, amg.solver_)));
334 }else
335 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, std::shared_ptr<InverseOperator<Domain,Range> >())));
336
337 std::ostringstream s;
338
339 if(matrix!=amg.matrices_->matrices().finest())
340 while(true) {
341 scalarproducts.push_back(std::shared_ptr<typename Amg::ScalarProduct>(Amg::ScalarProductChooser::construct(*pinfo)));
342 std::shared_ptr<InverseOperator<Domain,Range> > ks =
343 std::shared_ptr<InverseOperator<Domain,Range> >(new KrylovSolver(*matrix, *(scalarproducts.back()),
344 *(ksolvers.back()), levelDefectReduction,
345 maxLevelKrylovSteps, 0));
346 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, ks)));
347 --matrix;
348 --pinfo;
349 if(matrix==amg.matrices_->matrices().finest())
350 break;
351 }
352 }
353
354
355 template<class M, class X, class S, class P, class K, class A>
357 {
358 amg.post(x);
359
360 }
361
362 template<class M, class X, class S, class P, class K, class A>
364 {
365 if(ksolvers.size()==0)
366 {
367 Range td=d;
369 amg.solver_->apply(v,td,res);
370 }else
371 {
372 typedef typename Amg::LevelContext LevelContext;
373 std::shared_ptr<LevelContext> levelContext(new LevelContext);
374 amg.initIteratorsWithFineLevel(*levelContext);
375 typedef typename std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > >::iterator Iter;
376 for(Iter solver=ksolvers.begin(); solver!=ksolvers.end(); ++solver)
377 (*solver)->setLevelContext(levelContext);
378 ksolvers.back()->apply(v,d);
379 }
380 }
381
382 template<class M, class X, class S, class P, class K, class A>
383 std::size_t KAMG<M,X,S,P,K,A>::maxlevels()
384 {
385 return amg.maxlevels();
386 }
387
389 } // Amg
390} // Dune
391
392#endif
The AMG preconditioner.
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:56
LevelIterator< Hierarchy< ParallelInformation, Allocator >, ParallelInformation > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:257
an algebraic multigrid method using a Krylov-cycle.
Definition: kamg.hh:141
Amg::Domain Domain
the type of the domain.
Definition: kamg.hh:158
Amg::SmootherArgs SmootherArgs
The type of the arguments for construction of the smoothers.
Definition: kamg.hh:154
Amg::ParallelInformation ParallelInformation
the type of the parallelinformation to use.
Definition: kamg.hh:152
Amg::CoarseSolver CoarseSolver
The type of the coarse solver.
Definition: kamg.hh:150
Amg::OperatorHierarchy OperatorHierarchy
The type of the hierarchy of operators.
Definition: kamg.hh:148
Amg::Range Range
The type of the range.
Definition: kamg.hh:160
Amg::ScalarProduct ScalarProduct
The type of the scalar product.
Definition: kamg.hh:164
@ category
The solver category.
Definition: kamg.hh:168
AMG< M, X, S, PI, A > Amg
The type of the underlying AMG.
Definition: kamg.hh:144
Amg::Operator Operator
the type of the lineatr operator.
Definition: kamg.hh:156
Amg::ParallelInformationHierarchy ParallelInformationHierarchy
The type of the hierarchy of parallel information.
Definition: kamg.hh:162
K KrylovSolver
The type of the Krylov solver for the cycle.
Definition: kamg.hh:146
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:55
KAmgTwoGrid(AMG &amg, std::shared_ptr< InverseOperator< Domain, Range > > coarseSolver)
Constructor.
Definition: kamg.hh:50
~KAmgTwoGrid()
Destructor.
Definition: kamg.hh:111
void post(typename AMG::Domain &x)
Clean up.
Definition: kamg.hh:61
void setLevelContext(std::shared_ptr< typename AMG::LevelContext > p)
Set the level context pointer.
Definition: kamg.hh:105
InverseOperator< Domain, Range > * coarseSolver()
Get a pointer to the coarse grid solver.
Definition: kamg.hh:96
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:67
@ category
The solver category.
Definition: kamg.hh:40
The hierarchies build by the coarsening process.
Definition: hierarchy.hh:317
All parameters for AMG.
Definition: parameters.hh:391
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1504
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:26
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:44
#define DUNE_DEPRECATED
Mark some entity as deprecated.
Definition: deprecated.hh:84
void apply(Domain &v, const Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: kamg.hh:363
X Domain
The domain type.
Definition: amg.hh:78
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:91
M Operator
The matrix operator type.
Definition: amg.hh:64
void post(Domain &x)
Clean up.
Definition: kamg.hh:356
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) DUNE_DEPRECATED
Construct a new amg with a specific coarse solver.
Definition: kamg.hh:276
X Range
The range type.
Definition: amg.hh:80
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:408
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:430
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: kamg.hh:316
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:71
@ category
The solver category.
Definition: amg.hh:95
Dune namespace.
Definition: alignment.hh:11
STL namespace.
Define general preconditioner interface.
Statistics about the application of an inverse operator.
Definition: solver.hh:32
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)