3 #ifndef DUNE_AMG_KAMG_HH
4 #define DUNE_AMG_KAMG_HH
30 :
public Preconditioner<typename AMG::Domain,typename AMG::Range>
66 *levelContext_->update=0;
67 *levelContext_->rhs = d;
68 *levelContext_->lhs = v;
70 presmooth(*levelContext_, amg_.preSteps_);
71 bool processFineLevel =
72 amg_.moveToCoarseLevel(*levelContext_);
74 if(processFineLevel) {
78 coarseSolver_->apply(x, b, res);
79 *levelContext_->update=x;
82 amg_.moveToFineLevel(*levelContext_, processFineLevel);
85 v=*levelContext_->update;
114 Dune::shared_ptr<InverseOperator<Domain,Range> > coarseSolver_;
116 Dune::shared_ptr<typename AMG::LevelContext> levelContext_;
131 template<
class M,
class X,
class S,
class PI=SequentialInformation,
178 std::size_t preSmoothingSteps =1, std::size_t postSmoothingSteps = 1,
179 std::size_t maxLevelKrylovSteps = 3 ,
double minDefectReduction =1e-1);
200 std::size_t preSmoothingSteps=1, std::size_t postSmoothingSteps=1,
201 std::size_t maxLevelKrylovSteps=3,
double minDefectReduction=1e-1,
218 std::size_t maxLevelKrylovSteps;
221 double levelDefectReduction;
224 std::vector<shared_ptr<typename Amg::ScalarProduct> > scalarproducts;
227 std::vector<shared_ptr<KAmgTwoGrid<Amg> > > ksolvers;
230 template<
class M,
class X,
class S,
class P,
class K,
class A>
233 std::size_t gamma, std::size_t preSmoothingSteps,
234 std::size_t postSmoothingSteps,
235 std::size_t ksteps,
double reduction)
236 : amg(matrices, coarseSolver, smootherArgs, gamma, preSmoothingSteps,
237 postSmoothingSteps), maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
240 template<
class M,
class X,
class S,
class P,
class K,
class A>
244 std::size_t preSmoothingSteps, std::size_t postSmoothingSteps,
245 std::size_t ksteps,
double reduction,
247 : amg(fineOperator, criterion, smootherArgs, gamma, preSmoothingSteps,
248 postSmoothingSteps, false, pinfo), maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
252 template<
class M,
class X,
class S,
class P,
class K,
class A>
256 scalarproducts.reserve(amg.levels());
257 ksolvers.reserve(amg.levels());
259 typename OperatorHierarchy::ParallelMatrixHierarchy::Iterator
260 matrix = amg.matrices_->matrices().coarsest();
262 pinfo = amg.matrices_->parallelInformation().coarsest();
263 bool hasCoarsest=(amg.levels()==amg.maxlevels());
266 if(matrix==amg.matrices_->matrices().finest())
274 std::ostringstream s;
276 if(matrix!=amg.matrices_->matrices().finest())
278 scalarproducts.push_back(shared_ptr<typename Amg::ScalarProduct>(Amg::ScalarProductChooser::construct(*pinfo)));
279 shared_ptr<InverseOperator<Domain,Range> > ks =
280 shared_ptr<InverseOperator<Domain,Range> >(
new KrylovSolver(*matrix, *(scalarproducts.back()),
281 *(ksolvers.back()), levelDefectReduction,
282 maxLevelKrylovSteps, 0));
286 if(matrix==amg.matrices_->matrices().finest())
292 template<
class M,
class X,
class S,
class P,
class K,
class A>
299 template<
class M,
class X,
class S,
class P,
class K,
class A>
302 if(ksolvers.size()==0)
306 amg.solver_->apply(v,td,res);
309 typedef typename Amg::LevelContext LevelContext;
310 Dune::shared_ptr<LevelContext> levelContext(
new LevelContext);
311 amg.initIteratorsWithFineLevel(*levelContext);
312 typedef typename std::vector<shared_ptr<KAmgTwoGrid<Amg> > >::iterator Iter;
313 for(Iter solver=ksolvers.begin(); solver!=ksolvers.end(); ++solver)
314 (*solver)->setLevelContext(levelContext);
315 ksolvers.back()->apply(v,d);
319 template<
class M,
class X,
class S,
class P,
class K,
class A>
322 return amg.maxlevels();
The solver category.
Definition: amg.hh:96
Amg::ScalarProduct ScalarProduct
The type of the scalar product.
Definition: kamg.hh:157
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: kamg.hh:253
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:92
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:26
The solver category.
Definition: kamg.hh:40
InverseOperator< Domain, Range > * coarseSolver()
Get a pointer to the coarse grid solver.
Definition: kamg.hh:92
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:72
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1364
Amg::ParallelInformation ParallelInformation
the type of the parallelinformation to use.
Definition: kamg.hh:145
void apply(Domain &v, const Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: kamg.hh:300
Amg::SmootherArgs SmootherArgs
The type of the arguments for construction of the smoothers.
Definition: kamg.hh:147
Abstract base class for all solvers.
Definition: solver.hh:79
Amg::Range Range
The type of the range.
Definition: kamg.hh:153
Amg::Domain Domain
the type of the domain.
Definition: kamg.hh:151
Two grid operator for AMG with Krylov cycle.
Definition: amg.hh:45
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:63
void post(Domain &x)
Clean up.
Definition: kamg.hh:293
K KrylovSolver
The type of the Krylov solver for the cycle.
Definition: kamg.hh:139
~KAmgTwoGrid()
Destructor.
Definition: kamg.hh:107
Amg::ParallelInformationHierarchy ParallelInformationHierarchy
The type of the hierarchy of parallel information.
Definition: kamg.hh:155
AMG< M, X, S, PI, A > Amg
The type of the underlying AMG.
Definition: kamg.hh:137
void pre(typename AMG::Domain &x, typename AMG::Range &b)
Prepare the preconditioner.
Definition: kamg.hh:55
void setLevelContext(Dune::shared_ptr< typename AMG::LevelContext > p)
Set the level context pointer.
Definition: kamg.hh:101
The hierarchies build by the coarsening process.
Definition: hierarchy.hh:317
M Operator
The matrix operator type.
Definition: amg.hh:65
Matrix & A
Definition: matrixmatrix.hh:216
void post(typename AMG::Domain &x)
Clean up.
Definition: kamg.hh:59
std::size_t maxlevels()
Definition: kamg.hh:320
Amg::Operator Operator
the type of the lineatr operator.
Definition: kamg.hh:149
Statistics about the application of an inverse operator.
Definition: solver.hh:31
X Domain
The domain type.
Definition: amg.hh:79
KAmgTwoGrid(AMG &amg, shared_ptr< InverseOperator< Domain, Range > > coarseSolver)
Constructor.
Definition: kamg.hh:50
LevelIterator< Hierarchy< ParallelInformation, Allocator >, ParallelInformation > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:258
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:231
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:56
X Range
The range type.
Definition: amg.hh:81
Amg::CoarseSolver CoarseSolver
The type of the coarse solver.
Definition: kamg.hh:143
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
The solver category.
Definition: kamg.hh:161
Amg::OperatorHierarchy OperatorHierarchy
The type of the hierarchy of operators.
Definition: kamg.hh:141
Define general preconditioner interface.