DUNE PDELab (2.7)

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>
12#include <dune/common/propertymap.hh>
13#include <dune/common/unused.hh>
14
15namespace Dune
16{
17 namespace Amg
18 {
19
35 template<class T>
37 {
42
51
57 {}
58 };
59
63 template<class T>
65 {
67
68 };
69
70 template<class X, class Y>
71 struct SmootherTraits<Richardson<X,Y>>
72 {
74
75 };
76
77 template<class X, class Y, class C, class T>
78 struct SmootherTraits<BlockPreconditioner<X,Y,C,T> >
79 : public SmootherTraits<T>
80 {};
81
82 template<class C, class T>
83 struct SmootherTraits<NonoverlappingBlockPreconditioner<C,T> >
84 : public SmootherTraits<T>
85 {};
86
90 template<class T>
92 {
93 typedef typename T::matrix_type Matrix;
94
96
98
99 public:
101 {}
102
103 void setMatrix(const Matrix& matrix)
104 {
105 matrix_=&matrix;
106 }
107 virtual void setMatrix(const Matrix& matrix, const AggregatesMap& amap)
108 {
110 setMatrix(matrix);
111 }
112
113
114 const Matrix& getMatrix() const
115 {
116 return *matrix_;
117 }
118
119 void setArgs(const SmootherArgs& args)
120 {
121 args_=&args;
122 }
123
124 template<class T1>
125 void setComm(T1& comm)
126 {
128 }
129
130 const SequentialInformation& getComm()
131 {
132 return comm_;
133 }
134
135 const SmootherArgs getArgs() const
136 {
137 return *args_;
138 }
139
140 protected:
141 const Matrix* matrix_;
142 private:
143 const SmootherArgs* args_;
144 SequentialInformation comm_;
145 };
146
147 template<class T>
148 struct ConstructionArgs
149 : public DefaultConstructionArgs<T>
150 {};
151
152 template<class T, class C=SequentialInformation>
153 class DefaultParallelConstructionArgs
154 : public ConstructionArgs<T>
155 {
156 public:
157 virtual ~DefaultParallelConstructionArgs()
158 {}
159
160 void setComm(const C& comm)
161 {
162 comm_ = &comm;
163 }
164
165 const C& getComm() const
166 {
167 return *comm_;
168 }
169 private:
170 const C* comm_;
171 };
172
173
174 template<class X, class Y>
175 class DefaultConstructionArgs<Richardson<X,Y>>
176 {
177 typedef Richardson<X,Y> T;
178
179 typedef typename SmootherTraits<T>::Arguments SmootherArgs;
180
181 public:
182 virtual ~DefaultConstructionArgs()
183 {}
184
185 template <class... Args>
186 void setMatrix(const Args&...)
187 {}
188
189 void setArgs(const SmootherArgs& args)
190 {
191 args_=&args;
192 }
193
194 template<class T1>
195 void setComm(T1& comm)
196 {
198 }
199
200 const SequentialInformation& getComm()
201 {
202 return comm_;
203 }
204
205 const SmootherArgs getArgs() const
206 {
207 return *args_;
208 }
209
210 private:
211 const SmootherArgs* args_;
212 SequentialInformation comm_;
213 };
214
215
216
217 template<class T>
218 class ConstructionTraits;
219
223 template<class M, class X, class Y, int l>
224 struct ConstructionTraits<SeqSSOR<M,X,Y,l> >
225 {
227
228 static inline std::shared_ptr<SeqSSOR<M,X,Y,l>> construct(Arguments& args)
229 {
230 return std::make_shared<SeqSSOR<M,X,Y,l>>
231 (args.getMatrix(), args.getArgs().iterations, args.getArgs().relaxationFactor);
232 }
233 };
234
235
239 template<class M, class X, class Y, int l>
240 struct ConstructionTraits<SeqSOR<M,X,Y,l> >
241 {
243
244 static inline std::shared_ptr<SeqSOR<M,X,Y,l>> construct(Arguments& args)
245 {
246 return std::make_shared<SeqSOR<M,X,Y,l>>
247 (args.getMatrix(), args.getArgs().iterations, args.getArgs().relaxationFactor);
248 }
249 };
250
251
255 template<class M, class X, class Y, int l>
256 struct ConstructionTraits<SeqJac<M,X,Y,l> >
257 {
259
260 static inline std::shared_ptr<SeqJac<M,X,Y,l>> construct(Arguments& args)
261 {
262 return std::make_shared<SeqJac<M,X,Y,l>>
263 (args.getMatrix(), args.getArgs().iterations, args.getArgs().relaxationFactor);
264 }
265 };
266
270 template<class X, class Y>
272 {
274
275 static inline std::shared_ptr<Richardson<X,Y>> construct(Arguments& args)
276 {
277 return std::make_shared<Richardson<X,Y>>
278 (args.getArgs().relaxationFactor);
279 }
280 };
281
282
283DUNE_NO_DEPRECATED_BEGIN // for deprecated SeqILU0
287 template<class M, class X, class Y>
289 {
291
292 static inline std::shared_ptr<SeqILU0<M,X,Y>> construct(Arguments& args)
293 {
294 return std::make_shared<SeqILU0<M,X,Y>>
295 (args.getMatrix(), args.getArgs().relaxationFactor);
296 }
297 };
298DUNE_NO_DEPRECATED_END // for deprecated SeqILU0
299
300DUNE_NO_DEPRECATED_BEGIN // for deprecated SeqILUn
301 template<class M, class X, class Y>
302 class ConstructionArgs<SeqILUn<M,X,Y> >
303 : public DefaultConstructionArgs<SeqILUn<M,X,Y> >
304 {
305 public:
306 ConstructionArgs(int n=1)
307 : n_(n)
308 {}
309
310 void setN(int n)
311 {
312 n_ = n;
313 }
314 int getN()
315 {
316 return n_;
317 }
318
319 private:
320 int n_;
321 };
322
323
327 template<class M, class X, class Y>
329 {
330 typedef ConstructionArgs<SeqILUn<M,X,Y> > Arguments;
331
332 static inline std::shared_ptr<SeqILUn<M,X,Y>> construct(Arguments& args)
333 {
334 return std::make_shared<SeqILUn<M,X,Y>>
335 (args.getMatrix(), args.getN(), args.getArgs().relaxationFactor);
336 }
337 };
338DUNE_NO_DEPRECATED_END // for deprecated SeqILUn
339
340
341 template<class M, class X, class Y>
342 class ConstructionArgs<SeqILU<M,X,Y> >
343 : public DefaultConstructionArgs<SeqILU<M,X,Y> >
344 {
345 public:
346 ConstructionArgs(int n=0)
347 : n_(n)
348 {}
349
350 void setN(int n)
351 {
352 n_ = n;
353 }
354
355 int getN()
356 {
357 return n_;
358 }
359
360 private:
361 int n_;
362 };
363
364
368 template<class M, class X, class Y>
370 {
371 typedef ConstructionArgs<SeqILU<M,X,Y> > Arguments;
372
373 static inline std::shared_ptr<SeqILU<M,X,Y>> construct(Arguments& args)
374 {
375 return std::make_shared<SeqILU<M,X,Y>>
376 (args.getMatrix(), args.getN(), args.getArgs().relaxationFactor);
377 }
378 };
379
383 template<class M, class X, class Y, class C>
384 struct ConstructionTraits<ParSSOR<M,X,Y,C> >
385 {
386 typedef DefaultParallelConstructionArgs<M,C> Arguments;
387
388 static inline std::shared_ptr<ParSSOR<M,X,Y,C>> construct(Arguments& args)
389 {
390 return std::make_shared<ParSSOR<M,X,Y,C>>
391 (args.getMatrix(), args.getArgs().iterations,
392 args.getArgs().relaxationFactor, args.getComm());
393 }
394 };
395
396 template<class X, class Y, class C, class T>
398 {
399 typedef DefaultParallelConstructionArgs<T,C> Arguments;
400 typedef ConstructionTraits<T> SeqConstructionTraits;
401 static inline std::shared_ptr<BlockPreconditioner<X,Y,C,T>> construct(Arguments& args)
402 {
403 auto seqPrec = SeqConstructionTraits::construct(args);
404 return std::make_shared<BlockPreconditioner<X,Y,C,T>> (seqPrec, args.getComm());
405 }
406 };
407
408 template<class C, class T>
409 struct ConstructionTraits<NonoverlappingBlockPreconditioner<C,T> >
410 {
411 typedef DefaultParallelConstructionArgs<T,C> Arguments;
412 typedef ConstructionTraits<T> SeqConstructionTraits;
413 static inline std::shared_ptr<NonoverlappingBlockPreconditioner<C,T>> construct(Arguments& args)
414 {
415 auto seqPrec = SeqConstructionTraits::construct(args);
416 return std::make_shared<NonoverlappingBlockPreconditioner<C,T>> (seqPrec, args.getComm());
417 }
418 };
419
430 template<class T>
432 {
433 typedef T Smoother;
434 typedef typename Smoother::range_type Range;
435 typedef typename Smoother::domain_type Domain;
436
444 static void preSmooth(Smoother& smoother, Domain& v, const Range& d)
445 {
446 smoother.apply(v,d);
447 }
448
456 static void postSmooth(Smoother& smoother, Domain& v, const Range& d)
457 {
458 smoother.apply(v,d);
459 }
460 };
461
467 template<typename LevelContext>
468 void presmooth(LevelContext& levelContext, size_t steps)
469 {
470 for(std::size_t i=0; i < steps; ++i) {
471 *levelContext.lhs=0;
473 ::preSmooth(*levelContext.smoother, *levelContext.lhs,
474 *levelContext.rhs);
475 // Accumulate update
476 *levelContext.update += *levelContext.lhs;
477
478 // update defect
479 levelContext.matrix->applyscaleadd(-1, *levelContext.lhs, *levelContext.rhs);
480 levelContext.pinfo->project(*levelContext.rhs);
481 }
482 }
483
489 template<typename LevelContext>
490 void postsmooth(LevelContext& levelContext, size_t steps)
491 {
492 for(std::size_t i=0; i < steps; ++i) {
493 // update defect
494 levelContext.matrix->applyscaleadd(-1, *levelContext.lhs,
495 *levelContext.rhs);
496 *levelContext.lhs=0;
497 levelContext.pinfo->project(*levelContext.rhs);
499 ::postSmooth(*levelContext.smoother, *levelContext.lhs, *levelContext.rhs);
500 // Accumulate update
501 *levelContext.update += *levelContext.lhs;
502 }
503 }
504
505 template<class M, class X, class Y, int l>
506 struct SmootherApplier<SeqSOR<M,X,Y,l> >
507 {
508 typedef SeqSOR<M,X,Y,l> Smoother;
509 typedef typename Smoother::range_type Range;
510 typedef typename Smoother::domain_type Domain;
511
512 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
513 {
514 smoother.template apply<true>(v,d);
515 }
516
517
518 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
519 {
520 smoother.template apply<false>(v,d);
521 }
522 };
523
524 template<class M, class X, class Y, class C, int l>
525 struct SmootherApplier<BlockPreconditioner<X,Y,C,SeqSOR<M,X,Y,l> > >
526 {
527 typedef BlockPreconditioner<X,Y,C,SeqSOR<M,X,Y,l> > Smoother;
528 typedef typename Smoother::range_type Range;
529 typedef typename Smoother::domain_type Domain;
530
531 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
532 {
533 smoother.template apply<true>(v,d);
534 }
535
536
537 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
538 {
539 smoother.template apply<false>(v,d);
540 }
541 };
542
543 template<class M, class X, class Y, class C, int l>
544 struct SmootherApplier<NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y,l> > >
545 {
546 typedef NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y,l> > Smoother;
547 typedef typename Smoother::range_type Range;
548 typedef typename Smoother::domain_type Domain;
549
550 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
551 {
552 smoother.template apply<true>(v,d);
553 }
554
555
556 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
557 {
558 smoother.template apply<false>(v,d);
559 }
560 };
561
562 } // end namespace Amg
563
564 // forward declarations
565 template<class M, class X, class MO, class MS, class A>
566 class SeqOverlappingSchwarz;
567
568 struct MultiplicativeSchwarzMode;
569
570 namespace Amg
571 {
572 template<class M, class X, class MS, class TA>
573 struct SmootherApplier<SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,
574 MS,TA> >
575 {
576 typedef SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,MS,TA> Smoother;
577 typedef typename Smoother::range_type Range;
578 typedef typename Smoother::domain_type Domain;
579
580 static void preSmooth(Smoother& smoother, Domain& v, const Range& d)
581 {
582 smoother.template apply<true>(v,d);
583 }
584
585
586 static void postSmooth(Smoother& smoother, Domain& v, const Range& d)
587 {
588 smoother.template apply<false>(v,d);
589
590 }
591 };
592
593 // template<class M, class X, class TM, class TA>
594 // class SeqOverlappingSchwarz;
595
596 template<class T>
597 struct SeqOverlappingSchwarzSmootherArgs
598 : public DefaultSmootherArgs<T>
599 {
600 enum Overlap {vertex, aggregate, pairwise, none};
601
602 Overlap overlap;
603 bool onthefly;
604
605 SeqOverlappingSchwarzSmootherArgs(Overlap overlap_=vertex,
606 bool onthefly_=false)
607 : overlap(overlap_), onthefly(onthefly_)
608 {}
609 };
610
611 template<class M, class X, class TM, class TS, class TA>
612 struct SmootherTraits<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
613 {
614 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> Arguments;
615 };
616
617 template<class M, class X, class TM, class TS, class TA>
618 class ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
619 : public DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
620 {
621 typedef DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> > Father;
622
623 public:
624 typedef typename MatrixGraph<M>::VertexDescriptor VertexDescriptor;
625 typedef Dune::Amg::AggregatesMap<VertexDescriptor> AggregatesMap;
626 typedef typename AggregatesMap::AggregateDescriptor AggregateDescriptor;
628 typedef typename Vector::value_type Subdomain;
629
630 virtual void setMatrix(const M& matrix, const AggregatesMap& amap)
631 {
632 Father::setMatrix(matrix);
633
634 std::vector<bool> visited(amap.noVertices(), false);
635 typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
636 VisitedMapType visitedMap(visited.begin());
637
638 MatrixGraph<const M> graph(matrix);
639
640 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
641
642 switch(Father::getArgs().overlap) {
644 {
645 VertexAdder visitor(subdomains, amap);
646 createSubdomains(matrix, graph, amap, visitor, visitedMap);
647 }
648 break;
649 case SmootherArgs::pairwise :
650 {
651 createPairDomains(graph);
652 }
653 break;
654 case SmootherArgs::aggregate :
655 {
656 AggregateAdder<VisitedMapType> visitor(subdomains, amap, graph, visitedMap);
657 createSubdomains(matrix, graph, amap, visitor, visitedMap);
658 }
659 break;
660 case SmootherArgs::none :
661 NoneAdder visitor;
662 createSubdomains(matrix, graph, amap, visitor, visitedMap);
663 break;
664 default :
665 DUNE_THROW(NotImplemented, "This overlapping scheme is not supported!");
666 }
667 }
668 void setMatrix(const M& matrix)
669 {
670 Father::setMatrix(matrix);
671
672 /* Create aggregates map where each aggregate is just one vertex. */
673 AggregatesMap amap(matrix.N());
674 VertexDescriptor v=0;
675 for(typename AggregatesMap::iterator iter=amap.begin();
676 iter!=amap.end(); ++iter)
677 *iter=v++;
678
679 std::vector<bool> visited(amap.noVertices(), false);
680 typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
681 VisitedMapType visitedMap(visited.begin());
682
683 MatrixGraph<const M> graph(matrix);
684
685 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
686
687 switch(Father::getArgs().overlap) {
689 {
690 VertexAdder visitor(subdomains, amap);
691 createSubdomains(matrix, graph, amap, visitor, visitedMap);
692 }
693 break;
694 case SmootherArgs::aggregate :
695 {
696 DUNE_THROW(NotImplemented, "Aggregate overlap is not supported yet");
697 /*
698 AggregateAdder<VisitedMapType> visitor(subdomains, amap, graph, visitedMap);
699 createSubdomains(matrix, graph, amap, visitor, visitedMap);
700 */
701 }
702 break;
703 case SmootherArgs::pairwise :
704 {
705 createPairDomains(graph);
706 }
707 break;
708 case SmootherArgs::none :
709 NoneAdder visitor;
710 createSubdomains(matrix, graph, amap, visitor, visitedMap);
711
712 }
713 }
714
715 const Vector& getSubDomains()
716 {
717 return subdomains;
718 }
719
720 private:
721 struct VertexAdder
722 {
723 VertexAdder(Vector& subdomains_, const AggregatesMap& aggregates_)
724 : subdomains(subdomains_), max(-1), subdomain(-1), aggregates(aggregates_)
725 {}
726 template<class T>
727 void operator()(const T& edge)
728 {
729 if(aggregates[edge.target()]!=AggregatesMap::ISOLATED)
730 subdomains[subdomain].insert(edge.target());
731 }
732 int setAggregate(const AggregateDescriptor& aggregate_)
733 {
734 subdomain=aggregate_;
735 max = std::max(subdomain, aggregate_);
736 return subdomain;
737 }
738 int noSubdomains() const
739 {
740 return max+1;
741 }
742 private:
743 Vector& subdomains;
744 AggregateDescriptor max;
745 AggregateDescriptor subdomain;
746 const AggregatesMap& aggregates;
747 };
748 struct NoneAdder
749 {
750 template<class T>
751 void operator()(const T& edge)
752 {}
753 int setAggregate(const AggregateDescriptor& aggregate_)
754 {
755 return -1;
756 }
757 int noSubdomains() const
758 {
759 return -1;
760 }
761 };
762
763 template<class VM>
764 struct AggregateAdder
765 {
766 AggregateAdder(Vector& subdomains_, const AggregatesMap& aggregates_,
767 const MatrixGraph<const M>& graph_, VM& visitedMap_)
768 : subdomains(subdomains_), subdomain(-1), aggregates(aggregates_),
769 adder(subdomains_, aggregates_), graph(graph_), visitedMap(visitedMap_)
770 {}
771 template<class T>
772 void operator()(const T& edge)
773 {
774 subdomains[subdomain].insert(edge.target());
775 // If we (the neighbouring vertex of the aggregate)
776 // are not isolated, add the aggregate we belong to
777 // to the same subdomain using the OneOverlapAdder
778 if(aggregates[edge.target()]!=AggregatesMap::ISOLATED) {
779 assert(aggregates[edge.target()]!=aggregate);
780 typename AggregatesMap::VertexList vlist;
781 aggregates.template breadthFirstSearch<true,false>(edge.target(), aggregate,
782 graph, vlist, adder, adder,
783 visitedMap);
784 }
785 }
786
787 int setAggregate(const AggregateDescriptor& aggregate_)
788 {
789 adder.setAggregate(aggregate_);
790 aggregate=aggregate_;
791 return ++subdomain;
792 }
793 int noSubdomains() const
794 {
795 return subdomain+1;
796 }
797
798 private:
799 AggregateDescriptor aggregate;
800 Vector& subdomains;
801 int subdomain;
802 const AggregatesMap& aggregates;
803 VertexAdder adder;
804 const MatrixGraph<const M>& graph;
805 VM& visitedMap;
806 };
807
808 void createPairDomains(const MatrixGraph<const M>& graph)
809 {
810 typedef typename MatrixGraph<const M>::ConstVertexIterator VIter;
811 typedef typename MatrixGraph<const M>::ConstEdgeIterator EIter;
812 typedef typename M::size_type size_type;
813
814 std::set<std::pair<size_type,size_type> > pairs;
815 int total=0;
816 for(VIter v=graph.begin(), ve=graph.end(); ve != v; ++v)
817 for(EIter e = v.begin(), ee=v.end(); ee!=e; ++e)
818 {
819 ++total;
820 if(e.source()<e.target())
821 pairs.insert(std::make_pair(e.source(),e.target()));
822 else
823 pairs.insert(std::make_pair(e.target(),e.source()));
824 }
825
826
827 subdomains.resize(pairs.size());
828 Dune::dinfo <<std::endl<< "Created "<<pairs.size()<<" ("<<total<<") pair domains"<<std::endl<<std::endl;
829 typedef typename std::set<std::pair<size_type,size_type> >::const_iterator SIter;
830 typename Vector::iterator subdomain=subdomains.begin();
831
832 for(SIter s=pairs.begin(), se =pairs.end(); se!=s; ++s)
833 {
834 subdomain->insert(s->first);
835 subdomain->insert(s->second);
836 ++subdomain;
837 }
838 std::size_t minsize=10000;
839 std::size_t maxsize=0;
840 int sum=0;
841 for(typename Vector::size_type i=0; i < subdomains.size(); ++i) {
842 sum+=subdomains[i].size();
843 minsize=std::min(minsize, subdomains[i].size());
844 maxsize=std::max(maxsize, subdomains[i].size());
845 }
846 Dune::dinfo<<"Subdomain size: min="<<minsize<<" max="<<maxsize<<" avg="<<(sum/subdomains.size())
847 <<" no="<<subdomains.size()<<std::endl;
848 }
849
850 template<class Visitor>
851 void createSubdomains(const M& matrix, const MatrixGraph<const M>& graph,
852 const AggregatesMap& amap, Visitor& overlapVisitor,
853 IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap>& visitedMap )
854 {
855 // count number ag aggregates. We assume that the
856 // aggregates are numbered consecutively from 0 except
857 // for the isolated ones. All isolated vertices form
858 // one aggregate, here.
859 int isolated=0;
860 AggregateDescriptor maxAggregate=0;
861
862 for(std::size_t i=0; i < amap.noVertices(); ++i)
863 if(amap[i]==AggregatesMap::ISOLATED)
864 isolated++;
865 else
866 maxAggregate = std::max(maxAggregate, amap[i]);
867
868 subdomains.resize(maxAggregate+1+isolated);
869
870 // reset the subdomains
871 for(typename Vector::size_type i=0; i < subdomains.size(); ++i)
872 subdomains[i].clear();
873
874 // Create the subdomains from the aggregates mapping.
875 // For each aggregate we mark all entries and the
876 // neighbouring vertices as belonging to the same subdomain
877 VertexAdder aggregateVisitor(subdomains, amap);
878
879 for(VertexDescriptor i=0; i < amap.noVertices(); ++i)
880 if(!get(visitedMap, i)) {
881 AggregateDescriptor aggregate=amap[i];
882
883 if(amap[i]==AggregatesMap::ISOLATED) {
884 // isolated vertex gets its own aggregate
885 subdomains.push_back(Subdomain());
886 aggregate=subdomains.size()-1;
887 }
888 overlapVisitor.setAggregate(aggregate);
889 aggregateVisitor.setAggregate(aggregate);
890 subdomains[aggregate].insert(i);
891 typename AggregatesMap::VertexList vlist;
892 amap.template breadthFirstSearch<false,false>(i, aggregate, graph, vlist, aggregateVisitor,
893 overlapVisitor, visitedMap);
894 }
895
896 std::size_t minsize=10000;
897 std::size_t maxsize=0;
898 int sum=0;
899 for(typename Vector::size_type i=0; i < subdomains.size(); ++i) {
900 sum+=subdomains[i].size();
901 minsize=std::min(minsize, subdomains[i].size());
902 maxsize=std::max(maxsize, subdomains[i].size());
903 }
904 Dune::dinfo<<"Subdomain size: min="<<minsize<<" max="<<maxsize<<" avg="<<(sum/subdomains.size())
905 <<" no="<<subdomains.size()<<" isolated="<<isolated<<std::endl;
906
907
908
909 }
910 Vector subdomains;
911 };
912
913
914 template<class M, class X, class TM, class TS, class TA>
915 struct ConstructionTraits<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
916 {
917 typedef ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> > Arguments;
918
919 static inline std::shared_ptr<SeqOverlappingSchwarz<M,X,TM,TS,TA>> construct(Arguments& args)
920 {
921 return std::make_shared<SeqOverlappingSchwarz<M,X,TM,TS,TA>>
922 (args.getMatrix(),
923 args.getSubDomains(),
924 args.getArgs().relaxationFactor,
925 args.getArgs().onthefly);
926 }
927 };
928
929
930 } // namespace Amg
931} // namespace Dune
932
933
934
935#endif
Provides classes for the Coloring process of AMG.
Class providing information about the mapping of the vertices onto aggregates.
Definition: aggregates.hh:556
Traits class for generically constructing non default constructable types.
Definition: construction.hh:38
Construction Arguments for the default smoothers.
Definition: smoother.hh:92
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:273
A parallel SSOR preconditioner.
Definition: schwarz.hh:170
Richardson preconditioner.
Definition: preconditioners.hh:875
Sequential ILU0 preconditioner.
Definition: preconditioners.hh:670
Sequential ILU preconditioner.
Definition: preconditioners.hh:500
Sequential ILU(n) preconditioner.
Definition: preconditioners.hh:778
The sequential jacobian preconditioner.
Definition: preconditioners.hh:390
std::vector< subdomain_type, typename std::allocator_traits< TA >::template rebind_alloc< subdomain_type > > subdomain_vector
The vector type containing the subdomain to row index mapping.
Definition: overlappingschwarz.hh:794
Sequential SOR preconditioner.
Definition: preconditioners.hh:249
Sequential SSOR preconditioner.
Definition: preconditioners.hh:138
Helper classes for the construction of classes without empty constructor.
Definition of the DUNE_DEPRECATED macro for the case that config.h is not available.
#define DUNE_NO_DEPRECATED_END
Ignore deprecation warnings (end)
Definition: deprecated.hh:198
#define DUNE_NO_DEPRECATED_BEGIN
Ignore deprecation warnings (start)
Definition: deprecated.hh:192
#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:785
constexpr GeometryType vertex
GeometryType representing a vertex.
Definition: type.hh:798
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:468
V AggregateDescriptor
The aggregate descriptor type.
Definition: aggregates.hh:576
static const V ISOLATED
Identifier of isolated vertices.
Definition: aggregates.hh:567
DefaultSmootherArgs()
Default constructor.
Definition: smoother.hh:55
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:456
static std::shared_ptr< T > construct(Arguments &args)
Construct an object with the specified arguments.
Definition: construction.hh:52
SLList< VertexDescriptor, Allocator > VertexList
The type of a single linked list of vertex descriptors.
Definition: aggregates.hh:588
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:490
RelaxationFactor relaxationFactor
The relaxation factor to use.
Definition: smoother.hh:50
static void preSmooth(Smoother &smoother, Domain &v, const Range &d)
apply pre smoothing in forward direction
Definition: smoother.hh:444
int iterations
The numbe of iterations to perform.
Definition: smoother.hh:46
T RelaxationFactor
The type of the relaxation factor.
Definition: smoother.hh:41
auto min(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::min()
Definition: defaults.hh:87
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
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:14
Define general preconditioner interface.
The default class for the smoother arguments.
Definition: smoother.hh:37
Helper class for applying the smoothers.
Definition: smoother.hh:432
Traits class for getting the attribute class of a smoother.
Definition: smoother.hh:65
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.111.3 (Jul 15, 22:36, 2024)