Dune Core Modules (2.7.1)

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
184 KAMG(OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
185 const SmootherArgs& smootherArgs, const Parameters& parms,
186 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1);
187
201 template<class C>
202 KAMG(const Operator& fineOperator, const C& criterion,
203 const SmootherArgs& smootherArgs=SmootherArgs(),
204 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1,
206
208 void pre(Domain& x, Range& b);
210 void post(Domain& x);
212 void apply(Domain& v, const Range& d);
213
214 std::size_t maxlevels();
215
216 private:
218 Amg amg;
219
221 std::size_t maxLevelKrylovSteps;
222
224 double levelDefectReduction;
225
227 std::vector<std::shared_ptr<typename Amg::ScalarProduct> > scalarproducts;
228
230 std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > > ksolvers;
231 };
232
233
234 template<class M, class X, class S, class P, class K, class A>
236 const SmootherArgs& smootherArgs, const Parameters& params,
237 std::size_t ksteps, double reduction)
238 : amg(matrices, coarseSolver, smootherArgs, params),
239 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
240 {}
241
242
243 template<class M, class X, class S, class P, class K, class A>
244 template<class C>
245 KAMG<M,X,S,P,K,A>::KAMG(const Operator& fineOperator, const C& criterion,
246 const SmootherArgs& smootherArgs,
247 std::size_t ksteps, double reduction,
248 const ParallelInformation& pinfo)
249 : amg(fineOperator, criterion, smootherArgs, pinfo),
250 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
251 {}
252
253
254 template<class M, class X, class S, class P, class K, class A>
256 {
257 amg.pre(x,b);
258 scalarproducts.reserve(amg.levels());
259 ksolvers.reserve(amg.levels());
260
261 typename OperatorHierarchy::ParallelMatrixHierarchy::Iterator
262 matrix = amg.matrices_->matrices().coarsest();
264 pinfo = amg.matrices_->parallelInformation().coarsest();
265 bool hasCoarsest=(amg.levels()==amg.maxlevels());
266
267 if(hasCoarsest) {
268 if(matrix==amg.matrices_->matrices().finest())
269 return;
270 --matrix;
271 --pinfo;
272 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, amg.solver_)));
273 }else
274 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, std::shared_ptr<InverseOperator<Domain,Range> >())));
275
276 std::ostringstream s;
277
278 if(matrix!=amg.matrices_->matrices().finest())
279 while(true) {
280 scalarproducts.push_back(createScalarProduct<X>(*pinfo,category()));
281 std::shared_ptr<InverseOperator<Domain,Range> > ks =
282 std::shared_ptr<InverseOperator<Domain,Range> >(new KrylovSolver(*matrix, *(scalarproducts.back()),
283 *(ksolvers.back()), levelDefectReduction,
284 maxLevelKrylovSteps, 0));
285 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, ks)));
286 --matrix;
287 --pinfo;
288 if(matrix==amg.matrices_->matrices().finest())
289 break;
290 }
291 }
292
293
294 template<class M, class X, class S, class P, class K, class A>
296 {
297 amg.post(x);
298
299 }
300
301 template<class M, class X, class S, class P, class K, class A>
303 {
304 if(ksolvers.size()==0)
305 {
306 Range td=d;
308 amg.solver_->apply(v,td,res);
309 }else
310 {
311 typedef typename Amg::LevelContext LevelContext;
312 std::shared_ptr<LevelContext> levelContext(new LevelContext);
313 amg.initIteratorsWithFineLevel(*levelContext);
314 typedef typename std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > >::iterator Iter;
315 for(Iter solver=ksolvers.begin(); solver!=ksolvers.end(); ++solver)
316 (*solver)->setLevelContext(levelContext);
317 ksolvers.back()->apply(v,d);
318 }
319 }
320
321 template<class M, class X, class S, class P, class K, class A>
322 std::size_t KAMG<M,X,S,P,K,A>::maxlevels()
323 {
324 return amg.maxlevels();
325 }
326
328 } // Amg
329} // Dune
330
331#endif
The AMG preconditioner.
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:62
LevelIterator< Hierarchy< ParallelInformation, Allocator >, ParallelInformation > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:215
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: matrixhierarchy.hh:59
All parameters for AMG.
Definition: parameters.hh:391
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1285
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:30
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:48
#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:302
X Domain
The domain type.
Definition: amg.hh:84
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:235
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:97
M Operator
The matrix operator type.
Definition: amg.hh:70
void post(Domain &x)
Clean up.
Definition: kamg.hh:295
X Range
The range type.
Definition: amg.hh:86
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:468
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:490
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: kamg.hh:255
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: amg.hh:165
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:77
Dune namespace.
Definition: alignedallocator.hh:14
Define general preconditioner interface.
Statistics about the application of an inverse operator.
Definition: solver.hh:46
Category
Definition: solvercategory.hh:21
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)