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_;
138 template<
class M,
class X,
class S,
class PI=SequentialInformation,
185 std::size_t preSmoothingSteps =1, std::size_t postSmoothingSteps = 1,
186 std::size_t maxLevelKrylovSteps = 3 ,
double minDefectReduction =1e-1);
207 std::size_t preSmoothingSteps=1, std::size_t postSmoothingSteps=1,
208 std::size_t maxLevelKrylovSteps=3,
double minDefectReduction=1e-1,
225 std::size_t maxLevelKrylovSteps;
228 double levelDefectReduction;
231 std::vector<std::shared_ptr<typename Amg::ScalarProduct> > scalarproducts;
234 std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > > ksolvers;
237 template<
class M,
class X,
class S,
class P,
class K,
class A>
240 std::size_t gamma, std::size_t preSmoothingSteps,
241 std::size_t postSmoothingSteps,
242 std::size_t ksteps,
double reduction)
243 : amg(matrices, coarseSolver, smootherArgs, gamma, preSmoothingSteps,
244 postSmoothingSteps), maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
247 template<
class M,
class X,
class S,
class P,
class K,
class A>
251 std::size_t preSmoothingSteps, std::size_t postSmoothingSteps,
252 std::size_t ksteps,
double reduction,
254 : amg(fineOperator, criterion, smootherArgs, gamma, preSmoothingSteps,
255 postSmoothingSteps, false, pinfo), maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
259 template<
class M,
class X,
class S,
class P,
class K,
class A>
263 scalarproducts.reserve(amg.levels());
264 ksolvers.reserve(amg.levels());
266 typename OperatorHierarchy::ParallelMatrixHierarchy::Iterator
267 matrix = amg.matrices_->matrices().coarsest();
269 pinfo = amg.matrices_->parallelInformation().coarsest();
270 bool hasCoarsest=(amg.levels()==amg.maxlevels());
273 if(matrix==amg.matrices_->matrices().finest())
281 std::ostringstream s;
283 if(matrix!=amg.matrices_->matrices().finest())
285 scalarproducts.push_back(std::shared_ptr<typename Amg::ScalarProduct>(Amg::ScalarProductChooser::construct(*pinfo)));
286 std::shared_ptr<InverseOperator<Domain,Range> > ks =
287 std::shared_ptr<InverseOperator<Domain,Range> >(
new KrylovSolver(*matrix, *(scalarproducts.back()),
288 *(ksolvers.back()), levelDefectReduction,
289 maxLevelKrylovSteps, 0));
293 if(matrix==amg.matrices_->matrices().finest())
299 template<
class M,
class X,
class S,
class P,
class K,
class A>
306 template<
class M,
class X,
class S,
class P,
class K,
class A>
309 if(ksolvers.size()==0)
313 amg.solver_->apply(v,td,res);
316 typedef typename Amg::LevelContext LevelContext;
317 std::shared_ptr<LevelContext> levelContext(
new LevelContext);
318 amg.initIteratorsWithFineLevel(*levelContext);
319 typedef typename std::vector<std::shared_ptr<KAmgTwoGrid<Amg> > >::iterator Iter;
320 for(Iter solver=ksolvers.begin(); solver!=ksolvers.end(); ++solver)
321 (*solver)->setLevelContext(levelContext);
322 ksolvers.back()->apply(v,d);
326 template<
class M,
class X,
class S,
class P,
class K,
class A>
329 return amg.maxlevels();
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:71
Amg::OperatorHierarchy OperatorHierarchy
The type of the hierarchy of operators.
Definition: kamg.hh:148
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:430
void apply(Domain &v, const Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: kamg.hh:307
KAmgTwoGrid(AMG &amg, std::shared_ptr< InverseOperator< Domain, Range > > coarseSolver)
Constructor.
Definition: kamg.hh:50
void post(typename AMG::Domain &x)
Clean up.
Definition: kamg.hh:61
InverseOperator< Domain, Range > * coarseSolver()
Get a pointer to the coarse grid solver.
Definition: kamg.hh:96
Amg::Operator Operator
the type of the lineatr operator.
Definition: kamg.hh:156
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
Two grid operator for AMG with Krylov cycle.
Definition: amg.hh:44
~KAmgTwoGrid()
Destructor.
Definition: kamg.hh:111
Matrix & A
Definition: matrixmatrix.hh:216
Amg::SmootherArgs SmootherArgs
The type of the arguments for construction of the smoothers.
Definition: kamg.hh:154
X Range
The range type.
Definition: amg.hh:80
The solver category.
Definition: kamg.hh:168
void post(Domain &x)
Clean up.
Definition: kamg.hh:300
void setLevelContext(std::shared_ptr< typename AMG::LevelContext > p)
Set the level context pointer.
Definition: kamg.hh:105
The solver category.
Definition: amg.hh:95
K KrylovSolver
The type of the Krylov solver for the cycle.
Definition: kamg.hh:146
std::size_t maxlevels()
Definition: kamg.hh:327
X Domain
The domain type.
Definition: amg.hh:78
Amg::CoarseSolver CoarseSolver
The type of the coarse solver.
Definition: kamg.hh:150
Abstract base class for all solvers.
Definition: solver.hh:79
void pre(typename AMG::Domain &x, typename AMG::Range &b)
Prepare the preconditioner.
Definition: kamg.hh:55
M Operator
The matrix operator type.
Definition: amg.hh:64
Amg::ParallelInformation ParallelInformation
the type of the parallelinformation to use.
Definition: kamg.hh:152
Amg::ParallelInformationHierarchy ParallelInformationHierarchy
The type of the hierarchy of parallel information.
Definition: kamg.hh:162
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: kamg.hh:260
Define general preconditioner interface.
The hierarchies build by the coarsening process.
Definition: hierarchy.hh:316
Amg::Range Range
The type of the range.
Definition: kamg.hh:160
Amg::Domain Domain
the type of the domain.
Definition: kamg.hh:158
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:55
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:238
AMG< M, X, S, PI, A > Amg
The type of the underlying AMG.
Definition: kamg.hh:144
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:91
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
Generalized preconditioned conjugate gradient solver.
Definition: solvers.hh:1455
Statistics about the application of an inverse operator.
Definition: solver.hh:31
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:408
The solver category.
Definition: kamg.hh:40
LevelIterator< Hierarchy< ParallelInformation, Allocator >, ParallelInformation > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:257