Dune Core Modules (2.5.2)

smoother.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_AMGSMOOTHER_HH
4 #define DUNE_AMGSMOOTHER_HH
5 
9 #include <dune/istl/schwarz.hh>
10 #include <dune/istl/novlpschwarz.hh>
11 #include <dune/common/propertymap.hh>
12 #include <dune/common/unused.hh>
13 
14 namespace Dune
15 {
16  namespace Amg
17  {
18 
34  template<class T>
36  {
40  typedef T RelaxationFactor;
41 
50 
55  : iterations(1), relaxationFactor(1.0)
56  {}
57  };
58 
62  template<class T>
64  {
66 
67  };
68 
69  template<class X, class Y, class C, class T>
70  struct SmootherTraits<BlockPreconditioner<X,Y,C,T> >
71  {
73 
74  };
75 
76  template<class C, class T>
77  struct SmootherTraits<NonoverlappingBlockPreconditioner<C,T> >
78  {
79  typedef DefaultSmootherArgs<typename T::matrix_type::field_type> Arguments;
80 
81  };
82 
86  template<class T>
88  {
89  typedef typename T::matrix_type Matrix;
90 
92 
94 
95  public:
96  virtual ~DefaultConstructionArgs()
97  {}
98 
99  void setMatrix(const Matrix& matrix)
100  {
101  matrix_=&matrix;
102  }
103  virtual void setMatrix(const Matrix& matrix, const AggregatesMap& amap)
104  {
105  DUNE_UNUSED_PARAMETER(amap);
106  setMatrix(matrix);
107  }
108 
109 
110  const Matrix& getMatrix() const
111  {
112  return *matrix_;
113  }
114 
115  void setArgs(const SmootherArgs& args)
116  {
117  args_=&args;
118  }
119 
120  template<class T1>
121  void setComm(T1& comm)
122  {
123  DUNE_UNUSED_PARAMETER(comm);
124  }
125 
126  const SequentialInformation& getComm()
127  {
128  return comm_;
129  }
130 
131  const SmootherArgs getArgs() const
132  {
133  return *args_;
134  }
135 
136  protected:
137  const Matrix* matrix_;
138  private:
139  const SmootherArgs* args_;
140  SequentialInformation comm_;
141  };
142 
143  template<class T>
144  struct ConstructionArgs
145  : public DefaultConstructionArgs<T>
146  {};
147 
148  template<class T, class C=SequentialInformation>
149  class DefaultParallelConstructionArgs
150  : public ConstructionArgs<T>
151  {
152  public:
153  virtual ~DefaultParallelConstructionArgs()
154  {}
155 
156  void setComm(const C& comm)
157  {
158  comm_ = &comm;
159  }
160 
161  const C& getComm() const
162  {
163  return *comm_;
164  }
165  private:
166  const C* comm_;
167  };
168 
169 
170  template<class T>
171  class ConstructionTraits;
172 
176  template<class M, class X, class Y, int l>
177  struct ConstructionTraits<SeqSSOR<M,X,Y,l> >
178  {
180 
181  static inline SeqSSOR<M,X,Y,l>* construct(Arguments& args)
182  {
183  return new SeqSSOR<M,X,Y,l>(args.getMatrix(), args.getArgs().iterations,
184  args.getArgs().relaxationFactor);
185  }
186 
187  static inline void deconstruct(SeqSSOR<M,X,Y,l>* ssor)
188  {
189  delete ssor;
190  }
191 
192  };
193 
194 
198  template<class M, class X, class Y, int l>
199  struct ConstructionTraits<SeqSOR<M,X,Y,l> >
200  {
202 
203  static inline SeqSOR<M,X,Y,l>* construct(Arguments& args)
204  {
205  return new SeqSOR<M,X,Y,l>(args.getMatrix(), args.getArgs().iterations,
206  args.getArgs().relaxationFactor);
207  }
208 
209  static inline void deconstruct(SeqSOR<M,X,Y,l>* sor)
210  {
211  delete sor;
212  }
213 
214  };
218  template<class M, class X, class Y, int l>
219  struct ConstructionTraits<SeqJac<M,X,Y,l> >
220  {
222 
223  static inline SeqJac<M,X,Y,l>* construct(Arguments& args)
224  {
225  return new SeqJac<M,X,Y,l>(args.getMatrix(), args.getArgs().iterations,
226  args.getArgs().relaxationFactor);
227  }
228 
229  static void deconstruct(SeqJac<M,X,Y,l>* jac)
230  {
231  delete jac;
232  }
233 
234  };
235 
236 
240  template<class M, class X, class Y>
241  struct ConstructionTraits<SeqILU0<M,X,Y> >
242  {
244 
245  static inline SeqILU0<M,X,Y>* construct(Arguments& args)
246  {
247  return new SeqILU0<M,X,Y>(args.getMatrix(),
248  args.getArgs().relaxationFactor);
249  }
250 
251  static void deconstruct(SeqILU0<M,X,Y>* ilu)
252  {
253  delete ilu;
254  }
255 
256  };
257 
258  template<class M, class X, class Y>
259  class ConstructionArgs<SeqILUn<M,X,Y> >
260  : public DefaultConstructionArgs<SeqILUn<M,X,Y> >
261  {
262  public:
263  ConstructionArgs(int n=1)
264  : n_(n)
265  {}
266 
267  void setN(int n)
268  {
269  n_ = n;
270  }
271  int getN()
272  {
273  return n_;
274  }
275 
276  private:
277  int n_;
278  };
279 
280 
284  template<class M, class X, class Y>
285  struct ConstructionTraits<SeqILUn<M,X,Y> >
286  {
287  typedef ConstructionArgs<SeqILUn<M,X,Y> > Arguments;
288 
289  static inline SeqILUn<M,X,Y>* construct(Arguments& args)
290  {
291  return new SeqILUn<M,X,Y>(args.getMatrix(), args.getN(),
292  args.getArgs().relaxationFactor);
293  }
294 
295  static void deconstruct(SeqILUn<M,X,Y>* ilu)
296  {
297  delete ilu;
298  }
299 
300  };
301 
305  template<class M, class X, class Y, class C>
306  struct ConstructionTraits<ParSSOR<M,X,Y,C> >
307  {
308  typedef DefaultParallelConstructionArgs<M,C> Arguments;
309 
310  static inline ParSSOR<M,X,Y,C>* construct(Arguments& args)
311  {
312  return new ParSSOR<M,X,Y,C>(args.getMatrix(), args.getArgs().iterations,
313  args.getArgs().relaxationFactor,
314  args.getComm());
315  }
316  static inline void deconstruct(ParSSOR<M,X,Y,C>* ssor)
317  {
318  delete ssor;
319  }
320  };
321 
322  template<class X, class Y, class C, class T>
323  struct ConstructionTraits<BlockPreconditioner<X,Y,C,T> >
324  {
325  typedef DefaultParallelConstructionArgs<T,C> Arguments;
326  typedef ConstructionTraits<T> SeqConstructionTraits;
327  static inline BlockPreconditioner<X,Y,C,T>* construct(Arguments& args)
328  {
329  return new BlockPreconditioner<X,Y,C,T>(*SeqConstructionTraits::construct(args),
330  args.getComm());
331  }
332 
333  static inline void deconstruct(BlockPreconditioner<X,Y,C,T>* bp)
334  {
335  SeqConstructionTraits::deconstruct(static_cast<T*>(&bp->preconditioner));
336  delete bp;
337  }
338 
339  };
340 
341  template<class C, class T>
342  struct ConstructionTraits<NonoverlappingBlockPreconditioner<C,T> >
343  {
344  typedef DefaultParallelConstructionArgs<T,C> Arguments;
345  typedef ConstructionTraits<T> SeqConstructionTraits;
346  static inline NonoverlappingBlockPreconditioner<C,T>* construct(Arguments& args)
347  {
348  return new NonoverlappingBlockPreconditioner<C,T>(*SeqConstructionTraits::construct(args),
349  args.getComm());
350  }
351 
352  static inline void deconstruct(NonoverlappingBlockPreconditioner<C,T>* bp)
353  {
354  SeqConstructionTraits::deconstruct(static_cast<T*>(&bp->preconditioner));
355  delete bp;
356  }
357 
358  };
359 
370  template<class T>
372  {
373  typedef T Smoother;
374  typedef typename Smoother::range_type Range;
375  typedef typename Smoother::domain_type Domain;
376 
384  static void preSmooth(Smoother& smoother, Domain& v, const Range& d)
385  {
386  smoother.apply(v,d);
387  }
388 
396  static void postSmooth(Smoother& smoother, Domain& v, const Range& d)
397  {
398  smoother.apply(v,d);
399  }
400  };
401 
407  template<typename LevelContext>
408  void presmooth(LevelContext& levelContext, size_t steps)
409  {
410  for(std::size_t i=0; i < steps; ++i) {
411  *levelContext.lhs=0;
413  ::preSmooth(*levelContext.smoother, *levelContext.lhs,
414  *levelContext.rhs);
415  // Accumulate update
416  *levelContext.update += *levelContext.lhs;
417 
418  // update defect
419  levelContext.matrix->applyscaleadd(-1, *levelContext.lhs, *levelContext.rhs);
420  levelContext.pinfo->project(*levelContext.rhs);
421  }
422  }
423 
429  template<typename LevelContext>
430  void postsmooth(LevelContext& levelContext, size_t steps)
431  {
432  for(std::size_t i=0; i < steps; ++i) {
433  // update defect
434  levelContext.matrix->applyscaleadd(-1, *levelContext.lhs,
435  *levelContext.rhs);
436  *levelContext.lhs=0;
437  levelContext.pinfo->project(*levelContext.rhs);
439  ::postSmooth(*levelContext.smoother, *levelContext.lhs, *levelContext.rhs);
440  // Accumulate update
441  *levelContext.update += *levelContext.lhs;
442  }
443  }
444 
445  template<class M, class X, class Y, int l>
446  struct SmootherApplier<SeqSOR<M,X,Y,l> >
447  {
448  typedef SeqSOR<M,X,Y,l> Smoother;
449  typedef typename Smoother::range_type Range;
450  typedef typename Smoother::domain_type Domain;
451 
452  static void preSmooth(Smoother& smoother, Domain& v, Range& d)
453  {
454  smoother.template apply<true>(v,d);
455  }
456 
457 
458  static void postSmooth(Smoother& smoother, Domain& v, Range& d)
459  {
460  smoother.template apply<false>(v,d);
461  }
462  };
463 
464  template<class M, class X, class Y, class C, int l>
465  struct SmootherApplier<BlockPreconditioner<X,Y,C,SeqSOR<M,X,Y,l> > >
466  {
467  typedef BlockPreconditioner<X,Y,C,SeqSOR<M,X,Y,l> > Smoother;
468  typedef typename Smoother::range_type Range;
469  typedef typename Smoother::domain_type Domain;
470 
471  static void preSmooth(Smoother& smoother, Domain& v, Range& d)
472  {
473  smoother.template apply<true>(v,d);
474  }
475 
476 
477  static void postSmooth(Smoother& smoother, Domain& v, Range& d)
478  {
479  smoother.template apply<false>(v,d);
480  }
481  };
482 
483  template<class M, class X, class Y, class C, int l>
484  struct SmootherApplier<NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y,l> > >
485  {
486  typedef NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y,l> > Smoother;
487  typedef typename Smoother::range_type Range;
488  typedef typename Smoother::domain_type Domain;
489 
490  static void preSmooth(Smoother& smoother, Domain& v, Range& d)
491  {
492  smoother.template apply<true>(v,d);
493  }
494 
495 
496  static void postSmooth(Smoother& smoother, Domain& v, Range& d)
497  {
498  smoother.template apply<false>(v,d);
499  }
500  };
501 
502  } // end namespace Amg
503 
504  // forward declarations
505  template<class M, class X, class MO, class MS, class A>
506  class SeqOverlappingSchwarz;
507 
508  struct MultiplicativeSchwarzMode;
509 
510  namespace Amg
511  {
512  template<class M, class X, class MS, class TA>
513  struct SmootherApplier<SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,
514  MS,TA> >
515  {
516  typedef SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,MS,TA> Smoother;
517  typedef typename Smoother::range_type Range;
518  typedef typename Smoother::domain_type Domain;
519 
520  static void preSmooth(Smoother& smoother, Domain& v, const Range& d)
521  {
522  smoother.template apply<true>(v,d);
523  }
524 
525 
526  static void postSmooth(Smoother& smoother, Domain& v, const Range& d)
527  {
528  smoother.template apply<false>(v,d);
529 
530  }
531  };
532 
533  // template<class M, class X, class TM, class TA>
534  // class SeqOverlappingSchwarz;
535 
536  template<class T>
537  struct SeqOverlappingSchwarzSmootherArgs
538  : public DefaultSmootherArgs<T>
539  {
540  enum Overlap {vertex, aggregate, pairwise, none};
541 
542  Overlap overlap;
543  bool onthefly;
544 
545  SeqOverlappingSchwarzSmootherArgs(Overlap overlap_=vertex,
546  bool onthefly_=false)
547  : overlap(overlap_), onthefly(onthefly_)
548  {}
549  };
550 
551  template<class M, class X, class TM, class TS, class TA>
552  struct SmootherTraits<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
553  {
554  typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> Arguments;
555  };
556 
557  template<class M, class X, class TM, class TS, class TA>
558  class ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
559  : public DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
560  {
561  typedef DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> > Father;
562 
563  public:
564  typedef typename MatrixGraph<M>::VertexDescriptor VertexDescriptor;
565  typedef Dune::Amg::AggregatesMap<VertexDescriptor> AggregatesMap;
566  typedef typename AggregatesMap::AggregateDescriptor AggregateDescriptor;
568  typedef typename Vector::value_type Subdomain;
569 
570  virtual void setMatrix(const M& matrix, const AggregatesMap& amap)
571  {
572  Father::setMatrix(matrix);
573 
574  std::vector<bool> visited(amap.noVertices(), false);
575  typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
576  VisitedMapType visitedMap(visited.begin());
577 
578  MatrixGraph<const M> graph(matrix);
579 
580  typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
581 
582  switch(Father::getArgs().overlap) {
583  case SmootherArgs::vertex :
584  {
585  VertexAdder visitor(subdomains, amap);
586  createSubdomains(matrix, graph, amap, visitor, visitedMap);
587  }
588  break;
589  case SmootherArgs::pairwise :
590  {
591  createPairDomains(graph);
592  }
593  break;
594  case SmootherArgs::aggregate :
595  {
596  AggregateAdder<VisitedMapType> visitor(subdomains, amap, graph, visitedMap);
597  createSubdomains(matrix, graph, amap, visitor, visitedMap);
598  }
599  break;
600  case SmootherArgs::none :
601  NoneAdder visitor;
602  createSubdomains(matrix, graph, amap, visitor, visitedMap);
603  break;
604  default :
605  DUNE_THROW(NotImplemented, "This overlapping scheme is not supported!");
606  }
607  }
608  void setMatrix(const M& matrix)
609  {
610  Father::setMatrix(matrix);
611 
612  /* Create aggregates map where each aggregate is just one vertex. */
613  AggregatesMap amap(matrix.N());
614  VertexDescriptor v=0;
615  for(typename AggregatesMap::iterator iter=amap.begin();
616  iter!=amap.end(); ++iter)
617  *iter=v++;
618 
619  std::vector<bool> visited(amap.noVertices(), false);
620  typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
621  VisitedMapType visitedMap(visited.begin());
622 
623  MatrixGraph<const M> graph(matrix);
624 
625  typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
626 
627  switch(Father::getArgs().overlap) {
628  case SmootherArgs::vertex :
629  {
630  VertexAdder visitor(subdomains, amap);
631  createSubdomains(matrix, graph, amap, visitor, visitedMap);
632  }
633  break;
634  case SmootherArgs::aggregate :
635  {
636  DUNE_THROW(NotImplemented, "Aggregate overlap is not supported yet");
637  /*
638  AggregateAdder<VisitedMapType> visitor(subdomains, amap, graph, visitedMap);
639  createSubdomains(matrix, graph, amap, visitor, visitedMap);
640  */
641  }
642  break;
643  case SmootherArgs::pairwise :
644  {
645  createPairDomains(graph);
646  }
647  break;
648  case SmootherArgs::none :
649  NoneAdder visitor;
650  createSubdomains(matrix, graph, amap, visitor, visitedMap);
651 
652  }
653  }
654 
655  const Vector& getSubDomains()
656  {
657  return subdomains;
658  }
659 
660  private:
661  struct VertexAdder
662  {
663  VertexAdder(Vector& subdomains_, const AggregatesMap& aggregates_)
664  : subdomains(subdomains_), max(-1), subdomain(-1), aggregates(aggregates_)
665  {}
666  template<class T>
667  void operator()(const T& edge)
668  {
669  if(aggregates[edge.target()]!=AggregatesMap::ISOLATED)
670  subdomains[subdomain].insert(edge.target());
671  }
672  int setAggregate(const AggregateDescriptor& aggregate_)
673  {
674  subdomain=aggregate_;
675  max = std::max(subdomain, aggregate_);
676  return subdomain;
677  }
678  int noSubdomains() const
679  {
680  return max+1;
681  }
682  private:
683  Vector& subdomains;
684  AggregateDescriptor max;
685  AggregateDescriptor subdomain;
686  const AggregatesMap& aggregates;
687  };
688  struct NoneAdder
689  {
690  template<class T>
691  void operator()(const T& edge)
692  {}
693  int setAggregate(const AggregateDescriptor& aggregate_)
694  {
695  return -1;
696  }
697  int noSubdomains() const
698  {
699  return -1;
700  }
701  };
702 
703  template<class VM>
704  struct AggregateAdder
705  {
706  AggregateAdder(Vector& subdomains_, const AggregatesMap& aggregates_,
707  const MatrixGraph<const M>& graph_, VM& visitedMap_)
708  : subdomains(subdomains_), subdomain(-1), aggregates(aggregates_),
709  adder(subdomains_, aggregates_), graph(graph_), visitedMap(visitedMap_)
710  {}
711  template<class T>
712  void operator()(const T& edge)
713  {
714  subdomains[subdomain].insert(edge.target());
715  // If we (the neighbouring vertex of the aggregate)
716  // are not isolated, add the aggregate we belong to
717  // to the same subdomain using the OneOverlapAdder
718  if(aggregates[edge.target()]!=AggregatesMap::ISOLATED) {
719  assert(aggregates[edge.target()]!=aggregate);
720  typename AggregatesMap::VertexList vlist;
721  aggregates.template breadthFirstSearch<true,false>(edge.target(), aggregate,
722  graph, vlist, adder, adder,
723  visitedMap);
724  }
725  }
726 
727  int setAggregate(const AggregateDescriptor& aggregate_)
728  {
729  adder.setAggregate(aggregate_);
730  aggregate=aggregate_;
731  return ++subdomain;
732  }
733  int noSubdomains() const
734  {
735  return subdomain+1;
736  }
737 
738  private:
739  AggregateDescriptor aggregate;
740  Vector& subdomains;
741  int subdomain;
742  const AggregatesMap& aggregates;
743  VertexAdder adder;
744  const MatrixGraph<const M>& graph;
745  VM& visitedMap;
746  };
747 
748  void createPairDomains(const MatrixGraph<const M>& graph)
749  {
750  typedef typename MatrixGraph<const M>::ConstVertexIterator VIter;
751  typedef typename MatrixGraph<const M>::ConstEdgeIterator EIter;
752  typedef typename M::size_type size_type;
753 
754  std::set<std::pair<size_type,size_type> > pairs;
755  int total=0;
756  for(VIter v=graph.begin(), ve=graph.end(); ve != v; ++v)
757  for(EIter e = v.begin(), ee=v.end(); ee!=e; ++e)
758  {
759  ++total;
760  if(e.source()<e.target())
761  pairs.insert(std::make_pair(e.source(),e.target()));
762  else
763  pairs.insert(std::make_pair(e.target(),e.source()));
764  }
765 
766 
767  subdomains.resize(pairs.size());
768  Dune::dinfo <<std::endl<< "Created "<<pairs.size()<<" ("<<total<<") pair domains"<<std::endl<<std::endl;
769  typedef typename std::set<std::pair<size_type,size_type> >::const_iterator SIter;
770  typename Vector::iterator subdomain=subdomains.begin();
771 
772  for(SIter s=pairs.begin(), se =pairs.end(); se!=s; ++s)
773  {
774  subdomain->insert(s->first);
775  subdomain->insert(s->second);
776  ++subdomain;
777  }
778  std::size_t minsize=10000;
779  std::size_t maxsize=0;
780  int sum=0;
781  for(typename Vector::size_type i=0; i < subdomains.size(); ++i) {
782  sum+=subdomains[i].size();
783  minsize=std::min(minsize, subdomains[i].size());
784  maxsize=std::max(maxsize, subdomains[i].size());
785  }
786  Dune::dinfo<<"Subdomain size: min="<<minsize<<" max="<<maxsize<<" avg="<<(sum/subdomains.size())
787  <<" no="<<subdomains.size()<<std::endl;
788  }
789 
790  template<class Visitor>
791  void createSubdomains(const M& matrix, const MatrixGraph<const M>& graph,
792  const AggregatesMap& amap, Visitor& overlapVisitor,
793  IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap>& visitedMap )
794  {
795  // count number ag aggregates. We assume that the
796  // aggregates are numbered consecutively from 0 except
797  // for the isolated ones. All isolated vertices form
798  // one aggregate, here.
799  int isolated=0;
800  AggregateDescriptor maxAggregate=0;
801 
802  for(std::size_t i=0; i < amap.noVertices(); ++i)
803  if(amap[i]==AggregatesMap::ISOLATED)
804  isolated++;
805  else
806  maxAggregate = std::max(maxAggregate, amap[i]);
807 
808  subdomains.resize(maxAggregate+1+isolated);
809 
810  // reset the subdomains
811  for(typename Vector::size_type i=0; i < subdomains.size(); ++i)
812  subdomains[i].clear();
813 
814  // Create the subdomains from the aggregates mapping.
815  // For each aggregate we mark all entries and the
816  // neighbouring vertices as belonging to the same subdomain
817  VertexAdder aggregateVisitor(subdomains, amap);
818 
819  for(VertexDescriptor i=0; i < amap.noVertices(); ++i)
820  if(!get(visitedMap, i)) {
821  AggregateDescriptor aggregate=amap[i];
822 
823  if(amap[i]==AggregatesMap::ISOLATED) {
824  // isolated vertex gets its own aggregate
825  subdomains.push_back(Subdomain());
826  aggregate=subdomains.size()-1;
827  }
828  overlapVisitor.setAggregate(aggregate);
829  aggregateVisitor.setAggregate(aggregate);
830  subdomains[aggregate].insert(i);
831  typename AggregatesMap::VertexList vlist;
832  amap.template breadthFirstSearch<false,false>(i, aggregate, graph, vlist, aggregateVisitor,
833  overlapVisitor, visitedMap);
834  }
835 
836  std::size_t minsize=10000;
837  std::size_t maxsize=0;
838  int sum=0;
839  for(typename Vector::size_type i=0; i < subdomains.size(); ++i) {
840  sum+=subdomains[i].size();
841  minsize=std::min(minsize, subdomains[i].size());
842  maxsize=std::max(maxsize, subdomains[i].size());
843  }
844  Dune::dinfo<<"Subdomain size: min="<<minsize<<" max="<<maxsize<<" avg="<<(sum/subdomains.size())
845  <<" no="<<subdomains.size()<<" isolated="<<isolated<<std::endl;
846 
847 
848 
849  }
850  Vector subdomains;
851  };
852 
853 
854  template<class M, class X, class TM, class TS, class TA>
855  struct ConstructionTraits<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
856  {
857  typedef ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> > Arguments;
858 
859  static inline SeqOverlappingSchwarz<M,X,TM,TS,TA>* construct(Arguments& args)
860  {
861  return new SeqOverlappingSchwarz<M,X,TM,TS,TA>(args.getMatrix(),
862  args.getSubDomains(),
863  args.getArgs().relaxationFactor,
864  args.getArgs().onthefly);
865  }
866 
867  static void deconstruct(SeqOverlappingSchwarz<M,X,TM,TS,TA>* schwarz)
868  {
869  delete schwarz;
870  }
871  };
872 
873 
874  } // namespace Amg
875 } // namespace Dune
876 
877 
878 
879 #endif
Provides classes for the Coloring process of AMG.
Class providing information about the mapping of the vertices onto aggregates.
Definition: aggregates.hh:543
Traits class for generically constructing non default constructable types.
Definition: construction.hh:38
Construction Arguments for the default smoothers.
Definition: smoother.hh:88
VertexIteratorT< const MatrixGraph< Matrix > > ConstVertexIterator
The constant vertex iterator type.
Definition: graph.hh:306
M::size_type VertexDescriptor
The vertex descriptor.
Definition: graph.hh:71
EdgeIteratorT< const MatrixGraph< Matrix > > ConstEdgeIterator
The constant edge iterator type.
Definition: graph.hh:296
Block parallel preconditioner.
Definition: schwarz.hh:356
Nonoverlapping parallel preconditioner.
Definition: novlpschwarz.hh:342
A parallel SSOR preconditioner.
Definition: schwarz.hh:253
Sequential ILU0 preconditioner.
Definition: preconditioners.hh:488
Sequential ILU(n) preconditioner.
Definition: preconditioners.hh:572
The sequential jacobian preconditioner.
Definition: preconditioners.hh:402
std::vector< subdomain_type, typename TA::template rebind< subdomain_type >::other > subdomain_vector
The vector type containing the subdomain to row index mapping.
Definition: overlappingschwarz.hh:784
Sequential SOR preconditioner.
Definition: preconditioners.hh:215
Sequential SSOR preconditioner.
Definition: preconditioners.hh:127
Helper classes for the construction of classes without empty constructor.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
static T * construct(Arguments &args)
Construct an object with the specified arguments.
Definition: construction.hh:52
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:408
V AggregateDescriptor
The aggregate descriptor type.
Definition: aggregates.hh:563
static const V ISOLATED
Identifier of isolated vertices.
Definition: aggregates.hh:554
DefaultSmootherArgs()
Default constructor.
Definition: smoother.hh:54
const void * Arguments
A type holding all the arguments needed to call the constructor.
Definition: construction.hh:44
static void postSmooth(Smoother &smoother, Domain &v, const Range &d)
apply post smoothing in forward direction
Definition: smoother.hh:396
SLList< VertexDescriptor, Allocator > VertexList
The type of a single linked list of vertex descriptors.
Definition: aggregates.hh:575
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:430
RelaxationFactor relaxationFactor
The relaxation factor to use.
Definition: smoother.hh:49
static void preSmooth(Smoother &smoother, Domain &v, const Range &d)
apply pre smoothing in forward direction
Definition: smoother.hh:384
static void deconstruct(T *t)
Destroys an object.
Definition: construction.hh:61
int iterations
The numbe of iterations to perform.
Definition: smoother.hh:45
T RelaxationFactor
The type of the relaxation factor.
Definition: smoother.hh:40
DInfoType dinfo(std::cout)
Stream for informative output.
Definition: stdstreams.hh:138
PartitionSet<... > Overlap
Type of PartitionSet for the overlap partition.
Definition: partitionset.hh:198
Dune namespace.
Definition: alignment.hh:11
Define general preconditioner interface.
The default class for the smoother arguments.
Definition: smoother.hh:36
Helper class for applying the smoothers.
Definition: smoother.hh:372
Traits class for getting the attribute class of a smoother.
Definition: smoother.hh:64
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 9, 22:29, 2024)