Dune Core Modules (unstable)

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