smoother.hh

Go to the documentation of this file.
00001 #ifndef DUNE_AMGSMOOTHER_HH
00002 #define DUNE_AMGSMOOTHER_HH
00003 
00004 #include<dune/istl/paamg/construction.hh>
00005 #include<dune/istl/paamg/aggregates.hh>
00006 #include<dune/istl/preconditioners.hh>
00007 #include<dune/istl/schwarz.hh>
00008 #include<dune/common/propertymap.hh>
00009 
00010 namespace Dune
00011 {
00012   namespace Amg
00013   {
00014     
00030     template<class T>
00031     struct DefaultSmootherArgs
00032     {
00036       typedef T RelaxationFactor;
00037       
00041       int iterations;
00045       RelaxationFactor relaxationFactor;
00046       
00050       DefaultSmootherArgs()
00051         : iterations(1), relaxationFactor(1.0)
00052       {}
00053     };
00054     
00058     template<class T>
00059     struct SmootherTraits
00060     {
00061       typedef DefaultSmootherArgs<typename T::matrix_type::field_type> Arguments;
00062       
00063     };
00064     
00065     template<class X, class Y, class C, class T>
00066     struct SmootherTraits<BlockPreconditioner<X,Y,C,T> >
00067     {
00068       typedef DefaultSmootherArgs<typename T::matrix_type::field_type> Arguments;
00069       
00070     };
00071 
00072     
00076     template<class T>
00077     class DefaultConstructionArgs
00078     {      
00079       typedef typename T::matrix_type Matrix;
00080       
00081       typedef typename SmootherTraits<T>::Arguments SmootherArgs;
00082       
00083       typedef AggregatesMap<typename MatrixGraph<Matrix>::VertexDescriptor> AggregatesMap;
00084       
00085     public:
00086       virtual ~DefaultConstructionArgs()
00087       {}
00088       
00089       void setMatrix(const Matrix& matrix)
00090       {
00091         matrix_=&matrix;
00092       }
00093       virtual void setMatrix(const Matrix& matrix, const AggregatesMap& amap)
00094       {
00095         setMatrix(matrix);
00096       }
00097       
00098      
00099       const Matrix& getMatrix() const
00100       {
00101         return *matrix_;
00102       }
00103       
00104       void setArgs(const SmootherArgs& args)
00105       {
00106         args_=&args;
00107       }
00108 
00109       template<class T1>
00110       void setComm(T1& comm)
00111       {}
00112       
00113       const SmootherArgs getArgs() const
00114       {
00115         return *args_;
00116       }
00117       
00118     protected:
00119       const Matrix* matrix_;
00120     private:
00121       const SmootherArgs* args_;
00122     };
00123 
00124     template<class T>
00125     struct ConstructionArgs
00126       : public DefaultConstructionArgs<T>
00127     {};
00128 
00129     template<class T, class C=SequentialInformation>
00130     class DefaultParallelConstructionArgs 
00131       : public ConstructionArgs<T>
00132     {
00133     public:
00134       virtual ~DefaultParallelConstructionArgs()
00135       {}
00136       
00137       void setComm(const C& comm)
00138       {
00139         comm_ = &comm;
00140       }
00141 
00142       const C& getComm() const
00143       {
00144         return *comm_;
00145       }
00146     private:
00147       const C* comm_;
00148     };
00149     
00150 
00151     template<class T>
00152     class ConstructionTraits;
00153 
00157     template<class M, class X, class Y, int l>
00158     struct ConstructionTraits<SeqSSOR<M,X,Y,l> >
00159     {
00160       typedef DefaultConstructionArgs<SeqSSOR<M,X,Y,l> > Arguments;
00161       
00162       static inline SeqSSOR<M,X,Y,l>* construct(Arguments& args)
00163       {
00164         return new SeqSSOR<M,X,Y,l>(args.getMatrix(), args.getArgs().iterations,
00165                                     args.getArgs().relaxationFactor);
00166       }
00167       
00168       static inline void deconstruct(SeqSSOR<M,X,Y,l>* ssor)
00169       {
00170         delete ssor;
00171       }
00172       
00173     };
00174 
00175     
00179     template<class M, class X, class Y, int l>
00180     struct ConstructionTraits<SeqSOR<M,X,Y,l> >
00181     {
00182       typedef DefaultConstructionArgs<SeqSOR<M,X,Y,l> > Arguments;
00183       
00184       static inline SeqSOR<M,X,Y,l>* construct(Arguments& args)
00185       {
00186         return new SeqSOR<M,X,Y,l>(args.getMatrix(), args.getArgs().iterations,
00187                                     args.getArgs().relaxationFactor);
00188       }
00189       
00190       static inline void deconstruct(SeqSOR<M,X,Y,l>* sor)
00191       {
00192         delete sor;
00193       }
00194       
00195     };
00199     template<class M, class X, class Y, int l>
00200     struct ConstructionTraits<SeqJac<M,X,Y,l> >
00201     {
00202       typedef DefaultConstructionArgs<SeqJac<M,X,Y,l> > Arguments;
00203       
00204       static inline SeqJac<M,X,Y,l>* construct(Arguments& args)
00205       {
00206         return new SeqJac<M,X,Y,l>(args.getMatrix(), args.getArgs().iterations,
00207                                  args.getArgs().relaxationFactor);
00208       }
00209       
00210       static void deconstruct(SeqJac<M,X,Y,l>* jac)
00211       {
00212         delete jac;
00213       }
00214       
00215     };
00216 
00217     
00221     template<class M, class X, class Y>
00222     struct ConstructionTraits<SeqILU0<M,X,Y> >
00223     {
00224       typedef DefaultConstructionArgs<SeqILU0<M,X,Y> > Arguments;
00225       
00226       static inline SeqILU0<M,X,Y>* construct(Arguments& args)
00227       {
00228         return new SeqILU0<M,X,Y>(args.getMatrix(),
00229                                  args.getArgs().relaxationFactor);
00230       }
00231       
00232       static void deconstruct(SeqILU0<M,X,Y>* ilu)
00233       {
00234         delete ilu;
00235       }
00236       
00237     };
00238 
00239     template<class M, class X, class Y>
00240     class ConstructionArgs<SeqILUn<M,X,Y> >
00241       : public DefaultConstructionArgs<SeqILUn<M,X,Y> >
00242     {
00243     public:
00244       ConstructionArgs(int n=1)
00245         : n_(n)
00246       {}
00247       
00248       void setN(int n)
00249       {
00250         n_ = n;
00251       }
00252       int getN()
00253       {
00254         return n_;
00255       }
00256       
00257     private:
00258       int n_;
00259     };
00260     
00261     
00265     template<class M, class X, class Y>
00266     struct ConstructionTraits<SeqILUn<M,X,Y> >
00267     {
00268       typedef ConstructionArgs<SeqILUn<M,X,Y> > Arguments;
00269       
00270       static inline SeqILUn<M,X,Y>* construct(Arguments& args)
00271       {
00272         return new SeqILUn<M,X,Y>(args.getMatrix(), args.getN(),
00273                                  args.getArgs().relaxationFactor);
00274       }
00275       
00276       static void deconstruct(SeqILUn<M,X,Y>* ilu)
00277       {
00278         delete ilu;
00279       }
00280       
00281     };
00282 
00283 
00287     template<class M, class X, class Y, class C>
00288     struct ConstructionTraits<ParSSOR<M,X,Y,C> >
00289     {
00290       typedef DefaultParallelConstructionArgs<M,C> Arguments;
00291       
00292       static inline ParSSOR<M,X,Y,C>* construct(Arguments& args)
00293       {
00294         return new ParSSOR<M,X,Y,C>(args.getMatrix(), args.getArgs().iterations,
00295                                     args.getArgs().relaxationFactor,
00296                                     args.getComm());
00297       }      
00298       static inline void deconstruct(ParSSOR<M,X,Y,C>* ssor)
00299       {
00300         delete ssor;
00301       }
00302     };
00303 
00304     template<class X, class Y, class C, class T>
00305     struct ConstructionTraits<BlockPreconditioner<X,Y,C,T> >
00306     {
00307       typedef DefaultParallelConstructionArgs<T,C> Arguments;
00308       typedef ConstructionTraits<T> SeqConstructionTraits;
00309       static inline BlockPreconditioner<X,Y,C,T>* construct(Arguments& args)
00310       {
00311         return new BlockPreconditioner<X,Y,C,T>(*SeqConstructionTraits::construct(args),
00312                                                 args.getComm());
00313       }
00314 
00315       static inline void deconstruct(BlockPreconditioner<X,Y,C,T>* bp)
00316       {
00317         SeqConstructionTraits::deconstruct(static_cast<T*>(&bp->preconditioner));
00318         delete bp;
00319       }
00320       
00321     };
00322 
00333     template<class T>
00334     struct SmootherApplier
00335     {
00336       typedef T Smoother;
00337       typedef typename Smoother::range_type Range;
00338       typedef typename Smoother::domain_type Domain;
00339       
00347       static void preSmooth(Smoother& smoother, Domain& v, const Range& d)
00348       {
00349         smoother.apply(v,d);
00350       }
00351 
00359       static void postSmooth(Smoother& smoother, Domain& v, const Range& d)
00360       {
00361         smoother.apply(v,d);
00362       }
00363     };
00364     
00365     template<class M, class X, class Y, int l>
00366     struct SmootherApplier<SeqSOR<M,X,Y,l> >
00367     {
00368       typedef SeqSOR<M,X,Y,l> Smoother;
00369       typedef typename Smoother::range_type Range;
00370       typedef typename Smoother::domain_type Domain;
00371       
00372       static void preSmooth(Smoother& smoother, Domain& v, Range& d)
00373       {
00374         smoother.template apply<true>(v,d);
00375       }
00376 
00377        
00378       static void postSmooth(Smoother& smoother, Domain& v, Range& d)
00379       {
00380         smoother.template apply<false>(v,d);
00381       }
00382     };
00383 #ifdef HAVE_SUPERLU
00384     
00385     // forward dseclarations
00386     template<class M, class X, class MO, class A>
00387     class SeqOverlappingSchwarz;
00388 
00389     class MultiplicativeSchwarzMode;
00390     
00391     template<class M, class X, class TA>
00392     struct SmootherApplier<SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,TA> >
00393     {
00394       typedef SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,TA> Smoother;
00395       typedef typename Smoother::range_type Range;
00396       typedef typename Smoother::domain_type Domain;
00397       
00398       static void preSmooth(Smoother& smoother, Domain& v, const Range& d)
00399       {
00400         smoother.template apply<true>(v,d);
00401       }
00402 
00403        
00404       static void postSmooth(Smoother& smoother, Domain& v, const Range& d)
00405       {
00406         smoother.template apply<false>(v,d);
00407       
00408        }
00409     };
00410     
00411     //    template<class M, class X, class TM, class TA>
00412     //    class SeqOverlappingSchwarz;
00413 
00414     template<class T>
00415     struct SeqOverlappingSchwarzSmootherArgs
00416       : public DefaultSmootherArgs<T>
00417     {
00418       enum Overlap {vertex, aggregate, none};
00419       
00420       Overlap overlap;
00421       SeqOverlappingSchwarzSmootherArgs()
00422         : overlap(none)
00423       {}
00424     };  
00425       
00426     template<class M, class X, class TM, class TA>
00427     struct SmootherTraits<SeqOverlappingSchwarz<M,X,TM,TA> >
00428     {
00429       typedef  SeqOverlappingSchwarzSmootherArgs<typename M::field_type> Arguments;
00430     };
00431     
00432     template<class M, class X, class TM, class TA>
00433     class ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TA> >
00434       : public DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TA> >
00435     {
00436       typedef DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TA> > Father;
00437 
00438     public:
00439       typedef AggregatesMap<typename MatrixGraph<M>::VertexDescriptor> AggregatesMap;
00440       typedef typename AggregatesMap::AggregateDescriptor AggregateDescriptor;
00441       typedef typename AggregatesMap::VertexDescriptor VertexDescriptor;
00442       typedef typename SeqOverlappingSchwarz<M,X,TM,TA>::subdomain_vector Vector;
00443       typedef typename Vector::value_type Subdomain;
00444       
00445       virtual void setMatrix(const M& matrix, const AggregatesMap& amap)
00446       {
00447         Father::setMatrix(matrix);
00448         
00449         std::vector<bool> visited(amap.noVertices(), false);
00450         typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
00451         VisitedMapType visitedMap(visited.begin());
00452         
00453         MatrixGraph<const M> graph(matrix);
00454         
00455         typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
00456         
00457         switch(Father::getArgs().overlap){
00458         case SmootherArgs::vertex:
00459           {  
00460           VertexAdder visitor(subdomains, amap);
00461           createSubdomains(matrix, graph, amap, visitor,  visitedMap);
00462           }
00463           break;
00464         case SmootherArgs::aggregate:
00465           {
00466           AggregateAdder<VisitedMapType> visitor(subdomains, amap, graph, visitedMap);
00467           createSubdomains(matrix, graph, amap, visitor, visitedMap);
00468           }
00469           break;
00470         case SmootherArgs::none:
00471           NoneAdder visitor;
00472           createSubdomains(matrix, graph, amap, visitor, visitedMap);     
00473         }
00474       }
00475       void setMatrix(const M& matrix)
00476       {
00477         Father::setMatrix(matrix);
00478       }
00479       
00480       const Vector& getSubDomains()
00481       {
00482         return subdomains;
00483       }
00484       
00485     private:
00486       struct VertexAdder
00487       {
00488         VertexAdder(Vector& subdomains_, const AggregatesMap& aggregates_)
00489           : subdomains(subdomains_), max(-1), subdomain(-1), aggregates(aggregates_)
00490         {}
00491         template<class T>
00492         void operator()(const T& edge)
00493         {
00494           if(aggregates[edge.target()]!=AggregatesMap::ISOLATED)
00495             subdomains[subdomain].insert(edge.target());
00496         }
00497         int setAggregate(const AggregateDescriptor& aggregate_)
00498         {
00499           subdomain=aggregate_;
00500           max = std::max(subdomain, aggregate_);
00501           return subdomain;
00502         }
00503         int noSubdomains() const
00504         {
00505           return max+1;
00506         }
00507       private:
00508         Vector& subdomains;
00509         int max;
00510         int subdomain;
00511         const AggregatesMap& aggregates;
00512       };
00513       struct NoneAdder
00514       {
00515         template<class T>
00516         void operator()(const T& edge)
00517         {}
00518         int setAggregate(const AggregateDescriptor& aggregate_)
00519         {
00520         return -1;
00521         }
00522         int noSubdomains() const
00523         {
00524           return -1;
00525         }
00526       };
00527 
00528       template<class VM>
00529       struct AggregateAdder
00530       {
00531         AggregateAdder(Vector& subdomains_, const AggregatesMap& aggregates_, 
00532                        const MatrixGraph<const M>& graph_, VM& visitedMap_)
00533           : subdomains(subdomains_), subdomain(-1), aggregates(aggregates_),
00534             adder(subdomains_, aggregates_), graph(graph_), visitedMap(visitedMap_)
00535         {}
00536         template<class T>
00537         void operator()(const T& edge)
00538         {
00539           subdomains[subdomain].insert(edge.target());
00540           // If we (the neighbouring vertex of the aggregate)
00541           // are not isolated, add the aggregate we belong to 
00542           // to the same subdomain using the OneOverlapAdder
00543           if(aggregates[edge.target()]!=AggregatesMap::ISOLATED){
00544             assert(aggregates[edge.target()]!=aggregate);
00545             typename AggregatesMap::VertexList vlist;       
00546             aggregates.template breadthFirstSearch<true,false>(edge.target(), aggregate,
00547                                                                graph, vlist, adder, adder,
00548                                                          visitedMap);
00549           }
00550         }
00551         
00552         int setAggregate(const AggregateDescriptor& aggregate_)
00553         {
00554           adder.setAggregate(aggregate_);
00555           aggregate=aggregate_;
00556           return ++subdomain;
00557         }
00558         int noSubdomains() const
00559         {
00560           return subdomain+1;
00561         }
00562         
00563       private:
00564         AggregateDescriptor aggregate;
00565         Vector& subdomains;
00566         int subdomain;
00567         const AggregatesMap& aggregates;
00568         VertexAdder adder;
00569         const MatrixGraph<const M>& graph;
00570         VM& visitedMap;
00571       };
00572       
00573       template<class Visitor>
00574       void createSubdomains(const M& matrix, const MatrixGraph<const M>& graph, 
00575                             const AggregatesMap& amap, Visitor& overlapVisitor, 
00576                             IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap>& visitedMap )
00577       {
00578         // count  number ag aggregates. We asume that the
00579         // aggregates are numbered consecutively from 0 exept
00580         // for the isolated ones. All isolated vertices form 
00581         // one aggregate, here.
00582         int isolated=false;
00583         AggregateDescriptor maxAggregate=0;
00584         
00585         for(int i=0; i < amap.noVertices(); ++i)
00586           if(amap[i]==AggregatesMap::ISOLATED)
00587             isolated=true;
00588           else
00589             maxAggregate = std::max(maxAggregate, amap[i]);
00590         
00591         subdomains.resize(maxAggregate+1);
00592 
00593         // reset the subdomains
00594         for(int i=0; i < subdomains.size(); ++i)
00595           subdomains[i].clear();
00596         
00597         // Create the subdomains from the aggregates mapping.
00598         // For each aggregate we mark all entries and the 
00599         // neighbouring vertices as belonging to the same subdomain
00600         VertexAdder aggregateVisitor(subdomains, amap);
00601                 
00602         for(VertexDescriptor i=0; i < amap.noVertices(); ++i)
00603           if(!get(visitedMap, i)){
00604             AggregateDescriptor aggregate=amap[i];
00605 
00606             if(amap[i]==AggregatesMap::ISOLATED){
00607               // isolated vertex gets its own aggregate
00608               subdomains.push_back(Subdomain());
00609               aggregate=subdomains.size()-1;
00610             }
00611             overlapVisitor.setAggregate(aggregate);
00612             aggregateVisitor.setAggregate(aggregate);
00613             subdomains[aggregate].insert(i);
00614             typename AggregatesMap::VertexList vlist;       
00615             amap.template breadthFirstSearch<false,false>(i, aggregate, graph, vlist, aggregateVisitor, 
00616                             overlapVisitor, visitedMap);
00617           }
00618         subdomains.resize(aggregateVisitor.noSubdomains());
00619         
00620         std::size_t minsize=10000;
00621         std::size_t maxsize=0;
00622         int sum=0;
00623         for(int i=0; i < subdomains.size(); ++i){
00624           sum+=subdomains[i].size();
00625           minsize=std::min(minsize, subdomains[i].size());
00626           maxsize=std::max(maxsize, subdomains[i].size());
00627         }
00628         std::cout<<"Subdomain size: min="<<minsize<<" max="<<maxsize<<" avg="<<(sum/subdomains.size())
00629                  <<" no="<<subdomains.size()<<std::endl;
00630         
00631           
00632           
00633       }
00634       Vector subdomains;
00635     };
00636     
00637 
00638     template<class M, class X, class TM, class TA>
00639     struct ConstructionTraits<SeqOverlappingSchwarz<M,X,TM,TA> >
00640     {
00641       typedef ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TA> > Arguments;
00642       
00643       static inline SeqOverlappingSchwarz<M,X,TM,TA>* construct(Arguments& args)
00644       {
00645         return new SeqOverlappingSchwarz<M,X,TM,TA>(args.getMatrix(), 
00646                                                   args.getSubDomains(),
00647                                                   args.getArgs().relaxationFactor);
00648       }
00649       
00650       static void deconstruct(SeqOverlappingSchwarz<M,X,TM,TA>* schwarz)
00651       {
00652         delete schwarz;
00653       }
00654     };
00655 #endif
00656   } // namespace Amg
00657 } // namespace Dune
00658 
00659 
00660 
00661 #endif

Generated on 6 Nov 2008 with Doxygen (ver 1.5.6) [logfile].