3 #ifndef DUNE_AMG_KAMG_HH 4 #define DUNE_AMG_KAMG_HH 30 :
public Preconditioner<typename AMG::Domain,typename AMG::Range>
57 DUNE_UNUSED_PARAMETER(x); DUNE_UNUSED_PARAMETER(b);
63 DUNE_UNUSED_PARAMETER(x);
70 *levelContext_->update=0;
71 *levelContext_->rhs = d;
72 *levelContext_->lhs = v;
74 presmooth(*levelContext_, amg_.preSteps_);
75 bool processFineLevel =
76 amg_.moveToCoarseLevel(*levelContext_);
78 if(processFineLevel) {
82 coarseSolver_->apply(x, b, res);
83 *levelContext_->update=x;
86 amg_.moveToFineLevel(*levelContext_, processFineLevel);
89 v=*levelContext_->update;
118 std::shared_ptr<InverseOperator<Domain,Range> > coarseSolver_;
120 std::shared_ptr<typename AMG::LevelContext> levelContext_;
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;
199 KAMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
200 const SmootherArgs& smootherArgs, const
Parameters& parms,
201 std::
size_t maxLevelKrylovSteps=3,
double minDefectReduction=1e-1);
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,
228 const ParallelInformation& pinfo=ParallelInformation()) DUNE_DEPRECATED;
244 KAMG(const Operator& fineOperator, const C& criterion,
245 const SmootherArgs& smootherArgs=SmootherArgs(),
246 std::
size_t maxLevelKrylovSteps=3,
double minDefectReduction=1e-1,
247 const ParallelInformation& pinfo=ParallelInformation());
250 void pre(Domain& x, Range& b);
252 void post(Domain& x);
254 void apply(Domain& v, const Range& d);
256 std::
size_t maxlevels();
263 std::
size_t maxLevelKrylovSteps;
266 double levelDefectReduction;
269 std::vector<
std::shared_ptr<typename Amg::ScalarProduct> > scalarproducts;
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)
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)
293 template<
class M,
class X,
class S,
class P,
class K,
class A>
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)
304 template<
class M,
class X,
class S,
class P,
class K,
class A>
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)
315 template<
class M,
class X,
class S,
class P,
class K,
class A>
319 scalarproducts.reserve(amg.
levels());
320 ksolvers.reserve(amg.
levels());
323 matrix = amg.matrices_->matrices().coarsest();
325 pinfo = amg.matrices_->parallelInformation().coarsest();
329 if(matrix==amg.matrices_->matrices().finest())
337 std::ostringstream s;
339 if(matrix!=amg.matrices_->matrices().finest())
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));
349 if(matrix==amg.matrices_->matrices().finest())
355 template<
class M,
class X,
class S,
class P,
class K,
class A>
362 template<
class M,
class X,
class S,
class P,
class K,
class A>
365 if(ksolvers.size()==0)
369 amg.solver_->apply(v,td,res);
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);
382 template<
class M,
class X,
class S,
class P,
class K,
class A>
void apply(Domain &v, const Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: kamg.hh:363
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:430
The hierarchies build by the coarsening process.
Definition: hierarchy.hh:316
All parameters for AMG.
Definition: parameters.hh:390
AMG< M, X, S, PI, A > Amg
The type of the underlying AMG.
Definition: kamg.hh:144
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:276
X Range
The range type.
Definition: amg.hh:80
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
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:91
Statistics about the application of an inverse operator.
Definition: solver.hh:31
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:26
Amg::ScalarProduct ScalarProduct
The type of the scalar product.
Definition: kamg.hh:164
The solver category.
Definition: kamg.hh:40
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:55
void post(Domain &x)
Clean up.
Definition: kamg.hh:356
Define general preconditioner interface.
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: kamg.hh:316
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: amg.hh:455
std::size_t maxlevels()
Definition: kamg.hh:383
void post(Domain &x)
Clean up.
Definition: amg.hh:792
Amg::Range Range
The type of the range.
Definition: kamg.hh:160
Two grid operator for AMG with Krylov cycle.
Definition: amg.hh:44
Amg::ParallelInformation ParallelInformation
the type of the parallelinformation to use.
Definition: kamg.hh:152
Amg::Domain Domain
the type of the domain.
Definition: kamg.hh:158
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1503
InverseOperator< Domain, Range > * coarseSolver()
Get a pointer to the coarse grid solver.
Definition: kamg.hh:96
K KrylovSolver
The type of the Krylov solver for the cycle.
Definition: kamg.hh:146
void post(typename AMG::Domain &x)
Clean up.
Definition: kamg.hh:61
Definition: basearray.hh:19
Amg::ParallelInformationHierarchy ParallelInformationHierarchy
The type of the hierarchy of parallel information.
Definition: kamg.hh:162
KAmgTwoGrid(AMG &amg, std::shared_ptr< InverseOperator< Domain, Range > > coarseSolver)
Constructor.
Definition: kamg.hh:50
The solver category.
Definition: amg.hh:95
Amg::OperatorHierarchy OperatorHierarchy
The type of the hierarchy of operators.
Definition: kamg.hh:148
std::size_t maxlevels()
Definition: amg.hh:540
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:408
LevelIterator< Hierarchy< MatrixOperator, Allocator >, MatrixOperator > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:257
X Domain
The domain type.
Definition: amg.hh:78
M Operator
The matrix operator type.
Definition: amg.hh:64
Amg::SmootherArgs SmootherArgs
The type of the arguments for construction of the smoothers.
Definition: kamg.hh:154
Abstract base class for all solvers.
Definition: solver.hh:79
Amg::CoarseSolver CoarseSolver
The type of the coarse solver.
Definition: kamg.hh:150
~KAmgTwoGrid()
Destructor.
Definition: kamg.hh:111
void setLevelContext(std::shared_ptr< typename AMG::LevelContext > p)
Set the level context pointer.
Definition: kamg.hh:105
Amg::Operator Operator
the type of the lineatr operator.
Definition: kamg.hh:156
an algebraic multigrid method using a Krylov-cycle.
Definition: amg.hh:41
void pre(typename AMG::Domain &x, typename AMG::Range &b)
Prepare the preconditioner.
Definition: kamg.hh:55
std::size_t levels()
Definition: amg.hh:535
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:71