DUNE PDELab (2.8)

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([[maybe_unused]] typename AMG::Domain& x, [[maybe_unused]] typename AMG::Range& b)
57 {}
58
60 void post([[maybe_unused]] typename AMG::Domain& x)
61 {}
62
64 void apply(typename AMG::Domain& v, const typename AMG::Range& d)
65 {
66 // Copy data
67 *levelContext_->update=0;
68 *levelContext_->rhs = d;
69 *levelContext_->lhs = v;
70
71 presmooth(*levelContext_, amg_.preSteps_);
72 bool processFineLevel =
73 amg_.moveToCoarseLevel(*levelContext_);
74
75 if(processFineLevel) {
76 typename AMG::Range b=*levelContext_->rhs;
77 typename AMG::Domain x=*levelContext_->update;
79 coarseSolver_->apply(x, b, res);
80 *levelContext_->update=x;
81 }
82
83 amg_.moveToFineLevel(*levelContext_, processFineLevel);
84
85 postsmooth(*levelContext_, amg_.postSteps_);
86 v=*levelContext_->update;
87 }
88
94 {
95 return coarseSolver_;
96 }
97
102 void setLevelContext(std::shared_ptr<typename AMG::LevelContext> p)
103 {
104 levelContext_=p;
105 }
106
109 {}
110
111 private:
113 AMG& amg_;
115 std::shared_ptr<InverseOperator<Domain,Range> > coarseSolver_;
117 std::shared_ptr<typename AMG::LevelContext> levelContext_;
118 };
119
120
121
135 template<class M, class X, class S, class PI=SequentialInformation,
136 class K=GeneralizedPCGSolver<X>, class A=std::allocator<X> >
137 class KAMG : public Preconditioner<X,X>
138 {
139 public:
143 typedef K KrylovSolver;
153 typedef typename Amg::Operator Operator;
155 typedef typename Amg::Domain Domain;
157 typedef typename Amg::Range Range;
162
165 {
166 return amg.category();
167 };
168
180 KAMG(OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
181 const SmootherArgs& smootherArgs, const Parameters& parms,
182 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1);
183
197 template<class C>
198 KAMG(const Operator& fineOperator, const C& criterion,
199 const SmootherArgs& smootherArgs=SmootherArgs(),
200 std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1,
202
204 void pre(Domain& x, Range& b);
206 void post(Domain& x);
208 void apply(Domain& v, const Range& d);
209
210 std::size_t maxlevels();
211
212 private:
214 Amg amg;
215
217 std::size_t maxLevelKrylovSteps;
218
220 double levelDefectReduction;
221
223 std::vector<std::shared_ptr<typename Amg::ScalarProduct> > scalarproducts;
224
226 std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > > ksolvers;
227 };
228
229
230 template<class M, class X, class S, class P, class K, class A>
232 const SmootherArgs& smootherArgs, const Parameters& params,
233 std::size_t ksteps, double reduction)
234 : amg(matrices, coarseSolver, smootherArgs, params),
235 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
236 {}
237
238
239 template<class M, class X, class S, class P, class K, class A>
240 template<class C>
241 KAMG<M,X,S,P,K,A>::KAMG(const Operator& fineOperator, const C& criterion,
242 const SmootherArgs& smootherArgs,
243 std::size_t ksteps, double reduction,
244 const ParallelInformation& pinfo)
245 : amg(fineOperator, criterion, smootherArgs, pinfo),
246 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
247 {}
248
249
250 template<class M, class X, class S, class P, class K, class A>
252 {
253 amg.pre(x,b);
254 scalarproducts.reserve(amg.levels());
255 ksolvers.reserve(amg.levels());
256
257 typename OperatorHierarchy::ParallelMatrixHierarchy::Iterator
258 matrix = amg.matrices_->matrices().coarsest();
260 pinfo = amg.matrices_->parallelInformation().coarsest();
261 bool hasCoarsest=(amg.levels()==amg.maxlevels());
262
263 if(hasCoarsest) {
264 if(matrix==amg.matrices_->matrices().finest())
265 return;
266 --matrix;
267 --pinfo;
268 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, amg.solver_)));
269 }else
270 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, std::shared_ptr<InverseOperator<Domain,Range> >())));
271
272 std::ostringstream s;
273
274 if(matrix!=amg.matrices_->matrices().finest())
275 while(true) {
276 scalarproducts.push_back(createScalarProduct<X>(*pinfo,category()));
277 std::shared_ptr<InverseOperator<Domain,Range> > ks =
278 std::shared_ptr<InverseOperator<Domain,Range> >(new KrylovSolver(*matrix, *(scalarproducts.back()),
279 *(ksolvers.back()), levelDefectReduction,
280 maxLevelKrylovSteps, 0));
281 ksolvers.push_back(std::shared_ptr<KAmgTwoGrid<Amg> >(new KAmgTwoGrid<Amg>(amg, ks)));
282 --matrix;
283 --pinfo;
284 if(matrix==amg.matrices_->matrices().finest())
285 break;
286 }
287 }
288
289
290 template<class M, class X, class S, class P, class K, class A>
292 {
293 amg.post(x);
294
295 }
296
297 template<class M, class X, class S, class P, class K, class A>
299 {
300 if(ksolvers.size()==0)
301 {
302 Range td=d;
304 amg.solver_->apply(v,td,res);
305 }else
306 {
307 typedef typename Amg::LevelContext LevelContext;
308 std::shared_ptr<LevelContext> levelContext(new LevelContext);
309 amg.initIteratorsWithFineLevel(*levelContext);
310 typedef typename std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > >::iterator Iter;
311 for(Iter solver=ksolvers.begin(); solver!=ksolvers.end(); ++solver)
312 (*solver)->setLevelContext(levelContext);
313 ksolvers.back()->apply(v,d);
314 }
315 }
316
317 template<class M, class X, class S, class P, class K, class A>
318 std::size_t KAMG<M,X,S,P,K,A>::maxlevels()
319 {
320 return amg.maxlevels();
321 }
322
324 } // Amg
325} // Dune
326
327#endif
The AMG preconditioner.
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:63
LevelIterator< Hierarchy< ParallelInformation, Allocator >, ParallelInformation > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:214
an algebraic multigrid method using a Krylov-cycle.
Definition: kamg.hh:138
Amg::Domain Domain
the type of the domain.
Definition: kamg.hh:155
Amg::SmootherArgs SmootherArgs
The type of the arguments for construction of the smoothers.
Definition: kamg.hh:151
Amg::ParallelInformation ParallelInformation
the type of the parallelinformation to use.
Definition: kamg.hh:149
Amg::CoarseSolver CoarseSolver
The type of the coarse solver.
Definition: kamg.hh:147
Amg::OperatorHierarchy OperatorHierarchy
The type of the hierarchy of operators.
Definition: kamg.hh:145
Amg::Range Range
The type of the range.
Definition: kamg.hh:157
Amg::ScalarProduct ScalarProduct
The type of the scalar product.
Definition: kamg.hh:161
AMG< M, X, S, PI, A > Amg
The type of the underlying AMG.
Definition: kamg.hh:141
Amg::Operator Operator
the type of the lineatr operator.
Definition: kamg.hh:153
Amg::ParallelInformationHierarchy ParallelInformationHierarchy
The type of the hierarchy of parallel information.
Definition: kamg.hh:159
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: kamg.hh:164
K KrylovSolver
The type of the Krylov solver for the cycle.
Definition: kamg.hh:143
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:108
void post(typename AMG::Domain &x)
Clean up.
Definition: kamg.hh:60
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:102
InverseOperator< Domain, Range > * coarseSolver()
Get a pointer to the coarse grid solver.
Definition: kamg.hh:93
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:64
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:1284
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:30
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:50
void apply(Domain &v, const Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: kamg.hh:298
X Domain
The domain type.
Definition: amg.hh:85
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:231
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:98
M Operator
The matrix operator type.
Definition: amg.hh:71
void post(Domain &x)
Clean up.
Definition: kamg.hh:291
X Range
The range type.
Definition: amg.hh:87
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:404
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:426
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: kamg.hh:251
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: amg.hh:192
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:78
Dune namespace.
Definition: alignedallocator.hh:11
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 (Jul 27, 22:29, 2024)