dune-istl  2.4.1-rc2
twolevelmethod.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_ISTL_TWOLEVELMETHOD_HH
4 #define DUNE_ISTL_TWOLEVELMETHOD_HH
5 
7 #include"amg.hh"
8 #include"galerkin.hh"
9 #include<dune/istl/solver.hh>
10 
11 #include<dune/common/unused.hh>
12 
20 namespace Dune
21 {
22 namespace Amg
23 {
24 
34 template<class FO, class CO>
36 {
37 public:
42  typedef FO FineOperatorType;
46  typedef typename FineOperatorType::range_type FineRangeType;
50  typedef typename FineOperatorType::domain_type FineDomainType;
55  typedef CO CoarseOperatorType;
59  typedef typename CoarseOperatorType::range_type CoarseRangeType;
63  typedef typename CoarseOperatorType::domain_type CoarseDomainType;
68  std::shared_ptr<CoarseOperatorType>& getCoarseLevelOperator()
69  {
70  return operator_;
71  }
77  {
78  return rhs_;
79  }
80 
86  {
87  return lhs_;
88  }
98  virtual void moveToCoarseLevel(const FineRangeType& fineRhs)=0;
108  virtual void moveToFineLevel(FineDomainType& fineLhs)=0;
116  virtual void createCoarseLevelSystem(const FineOperatorType& fineOperator)=0;
117 
119  virtual LevelTransferPolicy* clone() const =0;
120 
123 
124  protected:
130  std::shared_ptr<CoarseOperatorType> operator_;
131 };
132 
138 template<class O, class C>
140  : public LevelTransferPolicy<O,O>
141 {
143 public:
145  typedef C Criterion;
147 
149  : criterion_(crit)
150  {}
151 
152  void createCoarseLevelSystem(const O& fineOperator)
153  {
154  prolongDamp_ = criterion_.getProlongationDampingFactor();
157  typedef typename Dune::Amg::PropertiesGraph<MatrixGraph,Dune::Amg::VertexProperties,
158  Dune::Amg::EdgeProperties,Dune::IdentityMap,Dune::IdentityMap> PropertiesGraph;
159  MatrixGraph mg(fineOperator.getmat());
160  PropertiesGraph pg(mg,Dune::IdentityMap(),Dune::IdentityMap());
161  typedef NegateSet<typename ParallelInformation::OwnerSet> OverlapFlags;
162 
163  aggregatesMap_.reset(new AggregatesMap(pg.maxVertex()+1));
164 
165  int noAggregates, isoAggregates, oneAggregates, skippedAggregates;
166 
167  tie(noAggregates, isoAggregates, oneAggregates, skippedAggregates) =
168  aggregatesMap_->buildAggregates(fineOperator.getmat(), pg, criterion_, true);
169  std::cout<<"no aggregates="<<noAggregates<<" iso="<<isoAggregates<<" one="<<oneAggregates<<" skipped="<<skippedAggregates<<std::endl;
170  // misuse coarsener to renumber aggregates
172  typedef std::vector<bool>::iterator Iterator;
173  typedef Dune::IteratorPropertyMap<Iterator, Dune::IdentityMap> VisitedMap;
174  std::vector<bool> excluded(fineOperator.getmat().N(), false);
175  VisitedMap vm(excluded.begin(), Dune::IdentityMap());
176  ParallelInformation pinfo;
177  std::size_t aggregates = renumberer.coarsen(pinfo, pg, vm,
178  *aggregatesMap_, pinfo,
179  noAggregates);
180  std::vector<bool>& visited=excluded;
181 
182  typedef std::vector<bool>::iterator Iterator;
183 
184  for(Iterator iter= visited.begin(), end=visited.end();
185  iter != end; ++iter)
186  *iter=false;
187  matrix_.reset(productBuilder.build(mg, vm,
189  *aggregatesMap_,
190  aggregates,
191  OverlapFlags()));
192  productBuilder.calculate(fineOperator.getmat(), *aggregatesMap_, *matrix_, pinfo, OverlapFlags());
193  this->lhs_.resize(this->matrix_->M());
194  this->rhs_.resize(this->matrix_->N());
195  this->operator_.reset(new O(*matrix_));
196  }
197 
198  void moveToCoarseLevel(const typename FatherType::FineRangeType& fineRhs)
199  {
201  ::restrictVector(*aggregatesMap_, this->rhs_, fineRhs, ParallelInformation());
202  this->lhs_=0;
203  }
204 
206  {
208  ::prolongateVector(*aggregatesMap_, this->lhs_, fineLhs,
209  prolongDamp_, ParallelInformation());
210  }
211 
213  {
214  return new AggregationLevelTransferPolicy(*this);
215  }
216 
217 private:
218  typename O::matrix_type::field_type prolongDamp_;
219  std::shared_ptr<AggregatesMap> aggregatesMap_;
220  Criterion criterion_;
221  std::shared_ptr<typename O::matrix_type> matrix_;
222 };
223 
230 template<class O, class S, class C>
232 {
233 public:
235  typedef O Operator;
237  typedef typename O::range_type X;
239  typedef C Criterion;
241  typedef S Smoother;
252  : smootherArgs_(args), criterion_(c)
253  {}
256  : coarseOperator_(other.coarseOperator_), smootherArgs_(other.smootherArgs_),
257  criterion_(other.criterion_)
258  {}
259 private:
266  struct AMGInverseOperator : public InverseOperator<X,X>
267  {
268  AMGInverseOperator(const typename AMGType::Operator& op,
269  const Criterion& crit,
270  const typename AMGType::SmootherArgs& args)
271  : amg_(op, crit,args), first_(true)
272  {}
273 
274  void apply(X& x, X& b, double reduction, InverseOperatorResult& res)
275  {
276  DUNE_UNUSED_PARAMETER(reduction);
277  DUNE_UNUSED_PARAMETER(res);
278  if(first_)
279  {
280  amg_.pre(x,b);
281  first_=false;
282  x_=x;
283  }
284  amg_.apply(x,b);
285  }
286 
287  void apply(X& x, X& b, InverseOperatorResult& res)
288  {
289  return apply(x,b,1e-8,res);
290  }
291 
292  ~AMGInverseOperator()
293  {
294  if(!first_)
295  amg_.post(x_);
296  }
297  AMGInverseOperator(const AMGInverseOperator& other)
298  : x_(other.x_), amg_(other.amg_), first_(other.first_)
299  {
300  }
301  private:
302  X x_;
303  AMGType amg_;
304  bool first_;
305  };
306 
307 public:
309  typedef AMGInverseOperator CoarseLevelSolver;
310 
318  template<class P>
320  {
321  coarseOperator_=transferPolicy.getCoarseLevelOperator();
322  AMGInverseOperator* inv = new AMGInverseOperator(*coarseOperator_,
323  criterion_,
324  smootherArgs_);
325 
326  return inv; //std::shared_ptr<InverseOperator<X,X> >(inv);
327 
328  }
329 
330 private:
332  std::shared_ptr<Operator> coarseOperator_;
334  SmootherArgs smootherArgs_;
336  Criterion criterion_;
337 };
338 
344 template<class FO, class CSP, class S>
346  public Preconditioner<typename FO::domain_type, typename FO::range_type>
347 {
348 public:
352  typedef typename CoarseLevelSolverPolicy::CoarseLevelSolver CoarseLevelSolver;
357  typedef FO FineOperatorType;
361  typedef typename FineOperatorType::range_type FineRangeType;
365  typedef typename FineOperatorType::domain_type FineDomainType;
370  typedef typename CSP::Operator CoarseOperatorType;
374  typedef typename CoarseOperatorType::range_type CoarseRangeType;
378  typedef typename CoarseOperatorType::domain_type CoarseDomainType;
382  typedef S SmootherType;
383  // define the category
384  enum {
387  };
388 
404  std::shared_ptr<SmootherType> smoother,
406  CoarseOperatorType>& policy,
407  CoarseLevelSolverPolicy& coarsePolicy,
408  std::size_t preSteps=1, std::size_t postSteps=1)
409  : operator_(&op), smoother_(smoother),
410  preSteps_(preSteps), postSteps_(postSteps)
411  {
412  policy_ = policy.clone();
413  policy_->createCoarseLevelSystem(*operator_);
414  coarseSolver_=coarsePolicy.createCoarseLevelSolver(*policy_);
415  }
416 
418  : operator_(other.operator_), coarseSolver_(new CoarseLevelSolver(*other.coarseSolver_)),
419  smoother_(other.smoother_), policy_(other.policy_->clone()),
420  preSteps_(other.preSteps_), postSteps_(other.postSteps_)
421  {}
422 
424  {
425  // Each instance has its own policy.
426  delete policy_;
427  delete coarseSolver_;
428  }
429 
431  {
432  smoother_->pre(x,b);
433  }
434 
436  {
437  DUNE_UNUSED_PARAMETER(x);
438  }
439 
440  void apply(FineDomainType& v, const FineRangeType& d)
441  {
442  FineDomainType u(v);
443  FineRangeType rhs(d);
444  LevelContext context;
446  context.pinfo=&info;
447  context.lhs=&u;
448  context.update=&v;
449  context.smoother=smoother_;
450  context.rhs=&rhs;
451  context.matrix=operator_;
452  // Presmoothing
453  presmooth(context, preSteps_);
454  //Coarse grid correction
455  policy_->moveToCoarseLevel(*context.rhs);
457  coarseSolver_->apply(policy_->getCoarseLevelLhs(), policy_->getCoarseLevelRhs(), res);
458  *context.lhs=0;
459  policy_->moveToFineLevel(*context.lhs);
460  *context.update += *context.lhs;
461  // Postsmoothing
462  postsmooth(context, postSteps_);
463 
464  }
465 
466 private:
470  struct LevelContext
471  {
473  typedef S SmootherType;
475  std::shared_ptr<SmootherType> smoother;
477  FineDomainType* lhs;
478  /*
479  * @brief The right hand side holding the current residual.
480  *
481  * This is passed to the smoother as the right hand side.
482  */
483  FineRangeType* rhs;
489  FineDomainType* update;
491  SequentialInformation* pinfo;
497  const FineOperatorType* matrix;
498  };
499  const FineOperatorType* operator_;
501  CoarseLevelSolver* coarseSolver_;
503  std::shared_ptr<S> smoother_;
505  LevelTransferPolicy<FO,typename CSP::Operator>* policy_;
507  std::size_t preSteps_;
509  std::size_t postSteps_;
510 };
511 }// end namespace Amg
512 }// end namespace Dune
513 
515 #endif
Definition: galerkin.hh:115
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:430
std::shared_ptr< CoarseOperatorType > operator_
the coarse level linear operator.
Definition: twolevelmethod.hh:130
Provides a class for building the galerkin product based on a aggregation scheme. ...
Definition: twolevelmethod.hh:345
Define general, extensible interface for operators. The available implementation wraps a matrix...
CSP::Operator CoarseOperatorType
The linear operator of the finel level system. Has to be derived from AssembledLinearOperator.
Definition: twolevelmethod.hh:370
CoarseOperatorType::range_type CoarseRangeType
The type of the range of the coarse level operator.
Definition: twolevelmethod.hh:59
G::MutableMatrix * build(G &fineGraph, V &visitedMap, const ParallelInformation &pinfo, AggregatesMap< typename G::VertexDescriptor > &aggregates, const typename G::Matrix::size_type &size, const Set &copy)
Calculates the coarse matrix via a Galerkin product.
Definition: galerkin.hh:564
virtual void apply(X &x, X &b, InverseOperatorResult &res)=0
Apply inverse operator,.
OneStepAMGCoarseSolverPolicy(const SmootherArgs &args, const Criterion &c)
Constructs the coarse solver policy.
Definition: twolevelmethod.hh:251
CoarseLevelSolver * createCoarseLevelSolver(P &transferPolicy)
Constructs a coarse level solver.
Definition: twolevelmethod.hh:319
FineOperatorType::domain_type FineDomainType
The type of the domain of the fine level operator.
Definition: twolevelmethod.hh:50
void post(FineDomainType &x)
Definition: twolevelmethod.hh:435
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:26
virtual void moveToFineLevel(FineDomainType &fineLhs)=0
Updates the fine level linear system after the correction of the coarse levels system.
void calculate(const M &fine, const AggregatesMap< V > &aggregates, M &coarse, const I &pinfo, const O &copy)
Calculate the galerkin product.
AMGInverseOperator CoarseLevelSolver
The type of solver constructed for the coarse level.
Definition: twolevelmethod.hh:309
CoarseDomainType lhs_
The coarse level lhs.
Definition: twolevelmethod.hh:128
static void prolongateVector(const AggregatesMap< Vertex > &aggregates, Vector &coarse, Vector &fine, Vector &fineRedist, T1 damp, R &redistributor=R())
The AMG preconditioner.
TwoLevelMethod(const FineOperatorType &op, std::shared_ptr< SmootherType > smoother, const LevelTransferPolicy< FineOperatorType, CoarseOperatorType > &policy, CoarseLevelSolverPolicy &coarsePolicy, std::size_t preSteps=1, std::size_t postSteps=1)
Constructs a two level method.
Definition: twolevelmethod.hh:403
AggregationLevelTransferPolicy(const Criterion &crit)
Definition: twolevelmethod.hh:148
Class representing the properties of an ede in the matrix graph.
Definition: dependency.hh:37
O::range_type X
The type of the range and domain of the operator.
Definition: twolevelmethod.hh:237
Define general, extensible interface for inverse operators.
CoarseOperatorType::range_type CoarseRangeType
The type of the range of the coarse level operator.
Definition: twolevelmethod.hh:374
void apply(FineDomainType &v, const FineRangeType &d)
Definition: twolevelmethod.hh:440
~TwoLevelMethod()
Definition: twolevelmethod.hh:423
void pre(FineDomainType &x, FineRangeType &b)
Definition: twolevelmethod.hh:430
Category for sequential solvers.
Definition: solvercategory.hh:21
A policy class for solving the coarse level system using one step of AMG.
Definition: twolevelmethod.hh:231
Class providing information about the mapping of the vertices onto aggregates.
Definition: aggregates.hh:542
CoarseOperatorType::domain_type CoarseDomainType
The type of the domain of the coarse level operator.
Definition: twolevelmethod.hh:378
OneStepAMGCoarseSolverPolicy(const OneStepAMGCoarseSolverPolicy &other)
Copy constructor.
Definition: twolevelmethod.hh:255
void moveToCoarseLevel(const typename FatherType::FineRangeType &fineRhs)
Definition: twolevelmethod.hh:198
CoarseRangeType rhs_
The coarse level rhs.
Definition: twolevelmethod.hh:126
virtual LevelTransferPolicy * clone() const =0
Clone the current object.
Definition: indicescoarsener.hh:34
CO CoarseOperatorType
The linear operator of the finel level system. Has to be derived from AssembledLinearOperator.
Definition: twolevelmethod.hh:55
The default class for the smoother arguments.
Definition: smoother.hh:35
FineOperatorType::range_type FineRangeType
The type of the range of the fine level operator.
Definition: twolevelmethod.hh:46
FineOperatorType::range_type FineRangeType
The type of the range of the fine level operator.
Definition: twolevelmethod.hh:361
O Operator
The type of the linear operator used.
Definition: twolevelmethod.hh:235
C Criterion
The type of the crition used for the aggregation within AMG.
Definition: twolevelmethod.hh:239
The category the preconditioner is part of.
Definition: twolevelmethod.hh:386
Abstract base class for all solvers.
Definition: solver.hh:79
Dune::Amg::SmootherTraits< S >::Arguments SmootherArgs
The type of the arguments used for constructing the smoother.
Definition: twolevelmethod.hh:243
FO FineOperatorType
The linear operator of the finel level system. Has to be derived from AssembledLinearOperator.
Definition: twolevelmethod.hh:42
static void restrictVector(const AggregatesMap< Vertex > &aggregates, Vector &coarse, const Vector &fine, T &comm)
Abstract base class for transfer between levels and creation of the coarse level system.
Definition: twolevelmethod.hh:35
Operator Operator
The matrix operator type.
Definition: amg.hh:64
LevelTransferPolicy< O, O > FatherType
Definition: twolevelmethod.hh:144
CoarseRangeType & getCoarseLevelRhs()
Get the coarse level right hand side.
Definition: twolevelmethod.hh:76
CoarseDomainType & getCoarseLevelLhs()
Get the coarse level left hand side.
Definition: twolevelmethod.hh:85
C Criterion
Definition: twolevelmethod.hh:145
virtual ~LevelTransferPolicy()
Destructor.
Definition: twolevelmethod.hh:122
void moveToFineLevel(typename FatherType::FineDomainType &fineLhs)
Updates the fine level linear system after the correction of the coarse levels system.
Definition: twolevelmethod.hh:205
AMG< Operator, X, Smoother > AMGType
The type of the AMG construct on the coarse level.
Definition: twolevelmethod.hh:245
std::shared_ptr< CoarseOperatorType > & getCoarseLevelOperator()
Get the coarse level operator.
Definition: twolevelmethod.hh:68
TwoLevelMethod(const TwoLevelMethod &other)
Definition: twolevelmethod.hh:417
S Smoother
The type of the smoother used in AMG.
Definition: twolevelmethod.hh:241
SequentialInformation ParallelInformation
Definition: twolevelmethod.hh:146
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:91
Class representing a node in the matrix graph.
Definition: dependency.hh:125
Statistics about the application of an inverse operator.
Definition: solver.hh:31
CSP CoarseLevelSolverPolicy
The type of the policy for constructing the coarse level solver.
Definition: twolevelmethod.hh:350
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:408
CoarseLevelSolverPolicy::CoarseLevelSolver CoarseLevelSolver
The type of the coarse level solver.
Definition: twolevelmethod.hh:352
S SmootherType
The type of the fine level smoother.
Definition: twolevelmethod.hh:382
Definition: pinfo.hh:25
virtual void moveToCoarseLevel(const FineRangeType &fineRhs)=0
Transfers the data to the coarse level.
FineOperatorType::domain_type FineDomainType
The type of the domain of the fine level operator.
Definition: twolevelmethod.hh:365
virtual void createCoarseLevelSystem(const FineOperatorType &fineOperator)=0
Algebraically creates the coarse level system.
The (undirected) graph of a matrix.
Definition: graph.hh:48
A LeveTransferPolicy that used aggregation to construct the coarse level system.
Definition: twolevelmethod.hh:139
Attaches properties to the edges and vertices of a graph.
Definition: graph.hh:975
AggregationLevelTransferPolicy * clone() const
Clone the current object.
Definition: twolevelmethod.hh:212
FO FineOperatorType
The linear operator of the finel level system. Has to be derived from AssembledLinearOperator.
Definition: twolevelmethod.hh:357
CoarseOperatorType::domain_type CoarseDomainType
The type of the domain of the coarse level operator.
Definition: twolevelmethod.hh:63
void createCoarseLevelSystem(const O &fineOperator)
Algebraically creates the coarse level system.
Definition: twolevelmethod.hh:152