3#ifndef DUNE_AMGSMOOTHER_HH
4#define DUNE_AMGSMOOTHER_HH
9#include <dune/istl/schwarz.hh>
10#include <dune/istl/novlpschwarz.hh>
12#include <dune/common/propertymap.hh>
70 template<
class X,
class Y>
77 template<
class X,
class Y,
class C,
class T>
79 :
public SmootherTraits<T>
82 template<
class C,
class T>
83 struct SmootherTraits<NonoverlappingBlockPreconditioner<C,T> >
84 :
public SmootherTraits<T>
93 typedef typename T::matrix_type Matrix;
103 void setMatrix(
const Matrix& matrix)
107 virtual void setMatrix(
const Matrix& matrix,
const AggregatesMap& amap)
114 const Matrix& getMatrix()
const
125 void setComm(T1& comm)
130 const SequentialInformation& getComm()
141 const Matrix* matrix_;
144 SequentialInformation comm_;
148 struct ConstructionArgs
152 template<
class T,
class C=SequentialInformation>
153 class DefaultParallelConstructionArgs
154 :
public ConstructionArgs<T>
157 virtual ~DefaultParallelConstructionArgs()
160 void setComm(
const C& comm)
165 const C& getComm()
const
174 template<
class X,
class Y>
175 class DefaultConstructionArgs<Richardson<X,Y>>
177 typedef Richardson<X,Y> T;
179 typedef typename SmootherTraits<T>::Arguments SmootherArgs;
182 virtual ~DefaultConstructionArgs()
185 template <
class... Args>
186 void setMatrix(
const Args&...)
189 void setArgs(
const SmootherArgs& args)
195 void setComm(T1& comm)
200 const SequentialInformation& getComm()
205 const SmootherArgs getArgs()
const
211 const SmootherArgs* args_;
212 SequentialInformation comm_;
218 class ConstructionTraits;
223 template<
class M,
class X,
class Y,
int l>
230 return std::make_shared<SeqSSOR<M,X,Y,l>>
239 template<
class M,
class X,
class Y,
int l>
246 return std::make_shared<SeqSOR<M,X,Y,l>>
255 template<
class M,
class X,
class Y,
int l>
262 return std::make_shared<SeqJac<M,X,Y,l>>
270 template<
class X,
class Y>
275 static inline std::shared_ptr<Richardson<X,Y>>
construct(Arguments& args)
277 return std::make_shared<Richardson<X,Y>>
278 (args.getArgs().relaxationFactor);
287 template<
class M,
class X,
class Y>
294 return std::make_shared<SeqILU0<M,X,Y>>
301 template<
class M,
class X,
class Y>
302 class ConstructionArgs<
SeqILUn<M,X,Y> >
306 ConstructionArgs(
int n=1)
327 template<
class M,
class X,
class Y>
330 typedef ConstructionArgs<SeqILUn<M,X,Y> >
Arguments;
332 static inline std::shared_ptr<SeqILUn<M,X,Y>>
construct(Arguments& args)
334 return std::make_shared<SeqILUn<M,X,Y>>
335 (args.getMatrix(), args.getN(), args.getArgs().relaxationFactor);
341 template<
class M,
class X,
class Y>
342 class ConstructionArgs<
SeqILU<M,X,Y> >
346 ConstructionArgs(
int n=0)
368 template<
class M,
class X,
class Y>
371 typedef ConstructionArgs<SeqILU<M,X,Y> >
Arguments;
373 static inline std::shared_ptr<SeqILU<M,X,Y>>
construct(Arguments& args)
375 return std::make_shared<SeqILU<M,X,Y>>
376 (args.getMatrix(), args.getN(), args.getArgs().relaxationFactor);
383 template<
class M,
class X,
class Y,
class C>
386 typedef DefaultParallelConstructionArgs<M,C>
Arguments;
388 static inline std::shared_ptr<ParSSOR<M,X,Y,C>>
construct(Arguments& args)
390 return std::make_shared<ParSSOR<M,X,Y,C>>
391 (args.getMatrix(), args.getArgs().iterations,
392 args.getArgs().relaxationFactor, args.getComm());
396 template<
class X,
class Y,
class C,
class T>
399 typedef DefaultParallelConstructionArgs<T,C>
Arguments;
401 static inline std::shared_ptr<BlockPreconditioner<X,Y,C,T>>
construct(
Arguments& args)
403 auto seqPrec = SeqConstructionTraits::construct(args);
404 return std::make_shared<BlockPreconditioner<X,Y,C,T>> (seqPrec, args.getComm());
408 template<
class C,
class T>
409 struct ConstructionTraits<NonoverlappingBlockPreconditioner<C,T> >
411 typedef DefaultParallelConstructionArgs<T,C>
Arguments;
412 typedef ConstructionTraits<T> SeqConstructionTraits;
413 static inline std::shared_ptr<NonoverlappingBlockPreconditioner<C,T>>
construct(
Arguments& args)
415 auto seqPrec = SeqConstructionTraits::construct(args);
416 return std::make_shared<NonoverlappingBlockPreconditioner<C,T>> (seqPrec, args.getComm());
434 typedef typename Smoother::range_type Range;
435 typedef typename Smoother::domain_type Domain;
444 static void preSmooth(Smoother& smoother, Domain& v,
const Range& d)
456 static void postSmooth(Smoother& smoother, Domain& v,
const Range& d)
467 template<
typename LevelContext>
468 void presmooth(LevelContext& levelContext,
size_t steps)
470 for(std::size_t i=0; i < steps; ++i) {
473 ::preSmooth(*levelContext.smoother, *levelContext.lhs,
476 *levelContext.update += *levelContext.lhs;
479 levelContext.matrix->applyscaleadd(-1, *levelContext.lhs, *levelContext.rhs);
480 levelContext.pinfo->project(*levelContext.rhs);
489 template<
typename LevelContext>
492 for(std::size_t i=0; i < steps; ++i) {
494 levelContext.matrix->applyscaleadd(-1, *levelContext.lhs,
497 levelContext.pinfo->project(*levelContext.rhs);
499 ::postSmooth(*levelContext.smoother, *levelContext.lhs, *levelContext.rhs);
501 *levelContext.update += *levelContext.lhs;
505 template<
class M,
class X,
class Y,
int l>
506 struct SmootherApplier<
SeqSOR<M,X,Y,l> >
509 typedef typename Smoother::range_type Range;
510 typedef typename Smoother::domain_type Domain;
512 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
514 smoother.template apply<true>(v,d);
518 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
520 smoother.template apply<false>(v,d);
524 template<
class M,
class X,
class Y,
class C,
int l>
525 struct SmootherApplier<BlockPreconditioner<X,Y,C,SeqSOR<M,X,Y,l> > >
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;
531 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
533 smoother.template apply<true>(v,d);
537 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
539 smoother.template apply<false>(v,d);
543 template<
class M,
class X,
class Y,
class C,
int l>
544 struct SmootherApplier<NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y,l> > >
546 typedef NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y,l> > Smoother;
547 typedef typename Smoother::range_type Range;
548 typedef typename Smoother::domain_type Domain;
550 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
552 smoother.template apply<true>(v,d);
556 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
558 smoother.template apply<false>(v,d);
565 template<
class M,
class X,
class MO,
class MS,
class A>
566 class SeqOverlappingSchwarz;
568 struct MultiplicativeSchwarzMode;
572 template<
class M,
class X,
class MS,
class TA>
573 struct SmootherApplier<SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,
576 typedef SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,MS,TA> Smoother;
577 typedef typename Smoother::range_type Range;
578 typedef typename Smoother::domain_type Domain;
580 static void preSmooth(Smoother& smoother, Domain& v,
const Range& d)
582 smoother.template apply<true>(v,d);
586 static void postSmooth(Smoother& smoother, Domain& v,
const Range& d)
588 smoother.template apply<false>(v,d);
597 struct SeqOverlappingSchwarzSmootherArgs
598 :
public DefaultSmootherArgs<T>
605 SeqOverlappingSchwarzSmootherArgs(Overlap overlap_=
vertex,
606 bool onthefly_=
false)
607 : overlap(overlap_), onthefly(onthefly_)
611 template<
class M,
class X,
class TM,
class TS,
class TA>
612 struct SmootherTraits<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
614 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> Arguments;
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> >
621 typedef DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> > Father;
628 typedef typename Vector::value_type Subdomain;
630 virtual void setMatrix(
const M& matrix,
const AggregatesMap& amap)
632 Father::setMatrix(matrix);
634 std::vector<bool> visited(amap.noVertices(),
false);
635 typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
636 VisitedMapType visitedMap(visited.begin());
638 MatrixGraph<const M> graph(matrix);
640 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
642 switch(Father::getArgs().overlap) {
645 VertexAdder visitor(subdomains, amap);
646 createSubdomains(matrix, graph, amap, visitor, visitedMap);
649 case SmootherArgs::pairwise :
651 createPairDomains(graph);
654 case SmootherArgs::aggregate :
656 AggregateAdder<VisitedMapType> visitor(subdomains, amap, graph, visitedMap);
657 createSubdomains(matrix, graph, amap, visitor, visitedMap);
662 createSubdomains(matrix, graph, amap, visitor, visitedMap);
665 DUNE_THROW(NotImplemented,
"This overlapping scheme is not supported!");
668 void setMatrix(
const M& matrix)
670 Father::setMatrix(matrix);
673 AggregatesMap amap(matrix.N());
674 VertexDescriptor v=0;
675 for(
typename AggregatesMap::iterator iter=amap.begin();
676 iter!=amap.end(); ++iter)
679 std::vector<bool> visited(amap.noVertices(),
false);
680 typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
681 VisitedMapType visitedMap(visited.begin());
683 MatrixGraph<const M> graph(matrix);
685 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
687 switch(Father::getArgs().overlap) {
690 VertexAdder visitor(subdomains, amap);
691 createSubdomains(matrix, graph, amap, visitor, visitedMap);
694 case SmootherArgs::aggregate :
696 DUNE_THROW(NotImplemented,
"Aggregate overlap is not supported yet");
703 case SmootherArgs::pairwise :
705 createPairDomains(graph);
710 createSubdomains(matrix, graph, amap, visitor, visitedMap);
715 const Vector& getSubDomains()
723 VertexAdder(Vector& subdomains_,
const AggregatesMap& aggregates_)
724 : subdomains(subdomains_),
max(-1), subdomain(-1), aggregates(aggregates_)
727 void operator()(
const T& edge)
730 subdomains[subdomain].insert(edge.target());
732 int setAggregate(
const AggregateDescriptor& aggregate_)
734 subdomain=aggregate_;
738 int noSubdomains()
const
744 AggregateDescriptor
max;
745 AggregateDescriptor subdomain;
746 const AggregatesMap& aggregates;
751 void operator()(
const T& edge)
753 int setAggregate(
const AggregateDescriptor& aggregate_)
757 int noSubdomains()
const
764 struct AggregateAdder
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_)
772 void operator()(
const T& edge)
774 subdomains[subdomain].insert(edge.target());
779 assert(aggregates[edge.target()]!=aggregate);
781 aggregates.template breadthFirstSearch<true,false>(edge.target(), aggregate,
782 graph, vlist, adder, adder,
787 int setAggregate(
const AggregateDescriptor& aggregate_)
789 adder.setAggregate(aggregate_);
790 aggregate=aggregate_;
793 int noSubdomains()
const
799 AggregateDescriptor aggregate;
802 const AggregatesMap& aggregates;
804 const MatrixGraph<const M>& graph;
808 void createPairDomains(
const MatrixGraph<const M>& graph)
812 typedef typename M::size_type size_type;
814 std::set<std::pair<size_type,size_type> > pairs;
816 for(VIter v=graph.begin(), ve=graph.end(); ve != v; ++v)
817 for(EIter e = v.begin(), ee=v.end(); ee!=e; ++e)
820 if(e.source()<e.target())
821 pairs.insert(std::make_pair(e.source(),e.target()));
823 pairs.insert(std::make_pair(e.target(),e.source()));
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();
832 for(SIter s=pairs.begin(), se =pairs.end(); se!=s; ++s)
834 subdomain->insert(s->first);
835 subdomain->insert(s->second);
838 std::size_t minsize=10000;
839 std::size_t maxsize=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());
846 Dune::dinfo<<
"Subdomain size: min="<<minsize<<
" max="<<maxsize<<
" avg="<<(sum/subdomains.size())
847 <<
" no="<<subdomains.size()<<std::endl;
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 )
860 AggregateDescriptor maxAggregate=0;
862 for(std::size_t i=0; i < amap.noVertices(); ++i)
866 maxAggregate =
std::max(maxAggregate, amap[i]);
868 subdomains.resize(maxAggregate+1+isolated);
871 for(
typename Vector::size_type i=0; i < subdomains.size(); ++i)
872 subdomains[i].clear();
877 VertexAdder aggregateVisitor(subdomains, amap);
879 for(VertexDescriptor i=0; i < amap.noVertices(); ++i)
880 if(!get(visitedMap, i)) {
881 AggregateDescriptor aggregate=amap[i];
885 subdomains.push_back(Subdomain());
886 aggregate=subdomains.size()-1;
888 overlapVisitor.setAggregate(aggregate);
889 aggregateVisitor.setAggregate(aggregate);
890 subdomains[aggregate].insert(i);
892 amap.template breadthFirstSearch<false,false>(i, aggregate, graph, vlist, aggregateVisitor,
893 overlapVisitor, visitedMap);
896 std::size_t minsize=10000;
897 std::size_t maxsize=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());
904 Dune::dinfo<<
"Subdomain size: min="<<minsize<<
" max="<<maxsize<<
" avg="<<(sum/subdomains.size())
905 <<
" no="<<subdomains.size()<<
" isolated="<<isolated<<std::endl;
914 template<
class M,
class X,
class TM,
class TS,
class TA>
915 struct ConstructionTraits<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
917 typedef ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
Arguments;
919 static inline std::shared_ptr<SeqOverlappingSchwarz<M,X,TM,TS,TA>>
construct(
Arguments& args)
921 return std::make_shared<SeqOverlappingSchwarz<M,X,TM,TS,TA>>
923 args.getSubDomains(),
924 args.getArgs().relaxationFactor,
925 args.getArgs().onthefly);
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
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:784
constexpr GeometryType vertex
GeometryType representing a vertex.
Definition: type.hh:797
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.