Dune Core Modules (2.6.0)

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 
302 
303  template<class M, class X, class Y>
304  class ConstructionArgs<SeqILU<M,X,Y> >
305  : public DefaultConstructionArgs<SeqILU<M,X,Y> >
306  {
307  public:
308  ConstructionArgs(int n=0)
309  : n_(n)
310  {}
311 
312  void setN(int n)
313  {
314  n_ = n;
315  }
316 
317  int getN()
318  {
319  return n_;
320  }
321 
322  private:
323  int n_;
324  };
325 
326 
330  template<class M, class X, class Y>
331  struct ConstructionTraits<SeqILU<M,X,Y> >
332  {
333  typedef ConstructionArgs<SeqILU<M,X,Y> > Arguments;
334 
335  static inline SeqILU<M,X,Y>* construct(Arguments& args)
336  {
337  return new SeqILU<M,X,Y>(args.getMatrix(), args.getN(),
338  args.getArgs().relaxationFactor);
339  }
340 
341  static void deconstruct(SeqILU<M,X,Y>* ilu)
342  {
343  delete ilu;
344  }
345 
346  };
347 
351  template<class M, class X, class Y, class C>
352  struct ConstructionTraits<ParSSOR<M,X,Y,C> >
353  {
354  typedef DefaultParallelConstructionArgs<M,C> Arguments;
355 
356  static inline ParSSOR<M,X,Y,C>* construct(Arguments& args)
357  {
358  return new ParSSOR<M,X,Y,C>(args.getMatrix(), args.getArgs().iterations,
359  args.getArgs().relaxationFactor,
360  args.getComm());
361  }
362  static inline void deconstruct(ParSSOR<M,X,Y,C>* ssor)
363  {
364  delete ssor;
365  }
366  };
367 
368  template<class X, class Y, class C, class T>
369  struct ConstructionTraits<BlockPreconditioner<X,Y,C,T> >
370  {
371  typedef DefaultParallelConstructionArgs<T,C> Arguments;
372  typedef ConstructionTraits<T> SeqConstructionTraits;
373  static inline BlockPreconditioner<X,Y,C,T>* construct(Arguments& args)
374  {
375  return new BlockPreconditioner<X,Y,C,T>(*SeqConstructionTraits::construct(args),
376  args.getComm());
377  }
378 
379  static inline void deconstruct(BlockPreconditioner<X,Y,C,T>* bp)
380  {
381  SeqConstructionTraits::deconstruct(static_cast<T*>(&bp->preconditioner));
382  delete bp;
383  }
384 
385  };
386 
387  template<class C, class T>
388  struct ConstructionTraits<NonoverlappingBlockPreconditioner<C,T> >
389  {
390  typedef DefaultParallelConstructionArgs<T,C> Arguments;
391  typedef ConstructionTraits<T> SeqConstructionTraits;
392  static inline NonoverlappingBlockPreconditioner<C,T>* construct(Arguments& args)
393  {
394  return new NonoverlappingBlockPreconditioner<C,T>(*SeqConstructionTraits::construct(args),
395  args.getComm());
396  }
397 
398  static inline void deconstruct(NonoverlappingBlockPreconditioner<C,T>* bp)
399  {
400  SeqConstructionTraits::deconstruct(static_cast<T*>(&bp->preconditioner));
401  delete bp;
402  }
403 
404  };
405 
416  template<class T>
418  {
419  typedef T Smoother;
420  typedef typename Smoother::range_type Range;
421  typedef typename Smoother::domain_type Domain;
422 
430  static void preSmooth(Smoother& smoother, Domain& v, const Range& d)
431  {
432  smoother.apply(v,d);
433  }
434 
442  static void postSmooth(Smoother& smoother, Domain& v, const Range& d)
443  {
444  smoother.apply(v,d);
445  }
446  };
447 
453  template<typename LevelContext>
454  void presmooth(LevelContext& levelContext, size_t steps)
455  {
456  for(std::size_t i=0; i < steps; ++i) {
457  *levelContext.lhs=0;
459  ::preSmooth(*levelContext.smoother, *levelContext.lhs,
460  *levelContext.rhs);
461  // Accumulate update
462  *levelContext.update += *levelContext.lhs;
463 
464  // update defect
465  levelContext.matrix->applyscaleadd(-1, *levelContext.lhs, *levelContext.rhs);
466  levelContext.pinfo->project(*levelContext.rhs);
467  }
468  }
469 
475  template<typename LevelContext>
476  void postsmooth(LevelContext& levelContext, size_t steps)
477  {
478  for(std::size_t i=0; i < steps; ++i) {
479  // update defect
480  levelContext.matrix->applyscaleadd(-1, *levelContext.lhs,
481  *levelContext.rhs);
482  *levelContext.lhs=0;
483  levelContext.pinfo->project(*levelContext.rhs);
485  ::postSmooth(*levelContext.smoother, *levelContext.lhs, *levelContext.rhs);
486  // Accumulate update
487  *levelContext.update += *levelContext.lhs;
488  }
489  }
490 
491  template<class M, class X, class Y, int l>
492  struct SmootherApplier<SeqSOR<M,X,Y,l> >
493  {
494  typedef SeqSOR<M,X,Y,l> Smoother;
495  typedef typename Smoother::range_type Range;
496  typedef typename Smoother::domain_type Domain;
497 
498  static void preSmooth(Smoother& smoother, Domain& v, Range& d)
499  {
500  smoother.template apply<true>(v,d);
501  }
502 
503 
504  static void postSmooth(Smoother& smoother, Domain& v, Range& d)
505  {
506  smoother.template apply<false>(v,d);
507  }
508  };
509 
510  template<class M, class X, class Y, class C, int l>
511  struct SmootherApplier<BlockPreconditioner<X,Y,C,SeqSOR<M,X,Y,l> > >
512  {
513  typedef BlockPreconditioner<X,Y,C,SeqSOR<M,X,Y,l> > Smoother;
514  typedef typename Smoother::range_type Range;
515  typedef typename Smoother::domain_type Domain;
516 
517  static void preSmooth(Smoother& smoother, Domain& v, Range& d)
518  {
519  smoother.template apply<true>(v,d);
520  }
521 
522 
523  static void postSmooth(Smoother& smoother, Domain& v, Range& d)
524  {
525  smoother.template apply<false>(v,d);
526  }
527  };
528 
529  template<class M, class X, class Y, class C, int l>
530  struct SmootherApplier<NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y,l> > >
531  {
532  typedef NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y,l> > Smoother;
533  typedef typename Smoother::range_type Range;
534  typedef typename Smoother::domain_type Domain;
535 
536  static void preSmooth(Smoother& smoother, Domain& v, Range& d)
537  {
538  smoother.template apply<true>(v,d);
539  }
540 
541 
542  static void postSmooth(Smoother& smoother, Domain& v, Range& d)
543  {
544  smoother.template apply<false>(v,d);
545  }
546  };
547 
548  } // end namespace Amg
549 
550  // forward declarations
551  template<class M, class X, class MO, class MS, class A>
552  class SeqOverlappingSchwarz;
553 
554  struct MultiplicativeSchwarzMode;
555 
556  namespace Amg
557  {
558  template<class M, class X, class MS, class TA>
559  struct SmootherApplier<SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,
560  MS,TA> >
561  {
562  typedef SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,MS,TA> Smoother;
563  typedef typename Smoother::range_type Range;
564  typedef typename Smoother::domain_type Domain;
565 
566  static void preSmooth(Smoother& smoother, Domain& v, const Range& d)
567  {
568  smoother.template apply<true>(v,d);
569  }
570 
571 
572  static void postSmooth(Smoother& smoother, Domain& v, const Range& d)
573  {
574  smoother.template apply<false>(v,d);
575 
576  }
577  };
578 
579  // template<class M, class X, class TM, class TA>
580  // class SeqOverlappingSchwarz;
581 
582  template<class T>
583  struct SeqOverlappingSchwarzSmootherArgs
584  : public DefaultSmootherArgs<T>
585  {
586  enum Overlap {vertex, aggregate, pairwise, none};
587 
588  Overlap overlap;
589  bool onthefly;
590 
591  SeqOverlappingSchwarzSmootherArgs(Overlap overlap_=vertex,
592  bool onthefly_=false)
593  : overlap(overlap_), onthefly(onthefly_)
594  {}
595  };
596 
597  template<class M, class X, class TM, class TS, class TA>
598  struct SmootherTraits<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
599  {
600  typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> Arguments;
601  };
602 
603  template<class M, class X, class TM, class TS, class TA>
604  class ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
605  : public DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
606  {
607  typedef DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> > Father;
608 
609  public:
610  typedef typename MatrixGraph<M>::VertexDescriptor VertexDescriptor;
611  typedef Dune::Amg::AggregatesMap<VertexDescriptor> AggregatesMap;
612  typedef typename AggregatesMap::AggregateDescriptor AggregateDescriptor;
614  typedef typename Vector::value_type Subdomain;
615 
616  virtual void setMatrix(const M& matrix, const AggregatesMap& amap)
617  {
618  Father::setMatrix(matrix);
619 
620  std::vector<bool> visited(amap.noVertices(), false);
621  typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
622  VisitedMapType visitedMap(visited.begin());
623 
624  MatrixGraph<const M> graph(matrix);
625 
626  typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
627 
628  switch(Father::getArgs().overlap) {
629  case SmootherArgs::vertex :
630  {
631  VertexAdder visitor(subdomains, amap);
632  createSubdomains(matrix, graph, amap, visitor, visitedMap);
633  }
634  break;
635  case SmootherArgs::pairwise :
636  {
637  createPairDomains(graph);
638  }
639  break;
640  case SmootherArgs::aggregate :
641  {
642  AggregateAdder<VisitedMapType> visitor(subdomains, amap, graph, visitedMap);
643  createSubdomains(matrix, graph, amap, visitor, visitedMap);
644  }
645  break;
646  case SmootherArgs::none :
647  NoneAdder visitor;
648  createSubdomains(matrix, graph, amap, visitor, visitedMap);
649  break;
650  default :
651  DUNE_THROW(NotImplemented, "This overlapping scheme is not supported!");
652  }
653  }
654  void setMatrix(const M& matrix)
655  {
656  Father::setMatrix(matrix);
657 
658  /* Create aggregates map where each aggregate is just one vertex. */
659  AggregatesMap amap(matrix.N());
660  VertexDescriptor v=0;
661  for(typename AggregatesMap::iterator iter=amap.begin();
662  iter!=amap.end(); ++iter)
663  *iter=v++;
664 
665  std::vector<bool> visited(amap.noVertices(), false);
666  typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
667  VisitedMapType visitedMap(visited.begin());
668 
669  MatrixGraph<const M> graph(matrix);
670 
671  typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
672 
673  switch(Father::getArgs().overlap) {
674  case SmootherArgs::vertex :
675  {
676  VertexAdder visitor(subdomains, amap);
677  createSubdomains(matrix, graph, amap, visitor, visitedMap);
678  }
679  break;
680  case SmootherArgs::aggregate :
681  {
682  DUNE_THROW(NotImplemented, "Aggregate overlap is not supported yet");
683  /*
684  AggregateAdder<VisitedMapType> visitor(subdomains, amap, graph, visitedMap);
685  createSubdomains(matrix, graph, amap, visitor, visitedMap);
686  */
687  }
688  break;
689  case SmootherArgs::pairwise :
690  {
691  createPairDomains(graph);
692  }
693  break;
694  case SmootherArgs::none :
695  NoneAdder visitor;
696  createSubdomains(matrix, graph, amap, visitor, visitedMap);
697 
698  }
699  }
700 
701  const Vector& getSubDomains()
702  {
703  return subdomains;
704  }
705 
706  private:
707  struct VertexAdder
708  {
709  VertexAdder(Vector& subdomains_, const AggregatesMap& aggregates_)
710  : subdomains(subdomains_), max(-1), subdomain(-1), aggregates(aggregates_)
711  {}
712  template<class T>
713  void operator()(const T& edge)
714  {
715  if(aggregates[edge.target()]!=AggregatesMap::ISOLATED)
716  subdomains[subdomain].insert(edge.target());
717  }
718  int setAggregate(const AggregateDescriptor& aggregate_)
719  {
720  subdomain=aggregate_;
721  max = std::max(subdomain, aggregate_);
722  return subdomain;
723  }
724  int noSubdomains() const
725  {
726  return max+1;
727  }
728  private:
729  Vector& subdomains;
730  AggregateDescriptor max;
731  AggregateDescriptor subdomain;
732  const AggregatesMap& aggregates;
733  };
734  struct NoneAdder
735  {
736  template<class T>
737  void operator()(const T& edge)
738  {}
739  int setAggregate(const AggregateDescriptor& aggregate_)
740  {
741  return -1;
742  }
743  int noSubdomains() const
744  {
745  return -1;
746  }
747  };
748 
749  template<class VM>
750  struct AggregateAdder
751  {
752  AggregateAdder(Vector& subdomains_, const AggregatesMap& aggregates_,
753  const MatrixGraph<const M>& graph_, VM& visitedMap_)
754  : subdomains(subdomains_), subdomain(-1), aggregates(aggregates_),
755  adder(subdomains_, aggregates_), graph(graph_), visitedMap(visitedMap_)
756  {}
757  template<class T>
758  void operator()(const T& edge)
759  {
760  subdomains[subdomain].insert(edge.target());
761  // If we (the neighbouring vertex of the aggregate)
762  // are not isolated, add the aggregate we belong to
763  // to the same subdomain using the OneOverlapAdder
764  if(aggregates[edge.target()]!=AggregatesMap::ISOLATED) {
765  assert(aggregates[edge.target()]!=aggregate);
766  typename AggregatesMap::VertexList vlist;
767  aggregates.template breadthFirstSearch<true,false>(edge.target(), aggregate,
768  graph, vlist, adder, adder,
769  visitedMap);
770  }
771  }
772 
773  int setAggregate(const AggregateDescriptor& aggregate_)
774  {
775  adder.setAggregate(aggregate_);
776  aggregate=aggregate_;
777  return ++subdomain;
778  }
779  int noSubdomains() const
780  {
781  return subdomain+1;
782  }
783 
784  private:
785  AggregateDescriptor aggregate;
786  Vector& subdomains;
787  int subdomain;
788  const AggregatesMap& aggregates;
789  VertexAdder adder;
790  const MatrixGraph<const M>& graph;
791  VM& visitedMap;
792  };
793 
794  void createPairDomains(const MatrixGraph<const M>& graph)
795  {
796  typedef typename MatrixGraph<const M>::ConstVertexIterator VIter;
797  typedef typename MatrixGraph<const M>::ConstEdgeIterator EIter;
798  typedef typename M::size_type size_type;
799 
800  std::set<std::pair<size_type,size_type> > pairs;
801  int total=0;
802  for(VIter v=graph.begin(), ve=graph.end(); ve != v; ++v)
803  for(EIter e = v.begin(), ee=v.end(); ee!=e; ++e)
804  {
805  ++total;
806  if(e.source()<e.target())
807  pairs.insert(std::make_pair(e.source(),e.target()));
808  else
809  pairs.insert(std::make_pair(e.target(),e.source()));
810  }
811 
812 
813  subdomains.resize(pairs.size());
814  Dune::dinfo <<std::endl<< "Created "<<pairs.size()<<" ("<<total<<") pair domains"<<std::endl<<std::endl;
815  typedef typename std::set<std::pair<size_type,size_type> >::const_iterator SIter;
816  typename Vector::iterator subdomain=subdomains.begin();
817 
818  for(SIter s=pairs.begin(), se =pairs.end(); se!=s; ++s)
819  {
820  subdomain->insert(s->first);
821  subdomain->insert(s->second);
822  ++subdomain;
823  }
824  std::size_t minsize=10000;
825  std::size_t maxsize=0;
826  int sum=0;
827  for(typename Vector::size_type i=0; i < subdomains.size(); ++i) {
828  sum+=subdomains[i].size();
829  minsize=std::min(minsize, subdomains[i].size());
830  maxsize=std::max(maxsize, subdomains[i].size());
831  }
832  Dune::dinfo<<"Subdomain size: min="<<minsize<<" max="<<maxsize<<" avg="<<(sum/subdomains.size())
833  <<" no="<<subdomains.size()<<std::endl;
834  }
835 
836  template<class Visitor>
837  void createSubdomains(const M& matrix, const MatrixGraph<const M>& graph,
838  const AggregatesMap& amap, Visitor& overlapVisitor,
839  IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap>& visitedMap )
840  {
841  // count number ag aggregates. We assume that the
842  // aggregates are numbered consecutively from 0 except
843  // for the isolated ones. All isolated vertices form
844  // one aggregate, here.
845  int isolated=0;
846  AggregateDescriptor maxAggregate=0;
847 
848  for(std::size_t i=0; i < amap.noVertices(); ++i)
849  if(amap[i]==AggregatesMap::ISOLATED)
850  isolated++;
851  else
852  maxAggregate = std::max(maxAggregate, amap[i]);
853 
854  subdomains.resize(maxAggregate+1+isolated);
855 
856  // reset the subdomains
857  for(typename Vector::size_type i=0; i < subdomains.size(); ++i)
858  subdomains[i].clear();
859 
860  // Create the subdomains from the aggregates mapping.
861  // For each aggregate we mark all entries and the
862  // neighbouring vertices as belonging to the same subdomain
863  VertexAdder aggregateVisitor(subdomains, amap);
864 
865  for(VertexDescriptor i=0; i < amap.noVertices(); ++i)
866  if(!get(visitedMap, i)) {
867  AggregateDescriptor aggregate=amap[i];
868 
869  if(amap[i]==AggregatesMap::ISOLATED) {
870  // isolated vertex gets its own aggregate
871  subdomains.push_back(Subdomain());
872  aggregate=subdomains.size()-1;
873  }
874  overlapVisitor.setAggregate(aggregate);
875  aggregateVisitor.setAggregate(aggregate);
876  subdomains[aggregate].insert(i);
877  typename AggregatesMap::VertexList vlist;
878  amap.template breadthFirstSearch<false,false>(i, aggregate, graph, vlist, aggregateVisitor,
879  overlapVisitor, visitedMap);
880  }
881 
882  std::size_t minsize=10000;
883  std::size_t maxsize=0;
884  int sum=0;
885  for(typename Vector::size_type i=0; i < subdomains.size(); ++i) {
886  sum+=subdomains[i].size();
887  minsize=std::min(minsize, subdomains[i].size());
888  maxsize=std::max(maxsize, subdomains[i].size());
889  }
890  Dune::dinfo<<"Subdomain size: min="<<minsize<<" max="<<maxsize<<" avg="<<(sum/subdomains.size())
891  <<" no="<<subdomains.size()<<" isolated="<<isolated<<std::endl;
892 
893 
894 
895  }
896  Vector subdomains;
897  };
898 
899 
900  template<class M, class X, class TM, class TS, class TA>
901  struct ConstructionTraits<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
902  {
903  typedef ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> > Arguments;
904 
905  static inline SeqOverlappingSchwarz<M,X,TM,TS,TA>* construct(Arguments& args)
906  {
907  return new SeqOverlappingSchwarz<M,X,TM,TS,TA>(args.getMatrix(),
908  args.getSubDomains(),
909  args.getArgs().relaxationFactor,
910  args.getArgs().onthefly);
911  }
912 
913  static void deconstruct(SeqOverlappingSchwarz<M,X,TM,TS,TA>* schwarz)
914  {
915  delete schwarz;
916  }
917  };
918 
919 
920  } // namespace Amg
921 } // namespace Dune
922 
923 
924 
925 #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:269
Nonoverlapping parallel preconditioner.
Definition: novlpschwarz.hh:266
A parallel SSOR preconditioner.
Definition: schwarz.hh:166
Sequential ILU0 preconditioner.
Definition: preconditioners.hh:652
Sequential ILU preconditioner.
Definition: preconditioners.hh:504
Sequential ILU(n) preconditioner.
Definition: preconditioners.hh:739
The sequential jacobian preconditioner.
Definition: preconditioners.hh:416
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:225
Sequential SSOR preconditioner.
Definition: preconditioners.hh:135
Helper classes for the construction of classes without empty constructor.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
constexpr GeometryType none(unsigned int dim)
Returns a GeometryType representing a singular of dimension dim.
Definition: type.hh:714
constexpr GeometryType vertex
GeometryType representing a vertex.
Definition: type.hh:727
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:454
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:442
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:476
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:430
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:248
Dune namespace.
Definition: alignedallocator.hh:10
Define general preconditioner interface.
The default class for the smoother arguments.
Definition: smoother.hh:36
Helper class for applying the smoothers.
Definition: smoother.hh:418
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.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)