3#ifndef DUNE_AMGSMOOTHER_HH
4#define DUNE_AMGSMOOTHER_HH
9#include <dune/istl/schwarz.hh>
10#include <dune/istl/novlpschwarz.hh>
11#include <dune/common/propertymap.hh>
69 template<
class X,
class Y,
class C,
class T>
76 template<
class C,
class T>
79 typedef DefaultSmootherArgs<typename T::matrix_type::field_type> Arguments;
89 typedef typename T::matrix_type Matrix;
99 void setMatrix(
const Matrix& matrix)
103 virtual void setMatrix(
const Matrix& matrix,
const AggregatesMap& amap)
110 const Matrix& getMatrix()
const
121 void setComm(T1& comm)
126 const SequentialInformation& getComm()
137 const Matrix* matrix_;
140 SequentialInformation comm_;
144 struct ConstructionArgs
148 template<
class T,
class C=SequentialInformation>
149 class DefaultParallelConstructionArgs
150 :
public ConstructionArgs<T>
153 virtual ~DefaultParallelConstructionArgs()
156 void setComm(
const C& comm)
161 const C& getComm()
const
171 class ConstructionTraits;
176 template<
class M,
class X,
class Y,
int l>
198 template<
class M,
class X,
class Y,
int l>
218 template<
class M,
class X,
class Y,
int l>
240 template<
class M,
class X,
class Y>
258 template<
class M,
class X,
class Y>
259 class ConstructionArgs<
SeqILUn<M,X,Y> >
263 ConstructionArgs(
int n=1)
284 template<
class M,
class X,
class Y>
287 typedef ConstructionArgs<SeqILUn<M,X,Y> >
Arguments;
292 args.getArgs().relaxationFactor);
305 template<
class M,
class X,
class Y,
class C>
308 typedef DefaultParallelConstructionArgs<M,C>
Arguments;
313 args.getArgs().relaxationFactor,
322 template<
class X,
class Y,
class C,
class T>
325 typedef DefaultParallelConstructionArgs<T,C>
Arguments;
335 SeqConstructionTraits::deconstruct(
static_cast<T*
>(&bp->preconditioner));
341 template<
class C,
class T>
342 struct ConstructionTraits<NonoverlappingBlockPreconditioner<C,T> >
344 typedef DefaultParallelConstructionArgs<T,C>
Arguments;
345 typedef ConstructionTraits<T> SeqConstructionTraits;
348 return new NonoverlappingBlockPreconditioner<C,T>(*SeqConstructionTraits::construct(args),
352 static inline void deconstruct(NonoverlappingBlockPreconditioner<C,T>* bp)
354 SeqConstructionTraits::deconstruct(
static_cast<T*
>(&bp->preconditioner));
374 typedef typename Smoother::range_type Range;
375 typedef typename Smoother::domain_type Domain;
384 static void preSmooth(Smoother& smoother, Domain& v,
const Range& d)
396 static void postSmooth(Smoother& smoother, Domain& v,
const Range& d)
407 template<
typename LevelContext>
408 void presmooth(LevelContext& levelContext,
size_t steps)
410 for(std::size_t i=0; i < steps; ++i) {
413 ::preSmooth(*levelContext.smoother, *levelContext.lhs,
416 *levelContext.update += *levelContext.lhs;
419 levelContext.matrix->applyscaleadd(-1, *levelContext.lhs, *levelContext.rhs);
420 levelContext.pinfo->project(*levelContext.rhs);
429 template<
typename LevelContext>
432 for(std::size_t i=0; i < steps; ++i) {
434 levelContext.matrix->applyscaleadd(-1, *levelContext.lhs,
437 levelContext.pinfo->project(*levelContext.rhs);
439 ::postSmooth(*levelContext.smoother, *levelContext.lhs, *levelContext.rhs);
441 *levelContext.update += *levelContext.lhs;
445 template<
class M,
class X,
class Y,
int l>
446 struct SmootherApplier<
SeqSOR<M,X,Y,l> >
449 typedef typename Smoother::range_type Range;
450 typedef typename Smoother::domain_type Domain;
452 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
454 smoother.template apply<true>(v,d);
458 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
460 smoother.template apply<false>(v,d);
464 template<
class M,
class X,
class Y,
class C,
int l>
465 struct SmootherApplier<BlockPreconditioner<X,Y,C,SeqSOR<M,X,Y,l> > >
467 typedef BlockPreconditioner<X,Y,C,SeqSOR<M,X,Y,l> > Smoother;
468 typedef typename Smoother::range_type Range;
469 typedef typename Smoother::domain_type Domain;
471 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
473 smoother.template apply<true>(v,d);
477 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
479 smoother.template apply<false>(v,d);
483 template<
class M,
class X,
class Y,
class C,
int l>
484 struct SmootherApplier<NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y,l> > >
486 typedef NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y,l> > Smoother;
487 typedef typename Smoother::range_type Range;
488 typedef typename Smoother::domain_type Domain;
490 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
492 smoother.template apply<true>(v,d);
496 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
498 smoother.template apply<false>(v,d);
505 template<
class M,
class X,
class MO,
class MS,
class A>
506 class SeqOverlappingSchwarz;
508 struct MultiplicativeSchwarzMode;
512 template<
class M,
class X,
class MS,
class TA>
513 struct SmootherApplier<SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,
516 typedef SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,MS,TA> Smoother;
517 typedef typename Smoother::range_type Range;
518 typedef typename Smoother::domain_type Domain;
520 static void preSmooth(Smoother& smoother, Domain& v,
const Range& d)
522 smoother.template apply<true>(v,d);
526 static void postSmooth(Smoother& smoother, Domain& v,
const Range& d)
528 smoother.template apply<false>(v,d);
537 struct SeqOverlappingSchwarzSmootherArgs
538 :
public DefaultSmootherArgs<T>
540 enum Overlap {vertex, aggregate, pairwise, none};
545 SeqOverlappingSchwarzSmootherArgs(Overlap overlap_=vertex,
546 bool onthefly_=
false)
547 : overlap(overlap_), onthefly(onthefly_)
551 template<
class M,
class X,
class TM,
class TS,
class TA>
552 struct SmootherTraits<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
554 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> Arguments;
557 template<
class M,
class X,
class TM,
class TS,
class TA>
558 class ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
559 :
public DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
561 typedef DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> > Father;
568 typedef typename Vector::value_type Subdomain;
570 virtual void setMatrix(
const M& matrix,
const AggregatesMap& amap)
572 Father::setMatrix(matrix);
574 std::vector<bool> visited(amap.noVertices(),
false);
575 typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
576 VisitedMapType visitedMap(visited.begin());
578 MatrixGraph<const M> graph(matrix);
580 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
582 switch(Father::getArgs().overlap) {
583 case SmootherArgs::vertex :
585 VertexAdder visitor(subdomains, amap);
586 createSubdomains(matrix, graph, amap, visitor, visitedMap);
589 case SmootherArgs::pairwise :
591 createPairDomains(graph);
594 case SmootherArgs::aggregate :
596 AggregateAdder<VisitedMapType> visitor(subdomains, amap, graph, visitedMap);
597 createSubdomains(matrix, graph, amap, visitor, visitedMap);
600 case SmootherArgs::none :
602 createSubdomains(matrix, graph, amap, visitor, visitedMap);
605 DUNE_THROW(NotImplemented,
"This overlapping scheme is not supported!");
608 void setMatrix(
const M& matrix)
610 Father::setMatrix(matrix);
613 AggregatesMap amap(matrix.N());
614 VertexDescriptor v=0;
615 for(
typename AggregatesMap::iterator iter=amap.begin();
616 iter!=amap.end(); ++iter)
619 std::vector<bool> visited(amap.noVertices(),
false);
620 typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
621 VisitedMapType visitedMap(visited.begin());
623 MatrixGraph<const M> graph(matrix);
625 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
627 switch(Father::getArgs().overlap) {
628 case SmootherArgs::vertex :
630 VertexAdder visitor(subdomains, amap);
631 createSubdomains(matrix, graph, amap, visitor, visitedMap);
634 case SmootherArgs::aggregate :
636 DUNE_THROW(NotImplemented,
"Aggregate overlap is not supported yet");
643 case SmootherArgs::pairwise :
645 createPairDomains(graph);
648 case SmootherArgs::none :
650 createSubdomains(matrix, graph, amap, visitor, visitedMap);
655 const Vector& getSubDomains()
663 VertexAdder(Vector& subdomains_,
const AggregatesMap& aggregates_)
664 : subdomains(subdomains_), max(-1), subdomain(-1), aggregates(aggregates_)
667 void operator()(
const T& edge)
670 subdomains[subdomain].insert(edge.target());
672 int setAggregate(
const AggregateDescriptor& aggregate_)
674 subdomain=aggregate_;
675 max = std::max(subdomain, aggregate_);
678 int noSubdomains()
const
684 AggregateDescriptor max;
685 AggregateDescriptor subdomain;
686 const AggregatesMap& aggregates;
691 void operator()(
const T& edge)
693 int setAggregate(
const AggregateDescriptor& aggregate_)
697 int noSubdomains()
const
704 struct AggregateAdder
706 AggregateAdder(Vector& subdomains_,
const AggregatesMap& aggregates_,
707 const MatrixGraph<const M>& graph_, VM& visitedMap_)
708 : subdomains(subdomains_), subdomain(-1), aggregates(aggregates_),
709 adder(subdomains_, aggregates_), graph(graph_), visitedMap(visitedMap_)
712 void operator()(
const T& edge)
714 subdomains[subdomain].insert(edge.target());
719 assert(aggregates[edge.target()]!=aggregate);
721 aggregates.template breadthFirstSearch<true,false>(edge.target(), aggregate,
722 graph, vlist, adder, adder,
727 int setAggregate(
const AggregateDescriptor& aggregate_)
729 adder.setAggregate(aggregate_);
730 aggregate=aggregate_;
733 int noSubdomains()
const
739 AggregateDescriptor aggregate;
742 const AggregatesMap& aggregates;
744 const MatrixGraph<const M>& graph;
748 void createPairDomains(
const MatrixGraph<const M>& graph)
752 typedef typename M::size_type size_type;
754 std::set<std::pair<size_type,size_type> > pairs;
756 for(VIter v=graph.begin(), ve=graph.end(); ve != v; ++v)
757 for(EIter e = v.begin(), ee=v.end(); ee!=e; ++e)
760 if(e.source()<e.target())
761 pairs.insert(std::make_pair(e.source(),e.target()));
763 pairs.insert(std::make_pair(e.target(),e.source()));
767 subdomains.resize(pairs.size());
768 Dune::dinfo <<std::endl<<
"Created "<<pairs.size()<<
" ("<<total<<
") pair domains"<<std::endl<<std::endl;
769 typedef typename std::set<std::pair<size_type,size_type> >::const_iterator SIter;
770 typename Vector::iterator subdomain=subdomains.begin();
772 for(SIter s=pairs.begin(), se =pairs.end(); se!=s; ++s)
774 subdomain->insert(s->first);
775 subdomain->insert(s->second);
778 std::size_t minsize=10000;
779 std::size_t maxsize=0;
781 for(
typename Vector::size_type i=0; i < subdomains.size(); ++i) {
782 sum+=subdomains[i].size();
783 minsize=std::min(minsize, subdomains[i].size());
784 maxsize=std::max(maxsize, subdomains[i].size());
786 Dune::dinfo<<
"Subdomain size: min="<<minsize<<
" max="<<maxsize<<
" avg="<<(sum/subdomains.size())
787 <<
" no="<<subdomains.size()<<std::endl;
790 template<
class Visitor>
791 void createSubdomains(
const M& matrix,
const MatrixGraph<const M>& graph,
792 const AggregatesMap& amap, Visitor& overlapVisitor,
793 IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap>& visitedMap )
800 AggregateDescriptor maxAggregate=0;
802 for(std::size_t i=0; i < amap.noVertices(); ++i)
806 maxAggregate = std::max(maxAggregate, amap[i]);
808 subdomains.resize(maxAggregate+1+isolated);
811 for(
typename Vector::size_type i=0; i < subdomains.size(); ++i)
812 subdomains[i].clear();
817 VertexAdder aggregateVisitor(subdomains, amap);
819 for(VertexDescriptor i=0; i < amap.noVertices(); ++i)
820 if(!get(visitedMap, i)) {
821 AggregateDescriptor aggregate=amap[i];
825 subdomains.push_back(Subdomain());
826 aggregate=subdomains.size()-1;
828 overlapVisitor.setAggregate(aggregate);
829 aggregateVisitor.setAggregate(aggregate);
830 subdomains[aggregate].insert(i);
832 amap.template breadthFirstSearch<false,false>(i, aggregate, graph, vlist, aggregateVisitor,
833 overlapVisitor, visitedMap);
836 std::size_t minsize=10000;
837 std::size_t maxsize=0;
839 for(
typename Vector::size_type i=0; i < subdomains.size(); ++i) {
840 sum+=subdomains[i].size();
841 minsize=std::min(minsize, subdomains[i].size());
842 maxsize=std::max(maxsize, subdomains[i].size());
844 Dune::dinfo<<
"Subdomain size: min="<<minsize<<
" max="<<maxsize<<
" avg="<<(sum/subdomains.size())
845 <<
" no="<<subdomains.size()<<
" isolated="<<isolated<<std::endl;
854 template<
class M,
class X,
class TM,
class TS,
class TA>
855 struct ConstructionTraits<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
857 typedef ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TS,TA> >
Arguments;
861 return new SeqOverlappingSchwarz<M,X,TM,TS,TA>(args.getMatrix(),
862 args.getSubDomains(),
863 args.getArgs().relaxationFactor,
864 args.getArgs().onthefly);
867 static void deconstruct(SeqOverlappingSchwarz<M,X,TM,TS,TA>* schwarz)
Provides classes for the Coloring process of AMG.
Class providing information about the mapping of the vertices onto aggregates.
Definition: aggregates.hh:498
Traits class for generically constructing non default constructable types.
Definition: construction.hh:39
Construction Arguments for the default smoothers.
Definition: smoother.hh:88
VertexIteratorT< const MatrixGraph< Matrix > > ConstVertexIterator
The constant vertex iterator type.
Definition: graph.hh:307
M::size_type VertexDescriptor
The vertex descriptor.
Definition: graph.hh:72
EdgeIteratorT< const MatrixGraph< Matrix > > ConstEdgeIterator
The constant edge iterator type.
Definition: graph.hh:297
Block parallel preconditioner.
Definition: schwarz.hh:288
Nonoverlapping parallel preconditioner.
Definition: novlpschwarz.hh:340
Sequential ILU0 preconditioner.
Definition: preconditioners.hh:488
Sequential ILU(n) preconditioner.
Definition: preconditioners.hh:572
The sequential jacobian preconditioner.
Definition: preconditioners.hh:402
std::vector< subdomain_type, typename TA::template rebind< subdomain_type >::other > subdomain_vector
The vector type containing the subdomain to row index mapping.
Definition: overlappingschwarz.hh:795
Sequential SOR preconditioner.
Definition: preconditioners.hh:215
Sequential SSOR preconditioner.
Definition: preconditioners.hh:127
Helper classes for the construction of classes without empty constructor.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:244
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:408
V AggregateDescriptor
The aggregate descriptor type.
Definition: aggregates.hh:518
static const V ISOLATED
Identifier of isolated vertices.
Definition: aggregates.hh:509
DefaultSmootherArgs()
Default constructor.
Definition: smoother.hh:54
const void * Arguments
A type holding all the arguments needed to call the constructor.
Definition: construction.hh:45
static void postSmooth(Smoother &smoother, Domain &v, const Range &d)
apply post smoothing in forward direction
Definition: smoother.hh:396
SLList< VertexDescriptor, Allocator > VertexList
The type of a single linked list of vertex descriptors.
Definition: aggregates.hh:530
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:430
static T * construct(Arguments &args)
Construct an object with the specified arguments.
Definition: construction.hh:53
RelaxationFactor relaxationFactor
The relaxation factor to use.
Definition: smoother.hh:49
static void preSmooth(Smoother &smoother, Domain &v, const Range &d)
apply pre smoothing in forward direction
Definition: smoother.hh:384
static void deconstruct(T *t)
Destroys an object.
Definition: construction.hh:62
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:139
Dune namespace.
Definition: alignment.hh:14
Define general preconditioner interface.
The default class for the smoother arguments.
Definition: smoother.hh:36
Helper class for applying the smoothers.
Definition: smoother.hh:372
Traits class for getting the attribute class of a smoother.
Definition: smoother.hh:64
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18