3#ifndef DUNE_AMG_KAMG_HH
4#define DUNE_AMG_KAMG_HH
67 *levelContext_->update=0;
68 *levelContext_->rhs = d;
69 *levelContext_->lhs = v;
71 presmooth(*levelContext_, amg_.preSteps_);
72 bool processFineLevel =
73 amg_.moveToCoarseLevel(*levelContext_);
75 if(processFineLevel) {
79 coarseSolver_->apply(x, b, res);
80 *levelContext_->update=x;
83 amg_.moveToFineLevel(*levelContext_, processFineLevel);
86 v=*levelContext_->update;
115 std::shared_ptr<InverseOperator<Domain,Range> > coarseSolver_;
117 std::shared_ptr<typename AMG::LevelContext> levelContext_;
135 template<
class M,
class X,
class S,
class PI=SequentialInformation,
182 std::size_t maxLevelKrylovSteps=3,
double minDefectReduction=1e-1);
200 std::size_t maxLevelKrylovSteps=3,
double minDefectReduction=1e-1,
217 std::size_t maxLevelKrylovSteps;
220 double levelDefectReduction;
223 std::vector<std::shared_ptr<typename Amg::ScalarProduct> > scalarproducts;
226 std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > > ksolvers;
230 template<
class M,
class X,
class S,
class P,
class K,
class A>
233 std::size_t ksteps,
double reduction)
234 : amg(matrices, coarseSolver, smootherArgs, params),
235 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
239 template<
class M,
class X,
class S,
class P,
class K,
class A>
243 std::size_t ksteps,
double reduction,
245 : amg(fineOperator, criterion, smootherArgs, pinfo),
246 maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
250 template<
class M,
class X,
class S,
class P,
class K,
class A>
254 scalarproducts.reserve(amg.levels());
255 ksolvers.reserve(amg.levels());
257 typename OperatorHierarchy::ParallelMatrixHierarchy::Iterator
258 matrix = amg.matrices_->matrices().coarsest();
260 pinfo = amg.matrices_->parallelInformation().coarsest();
261 bool hasCoarsest=(amg.levels()==amg.maxlevels());
264 if(matrix==amg.matrices_->matrices().finest())
272 std::ostringstream s;
274 if(matrix!=amg.matrices_->matrices().finest())
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));
284 if(matrix==amg.matrices_->matrices().finest())
290 template<
class M,
class X,
class S,
class P,
class K,
class A>
297 template<
class M,
class X,
class S,
class P,
class K,
class A>
300 if(ksolvers.size()==0)
304 amg.solver_->apply(v,td,res);
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);
317 template<
class M,
class X,
class S,
class P,
class K,
class A>
320 return amg.maxlevels();
Define general preconditioner interface.
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
std::size_t maxlevels()
Definition: kamg.hh:318
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
Definition: allocator.hh:9
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
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
The hierarchies build by the coarsening process.
Definition: matrixhierarchy.hh:59
All parameters for AMG.
Definition: parameters.hh:391
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:30
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:50
Statistics about the application of an inverse operator.
Definition: solver.hh:46
Abstract base class for all solvers.
Definition: solver.hh:97
Category
Definition: solvercategory.hh:21
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1284