dune-istl  2.2.0
kamg.hh
Go to the documentation of this file.
1 #ifndef DUNE_AMG_KAMG_HH
2 #define DUNE_AMG_KAMG_HH
3 
5 #include"amg.hh"
6 
7 namespace Dune
8 {
9  namespace Amg
10  {
11 
26  template<class AMG>
28  : public Preconditioner<typename AMG::Domain,typename AMG::Range>
29  {
31  typedef typename AMG::Domain Domain;
33  typedef typename AMG::Range Range;
34  public:
35 
36  enum {
39  };
40 
49  : amg_(amg), coarseSolver_(coarseSolver)
50  {}
51 
53  void pre(typename AMG::Domain& x, typename AMG::Range& b)
54  {}
55 
57  void post(typename AMG::Domain& x)
58  {}
59 
61  void apply(typename AMG::Domain& v, const typename AMG::Range& d)
62  {
63  *amg_.rhs = d;
64  *amg_.lhs = v;
65  // Copy data
66  *amg_.update=0;
67 
68  amg_.presmooth();
69  bool processFineLevel =
70  amg_.moveToCoarseLevel();
71 
72  if(processFineLevel){
73  typename AMG::Range b=*amg_.rhs;
74  typename AMG::Domain x=*amg_.update;
76  coarseSolver_->apply(x, b, res);
77  *amg_.update=x;
78  }
79 
80  amg_.moveToFineLevel(processFineLevel);
81 
82  amg_.postsmooth();
83  v=*amg_.update;
84  }
85 
91  {
92  return coarseSolver_;
93  }
94 
97  {}
98 
99  private:
101  AMG& amg_;
103  InverseOperator<Domain,Range>* coarseSolver_;
104  };
105 
106 
107 
118  template<class M, class X, class S, class PI=SequentialInformation,
119  class K=BiCGSTABSolver<X>, class A=std::allocator<X> >
120  class KAMG : public Preconditioner<X,X>
121  {
122  public:
126  typedef K KrylovSolver;
130  typedef typename Amg::CoarseSolver CoarseSolver;
134  typedef typename Amg::SmootherArgs SmootherArgs;
136  typedef typename Amg::Operator Operator;
138  typedef typename Amg::Domain Domain;
140  typedef typename Amg::Range Range;
144  typedef typename Amg::ScalarProduct ScalarProduct;
145 
146  enum {
149  };
163  KAMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
164  const SmootherArgs& smootherArgs, std::size_t gamma,
165  std::size_t preSmoothingSteps =1, std::size_t postSmoothingSteps = 1,
166  std::size_t maxLevelKrylovSteps = 3 , double minDefectReduction =1e-1);
167 
184  template<class C>
185  KAMG(const Operator& fineOperator, const C& criterion,
186  const SmootherArgs& smootherArgs, std::size_t gamma=1,
187  std::size_t preSmoothingSteps=1, std::size_t postSmoothingSteps=1,
188  std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1,
190 
192  void pre(Domain& x, Range& b);
194  void post(Domain& x);
196  void apply(Domain& x, const Range& b);
197 
198  std::size_t maxlevels();
199 
200  private:
202  Amg amg;
203 
205  std::size_t maxLevelKrylovSteps;
206 
208  double levelDefectReduction;
209 
211  std::vector<typename Amg::ScalarProduct*> scalarproducts;
212 
214  std::vector<KAmgTwoGrid<Amg>*> ksolvers;
215  };
216 
217  template<class M, class X, class S, class P, class K, class A>
218  KAMG<M,X,S,P,K,A>::KAMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
219  const SmootherArgs& smootherArgs,
220  std::size_t gamma, std::size_t preSmoothingSteps,
221  std::size_t postSmoothingSteps,
222  std::size_t ksteps, double reduction)
223  : amg(matrices, coarseSolver, smootherArgs, gamma, preSmoothingSteps,
224  postSmoothingSteps), maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
225  {}
226 
227  template<class M, class X, class S, class P, class K, class A>
228  template<class C>
229  KAMG<M,X,S,P,K,A>::KAMG(const Operator& fineOperator, const C& criterion,
230  const SmootherArgs& smootherArgs, std::size_t gamma,
231  std::size_t preSmoothingSteps, std::size_t postSmoothingSteps,
232  std::size_t ksteps, double reduction,
233  const ParallelInformation& pinfo)
234  : amg(fineOperator, criterion, smootherArgs, gamma, preSmoothingSteps,
235  postSmoothingSteps, false, pinfo), maxLevelKrylovSteps(ksteps), levelDefectReduction(reduction)
236  {}
237 
238 
239  template<class M, class X, class S, class P, class K, class A>
241  {
242  amg.pre(x,b);
243  scalarproducts.reserve(amg.levels());
244  ksolvers.reserve(amg.levels());
245 
246  typename OperatorHierarchy::ParallelMatrixHierarchy::Iterator
247  matrix = amg.matrices_->matrices().coarsest();
249  pinfo = amg.matrices_->parallelInformation().coarsest();
250  bool hasCoarsest=(amg.levels()==amg.maxlevels());
251  KrylovSolver* ks=0;
252 
253  if(hasCoarsest){
254  if(matrix==amg.matrices_->matrices().finest())
255  return;
256  --matrix;
257  --pinfo;
258  ksolvers.push_back(new KAmgTwoGrid<Amg>(amg, amg.solver_));
259  }else
260  ksolvers.push_back(new KAmgTwoGrid<Amg>(amg, ks));
261 
262  std::ostringstream s;
263 
264  if(matrix!=amg.matrices_->matrices().finest())
265  while(true){
266  scalarproducts.push_back(Amg::ScalarProductChooser::construct(*pinfo));
267  ks = new KrylovSolver(*matrix, *(scalarproducts.back()),
268  *(ksolvers.back()), levelDefectReduction,
269  maxLevelKrylovSteps, 0);
270  ksolvers.push_back(new KAmgTwoGrid<Amg>(amg, ks));
271  --matrix;
272  --pinfo;
273  if(matrix==amg.matrices_->matrices().finest())
274  break;
275  }
276  }
277 
278 
279  template<class M, class X, class S, class P, class K, class A>
281  {
282  typedef typename std::vector<KAmgTwoGrid<Amg>*>::reverse_iterator KIterator;
283  typedef typename std::vector<ScalarProduct*>::iterator SIterator;
284  for(KIterator kiter = ksolvers.rbegin(); kiter != ksolvers.rend();
285  ++kiter)
286  {
287  if((*kiter)->coarseSolver()!=amg.solver_)
288  delete (*kiter)->coarseSolver();
289  delete *kiter;
290  }
291  for(SIterator siter = scalarproducts.begin(); siter!=scalarproducts.end();
292  ++siter)
293  delete *siter;
294  amg.post(x);
295 
296  }
297 
298  template<class M, class X, class S, class P, class K, class A>
300  {
301  amg.initIteratorsWithFineLevel();
302  if(ksolvers.size()==0)
303  {
304  Range tb=b;
306  amg.solver_->apply(x,tb,res);
307  }else
308  ksolvers.back()->apply(x,b);
309  }
310 
311  template<class M, class X, class S, class P, class K, class A>
313  {
314  return amg.maxlevels();
315  }
316 
318  } // Amg
319 } // Dune
320 
321 #endif