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 Dune::Amg::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
00286 template<class M, class X, class Y, class C>
00287 struct ConstructionTraits<ParSSOR<M,X,Y,C> >
00288 {
00289 typedef DefaultParallelConstructionArgs<M,C> Arguments;
00290
00291 static inline ParSSOR<M,X,Y,C>* construct(Arguments& args)
00292 {
00293 return new ParSSOR<M,X,Y,C>(args.getMatrix(), args.getArgs().iterations,
00294 args.getArgs().relaxationFactor,
00295 args.getComm());
00296 }
00297 static inline void deconstruct(ParSSOR<M,X,Y,C>* ssor)
00298 {
00299 delete ssor;
00300 }
00301 };
00302
00303 template<class X, class Y, class C, class T>
00304 struct ConstructionTraits<BlockPreconditioner<X,Y,C,T> >
00305 {
00306 typedef DefaultParallelConstructionArgs<T,C> Arguments;
00307 typedef ConstructionTraits<T> SeqConstructionTraits;
00308 static inline BlockPreconditioner<X,Y,C,T>* construct(Arguments& args)
00309 {
00310 return new BlockPreconditioner<X,Y,C,T>(*SeqConstructionTraits::construct(args),
00311 args.getComm());
00312 }
00313
00314 static inline void deconstruct(BlockPreconditioner<X,Y,C,T>* bp)
00315 {
00316 SeqConstructionTraits::deconstruct(static_cast<T*>(&bp->preconditioner));
00317 delete bp;
00318 }
00319
00320 };
00321
00332 template<class T>
00333 struct SmootherApplier
00334 {
00335 typedef T Smoother;
00336 typedef typename Smoother::range_type Range;
00337 typedef typename Smoother::domain_type Domain;
00338
00346 static void preSmooth(Smoother& smoother, Domain& v, const Range& d)
00347 {
00348 smoother.apply(v,d);
00349 }
00350
00358 static void postSmooth(Smoother& smoother, Domain& v, const Range& d)
00359 {
00360 smoother.apply(v,d);
00361 }
00362 };
00363
00364 template<class M, class X, class Y, int l>
00365 struct SmootherApplier<SeqSOR<M,X,Y,l> >
00366 {
00367 typedef SeqSOR<M,X,Y,l> Smoother;
00368 typedef typename Smoother::range_type Range;
00369 typedef typename Smoother::domain_type Domain;
00370
00371 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
00372 {
00373 smoother.template apply<true>(v,d);
00374 }
00375
00376
00377 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
00378 {
00379 smoother.template apply<false>(v,d);
00380 }
00381 };
00382 #ifdef HAVE_SUPERLU
00383 }
00384
00385
00386 template<class M, class X, class MO, class A>
00387 class SeqOverlappingSchwarz;
00388
00389 class MultiplicativeSchwarzMode;
00390
00391 namespace Amg
00392 {
00393 template<class M, class X, class TA>
00394 struct SmootherApplier<SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,TA> >
00395 {
00396 typedef SeqOverlappingSchwarz<M,X,MultiplicativeSchwarzMode,TA> Smoother;
00397 typedef typename Smoother::range_type Range;
00398 typedef typename Smoother::domain_type Domain;
00399
00400 static void preSmooth(Smoother& smoother, Domain& v, const Range& d)
00401 {
00402 smoother.template apply<true>(v,d);
00403 }
00404
00405
00406 static void postSmooth(Smoother& smoother, Domain& v, const Range& d)
00407 {
00408 smoother.template apply<false>(v,d);
00409
00410 }
00411 };
00412
00413
00414
00415
00416 template<class T>
00417 struct SeqOverlappingSchwarzSmootherArgs
00418 : public DefaultSmootherArgs<T>
00419 {
00420 enum Overlap {vertex, aggregate, none};
00421
00422 Overlap overlap;
00423 SeqOverlappingSchwarzSmootherArgs(Overlap overlap_=vertex)
00424 : overlap(overlap)
00425 {}
00426 };
00427
00428 template<class M, class X, class TM, class TA>
00429 struct SmootherTraits<SeqOverlappingSchwarz<M,X,TM,TA> >
00430 {
00431 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> Arguments;
00432 };
00433
00434 template<class M, class X, class TM, class TA>
00435 class ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TA> >
00436 : public DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TA> >
00437 {
00438 typedef DefaultConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TA> > Father;
00439
00440 public:
00441 typedef typename MatrixGraph<M>::VertexDescriptor VertexDescriptor;
00442 typedef Dune::Amg::AggregatesMap<VertexDescriptor> AggregatesMap;
00443 typedef typename AggregatesMap::AggregateDescriptor AggregateDescriptor;
00444 typedef typename SeqOverlappingSchwarz<M,X,TM,TA>::subdomain_vector Vector;
00445 typedef typename Vector::value_type Subdomain;
00446
00447 virtual void setMatrix(const M& matrix, const AggregatesMap& amap)
00448 {
00449 Father::setMatrix(matrix);
00450
00451 std::vector<bool> visited(amap.noVertices(), false);
00452 typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
00453 VisitedMapType visitedMap(visited.begin());
00454
00455 MatrixGraph<const M> graph(matrix);
00456
00457 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
00458
00459 switch(Father::getArgs().overlap){
00460 case SmootherArgs::vertex:
00461 {
00462 VertexAdder visitor(subdomains, amap);
00463 createSubdomains(matrix, graph, amap, visitor, visitedMap);
00464 }
00465 break;
00466 case SmootherArgs::aggregate:
00467 {
00468 AggregateAdder<VisitedMapType> visitor(subdomains, amap, graph, visitedMap);
00469 createSubdomains(matrix, graph, amap, visitor, visitedMap);
00470 }
00471 break;
00472 case SmootherArgs::none:
00473 NoneAdder visitor;
00474 createSubdomains(matrix, graph, amap, visitor, visitedMap);
00475 }
00476 }
00477 void setMatrix(const M& matrix)
00478 {
00479 Father::setMatrix(matrix);
00480
00481
00482 AggregatesMap amap(matrix.N());
00483 VertexDescriptor v=0;
00484 for(typename AggregatesMap::iterator iter=amap.begin();
00485 iter!=amap.end(); ++iter)
00486 *iter=v++;
00487
00488 std::vector<bool> visited(amap.noVertices(), false);
00489 typedef IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap> VisitedMapType;
00490 VisitedMapType visitedMap(visited.begin());
00491
00492 MatrixGraph<const M> graph(matrix);
00493
00494 typedef SeqOverlappingSchwarzSmootherArgs<typename M::field_type> SmootherArgs;
00495
00496 switch(Father::getArgs().overlap){
00497 case SmootherArgs::vertex:
00498 {
00499 VertexAdder visitor(subdomains, amap);
00500 createSubdomains(matrix, graph, amap, visitor, visitedMap);
00501 }
00502 break;
00503 case SmootherArgs::aggregate:
00504 {
00505 DUNE_THROW(NotImplemented, "Aggregate overlap is not supported yet");
00506
00507
00508
00509
00510 }
00511 break;
00512 case SmootherArgs::none:
00513 NoneAdder visitor;
00514 createSubdomains(matrix, graph, amap, visitor, visitedMap);
00515
00516 }
00517 }
00518
00519 const Vector& getSubDomains()
00520 {
00521 return subdomains;
00522 }
00523
00524 private:
00525 struct VertexAdder
00526 {
00527 VertexAdder(Vector& subdomains_, const AggregatesMap& aggregates_)
00528 : subdomains(subdomains_), max(-1), subdomain(-1), aggregates(aggregates_)
00529 {}
00530 template<class T>
00531 void operator()(const T& edge)
00532 {
00533 if(aggregates[edge.target()]!=AggregatesMap::ISOLATED)
00534 subdomains[subdomain].insert(edge.target());
00535 }
00536 int setAggregate(const AggregateDescriptor& aggregate_)
00537 {
00538 subdomain=aggregate_;
00539 max = std::max(subdomain, aggregate_);
00540 return subdomain;
00541 }
00542 int noSubdomains() const
00543 {
00544 return max+1;
00545 }
00546 private:
00547 Vector& subdomains;
00548 AggregateDescriptor max;
00549 AggregateDescriptor subdomain;
00550 const AggregatesMap& aggregates;
00551 };
00552 struct NoneAdder
00553 {
00554 template<class T>
00555 void operator()(const T& edge)
00556 {}
00557 int setAggregate(const AggregateDescriptor& aggregate_)
00558 {
00559 return -1;
00560 }
00561 int noSubdomains() const
00562 {
00563 return -1;
00564 }
00565 };
00566
00567 template<class VM>
00568 struct AggregateAdder
00569 {
00570 AggregateAdder(Vector& subdomains_, const AggregatesMap& aggregates_,
00571 const MatrixGraph<const M>& graph_, VM& visitedMap_)
00572 : subdomains(subdomains_), subdomain(-1), aggregates(aggregates_),
00573 adder(subdomains_, aggregates_), graph(graph_), visitedMap(visitedMap_)
00574 {}
00575 template<class T>
00576 void operator()(const T& edge)
00577 {
00578 subdomains[subdomain].insert(edge.target());
00579
00580
00581
00582 if(aggregates[edge.target()]!=AggregatesMap::ISOLATED){
00583 assert(aggregates[edge.target()]!=aggregate);
00584 typename AggregatesMap::VertexList vlist;
00585 aggregates.template breadthFirstSearch<true,false>(edge.target(), aggregate,
00586 graph, vlist, adder, adder,
00587 visitedMap);
00588 }
00589 }
00590
00591 int setAggregate(const AggregateDescriptor& aggregate_)
00592 {
00593 adder.setAggregate(aggregate_);
00594 aggregate=aggregate_;
00595 return ++subdomain;
00596 }
00597 int noSubdomains() const
00598 {
00599 return subdomain+1;
00600 }
00601
00602 private:
00603 AggregateDescriptor aggregate;
00604 Vector& subdomains;
00605 int subdomain;
00606 const AggregatesMap& aggregates;
00607 VertexAdder adder;
00608 const MatrixGraph<const M>& graph;
00609 VM& visitedMap;
00610 };
00611
00612 template<class Visitor>
00613 void createSubdomains(const M& matrix, const MatrixGraph<const M>& graph,
00614 const AggregatesMap& amap, Visitor& overlapVisitor,
00615 IteratorPropertyMap<std::vector<bool>::iterator,IdentityMap>& visitedMap )
00616 {
00617
00618
00619
00620
00621 int isolated=false;
00622 AggregateDescriptor maxAggregate=0;
00623
00624 for(int i=0; i < amap.noVertices(); ++i)
00625 if(amap[i]==AggregatesMap::ISOLATED)
00626 isolated=true;
00627 else
00628 maxAggregate = std::max(maxAggregate, amap[i]);
00629
00630 subdomains.resize(maxAggregate+1);
00631
00632
00633 for(int i=0; i < subdomains.size(); ++i)
00634 subdomains[i].clear();
00635
00636
00637
00638
00639 VertexAdder aggregateVisitor(subdomains, amap);
00640
00641 for(VertexDescriptor i=0; i < amap.noVertices(); ++i)
00642 if(!get(visitedMap, i)){
00643 AggregateDescriptor aggregate=amap[i];
00644
00645 if(amap[i]==AggregatesMap::ISOLATED){
00646
00647 subdomains.push_back(Subdomain());
00648 aggregate=subdomains.size()-1;
00649 }
00650 overlapVisitor.setAggregate(aggregate);
00651 aggregateVisitor.setAggregate(aggregate);
00652 subdomains[aggregate].insert(i);
00653 typename AggregatesMap::VertexList vlist;
00654 amap.template breadthFirstSearch<false,false>(i, aggregate, graph, vlist, aggregateVisitor,
00655 overlapVisitor, visitedMap);
00656 }
00657 subdomains.resize(aggregateVisitor.noSubdomains());
00658
00659 std::size_t minsize=10000;
00660 std::size_t maxsize=0;
00661 int sum=0;
00662 for(int i=0; i < subdomains.size(); ++i){
00663 sum+=subdomains[i].size();
00664 minsize=std::min(minsize, subdomains[i].size());
00665 maxsize=std::max(maxsize, subdomains[i].size());
00666 }
00667 std::cout<<"Subdomain size: min="<<minsize<<" max="<<maxsize<<" avg="<<(sum/subdomains.size())
00668 <<" no="<<subdomains.size()<<std::endl;
00669
00670
00671
00672 }
00673 Vector subdomains;
00674 };
00675
00676
00677 template<class M, class X, class TM, class TA>
00678 struct ConstructionTraits<SeqOverlappingSchwarz<M,X,TM,TA> >
00679 {
00680 typedef ConstructionArgs<SeqOverlappingSchwarz<M,X,TM,TA> > Arguments;
00681
00682 static inline SeqOverlappingSchwarz<M,X,TM,TA>* construct(Arguments& args)
00683 {
00684 return new SeqOverlappingSchwarz<M,X,TM,TA>(args.getMatrix(),
00685 args.getSubDomains(),
00686 args.getArgs().relaxationFactor);
00687 }
00688
00689 static void deconstruct(SeqOverlappingSchwarz<M,X,TM,TA>* schwarz)
00690 {
00691 delete schwarz;
00692 }
00693 };
00694
00695 #endif
00696 }
00697 }
00698
00699
00700
00701 #endif