amg.hh

Go to the documentation of this file.
00001 // $Id: amg.hh 924 2008-07-14 12:41:33Z mblatt $
00002 #ifndef DUNE_AMG_AMG_HH
00003 #define DUNE_AMG_AMG_HH
00004 
00005 #include<memory>
00006 #include<dune/common/exceptions.hh>
00007 #include<dune/istl/paamg/smoother.hh>
00008 #include<dune/istl/paamg/transfer.hh>
00009 #include<dune/istl/paamg/hierarchy.hh>
00010 #include<dune/istl/solvers.hh>
00011 #include<dune/istl/scalarproducts.hh>
00012 #include<dune/istl/superlu.hh>
00013 #include<dune/common/typetraits.hh>
00014 namespace Dune
00015 {
00016   namespace Amg
00017   {
00042     template<class M, class X, class S, class PI=SequentialInformation,
00043              class A=std::allocator<X> >
00044     class AMG : public Preconditioner<X,X>
00045     {
00046     public:
00048       typedef M Operator;
00055       typedef PI ParallelInformation;
00057       typedef MatrixHierarchy<M, ParallelInformation, A> OperatorHierarchy;
00059       typedef typename OperatorHierarchy::ParallelInformationHierarchy ParallelInformationHierarchy;
00060 
00062       typedef X Domain;
00064       typedef X Range;
00066       typedef InverseOperator<X,X> CoarseSolver;
00072       typedef S Smoother;
00073   
00075       typedef typename SmootherTraits<Smoother>::Arguments SmootherArgs;
00076       
00077       enum {
00079         category = S::category
00080       };
00081 
00093       AMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver, 
00094           const SmootherArgs& smootherArgs, std::size_t gamma,
00095           std::size_t preSmoothingSteps,
00096           std::size_t postSmoothingSteps, bool additive=false);
00097 
00112       template<class C>
00113       AMG(const Operator& fineOperator, const C& criterion,
00114           const SmootherArgs& smootherArgs, std::size_t gamma=1,
00115           std::size_t preSmoothingSteps=2, std::size_t postSmoothingSteps=2,
00116           bool additive=false, const ParallelInformation& pinfo=ParallelInformation());
00117 
00118       ~AMG();
00119 
00121       void pre(Domain& x, Range& b);
00122 
00124       void apply(Domain& v, const Range& d);
00125       
00127       void post(Domain& x);
00128       
00129     private:
00131       void mgc(typename Hierarchy<Smoother,A>::Iterator& smoother,
00132                typename OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator& matrix,
00133                typename ParallelInformationHierarchy::Iterator& pinfo,
00134                typename OperatorHierarchy::AggregatesMapList::const_iterator& aggregates,
00135                typename Hierarchy<Domain,A>::Iterator& lhs,
00136                typename Hierarchy<Domain,A>::Iterator& update,
00137                typename Hierarchy<Range,A>::Iterator& rhs);
00138 
00139       void additiveMgc();
00140       
00141       //      void setupIndices(typename Matrix::ParallelIndexSet& indices, const Matrix& matrix);
00142       
00144       const OperatorHierarchy* matrices_;
00146       SmootherArgs smootherArgs_;
00148       Hierarchy<Smoother,A> smoothers_;
00150       CoarseSolver* solver_;
00152       Hierarchy<Range,A>* rhs_;
00154       Hierarchy<Domain,A>* lhs_;
00156       Hierarchy<Domain,A>* update_;
00158       typedef Dune::ScalarProductChooser<X,PI,M::category> ScalarProductChooser;
00160       typedef typename ScalarProductChooser::ScalarProduct ScalarProduct;
00162       ScalarProduct* scalarProduct_;
00164       std::size_t gamma_;
00166       std::size_t preSteps_;
00168       std::size_t postSteps_;
00169       std::size_t level;
00170       bool buildHierarchy_;
00171       bool additive;
00172       Smoother *coarseSmoother_;
00173     };
00174 
00175     template<class M, class X, class S, class P, class A>
00176     AMG<M,X,S,P,A>::AMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver, 
00177                         const SmootherArgs& smootherArgs,
00178                         std::size_t gamma, std::size_t preSmoothingSteps,
00179                         std::size_t postSmoothingSteps, bool additive_)
00180       : matrices_(&matrices), smootherArgs_(smootherArgs),
00181         smoothers_(), solver_(&coarseSolver), scalarProduct_(0),
00182         gamma_(gamma), preSteps_(preSmoothingSteps), postSteps_(postSmoothingSteps), buildHierarchy_(false),
00183         additive(additive_)
00184     {
00185       assert(matrices_->isBuilt());
00186       //printMatrix(matrices_->finest());
00187       
00188       // build the necessary smoother hierarchies
00189       matrices_->coarsenSmoother(smoothers_, smootherArgs_);
00190       
00191     }
00192 
00193     template<class M, class X, class S, class P, class A>
00194     template<class C>
00195     AMG<M,X,S,P,A>::AMG(const Operator& matrix,
00196                         const C& criterion,
00197                         const SmootherArgs& smootherArgs,
00198                         std::size_t gamma, std::size_t preSmoothingSteps,
00199                         std::size_t postSmoothingSteps,
00200                         bool additive_,
00201                         const P& pinfo)
00202       : smootherArgs_(smootherArgs),
00203         smoothers_(), scalarProduct_(0), gamma_(gamma),
00204         preSteps_(preSmoothingSteps), postSteps_(postSmoothingSteps), buildHierarchy_(true),
00205         additive(additive_)
00206     {
00207       dune_static_assert(static_cast<int>(M::category)==static_cast<int>(S::category),
00208                          "Matrix and Solver must match in terms of category!");
00209       dune_static_assert(static_cast<int>(P::category)==static_cast<int>(S::category),
00210                          "Matrix and Solver must match in terms of category!");
00211       OperatorHierarchy* matrices = new OperatorHierarchy(const_cast<Operator&>(matrix), pinfo);
00212             
00213       matrices->template build<typename P::CopyFlags>(criterion);
00214       
00215       matrices_ = matrices;
00216       // build the necessary smoother hierarchies
00217       matrices_->coarsenSmoother(smoothers_, smootherArgs_);
00218             
00219     }
00220     
00221     template<class M, class X, class S, class P, class A>
00222     AMG<M,X,S,P,A>::~AMG()
00223     {
00224       if(buildHierarchy_){
00225         delete matrices_;
00226       }
00227       if(scalarProduct_)
00228         delete scalarProduct_;
00229     }
00230     
00231     template<class M, class X, class S, class P, class A>
00232     void AMG<M,X,S,P,A>::pre(Domain& x, Range& b)
00233     {
00234       Range* copy = new Range(b);
00235       rhs_ = new Hierarchy<Range,A>(*copy);
00236       Domain* dcopy = new Domain(x);
00237       lhs_ = new Hierarchy<Domain,A>(*dcopy);
00238       dcopy = new Domain(x);
00239       update_ = new Hierarchy<Domain,A>(*dcopy);
00240       matrices_->coarsenVector(*rhs_);
00241       matrices_->coarsenVector(*lhs_);
00242       matrices_->coarsenVector(*update_);
00243       
00244       // Preprocess all smoothers
00245       typedef typename Hierarchy<Smoother,A>::Iterator Iterator;
00246       typedef typename Hierarchy<Range,A>::Iterator RIterator;
00247       typedef typename Hierarchy<Domain,A>::Iterator DIterator;
00248       Iterator coarsest = smoothers_.coarsest();
00249       Iterator smoother = smoothers_.finest();
00250       RIterator rhs = rhs_->finest();
00251       DIterator lhs = lhs_->finest();
00252       
00253       if(rhs!=rhs_->coarsest()){
00254         for(; smoother != coarsest; ++smoother, ++lhs, ++rhs)
00255           smoother->pre(*lhs,*rhs);
00256         smoother->pre(*lhs,*rhs);
00257       }
00258       
00259       // The precondtioner might change x and b. So we have to 
00260       // copy the changes to the original vectors.
00261       x = *lhs_->finest();
00262       b = *rhs_->finest();
00263       
00264       if(buildHierarchy_){
00265         // Create the coarse Solver
00266         SmootherArgs sargs(smootherArgs_);
00267         sargs.iterations = 1;
00268         
00269         typename ConstructionTraits<Smoother>::Arguments cargs;
00270         cargs.setArgs(sargs);
00271         cargs.setMatrix(matrices_->matrices().coarsest()->getmat());
00272         cargs.setComm(*matrices_->parallelInformation().coarsest());
00273         
00274         coarseSmoother_ = ConstructionTraits<Smoother>::construct(cargs);
00275         scalarProduct_ = ScalarProductChooser::construct(*matrices_->parallelInformation().coarsest());
00276 
00277 #ifdef HAVE_SUPERLU
00278         if(is_same<ParallelInformation,SequentialInformation>::value){
00279           std::cout<<"Using superlu"<<std::endl;
00280           solver_  = new SuperLU<typename M::matrix_type>(matrices_->matrices().coarsest()->getmat());
00281         }else
00282 #endif
00283         solver_ = new BiCGSTABSolver<X>(const_cast<M&>(*matrices_->matrices().coarsest()), 
00284                                   *scalarProduct_, 
00285                                   *coarseSmoother_, 1E-2, 10000, 0);
00286       }
00287     }
00288 
00290     template<class M, class X, class S, class P, class A>
00291     void AMG<M,X,S,P,A>::apply(Domain& v, const Range& d)
00292     {
00293       if(additive){
00294         *(rhs_->finest())=d;
00295         additiveMgc();
00296         v=*lhs_->finest();
00297       }else{
00298         typename Hierarchy<Smoother,A>::Iterator smoother = smoothers_.finest();
00299         typename OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator matrix = matrices_->matrices().finest();
00300         typename ParallelInformationHierarchy::Iterator pinfo = matrices_->parallelInformation().finest();
00301         typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates = matrices_->aggregatesMaps().begin();
00302         typename Hierarchy<Domain,A>::Iterator lhs = lhs_->finest();
00303         typename Hierarchy<Domain,A>::Iterator update = update_->finest();
00304         typename Hierarchy<Range,A>::Iterator rhs = rhs_->finest();
00305         
00306         *lhs = v;
00307         *rhs = d;
00308         *update=0;
00309         level=0;
00310         mgc(smoother, matrix, pinfo, aggregates, lhs, update, rhs);
00311         v=*update;
00312       }
00313       
00314     }
00315     
00316     template<class M, class X, class S, class P, class A>
00317     void AMG<M,X,S,P,A>::mgc(typename Hierarchy<Smoother,A>::Iterator& smoother,
00318                              typename OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator& matrix,
00319                              typename ParallelInformationHierarchy::Iterator& pinfo,
00320                              typename OperatorHierarchy::AggregatesMapList::const_iterator& aggregates,
00321                              typename Hierarchy<Domain,A>::Iterator& lhs,
00322                              typename Hierarchy<Domain,A>::Iterator& update,
00323                              typename Hierarchy<Range,A>::Iterator& rhs){
00324       if(matrix == matrices_->matrices().coarsest()){
00325         // Solve directly
00326         InverseOperatorResult res;
00327         pinfo->copyOwnerToAll(*rhs, *rhs);
00328         solver_->apply(*update, *rhs, res);
00329         if(!res.converged)
00330           DUNE_THROW(MathError, "Coarse solver did not converge");
00331       }else{
00332         // presmoothing
00333           *lhs=0;
00334         for(std::size_t i=0; i < preSteps_; ++i)
00335           SmootherApplier<S>::preSmooth(*smoother, *lhs, *rhs);
00336       
00337         // Accumulate update
00338         *update += *lhs;
00339         // update defect
00340         matrix->applyscaleadd(-1,static_cast<const Domain&>(*lhs), *rhs);
00341 
00342 #ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION                
00343         //restrict defect to coarse level right hand side.
00344         typename Hierarchy<Range,A>::Iterator fineRhs = rhs++;
00345         ++pinfo;
00346         
00347         Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
00348           ::restrict(*(*aggregates), *rhs, static_cast<const Range&>(*fineRhs), *pinfo);
00349         
00350         // prepare coarse system
00351         ++lhs;
00352         ++update;
00353         ++matrix;
00354         ++level;
00355 
00356         if(matrix != matrices_->matrices().coarsest()){
00357           ++smoother;
00358           ++aggregates;
00359         }
00360         
00361         // prepare the update on the next level
00362         *update=0;
00363         
00364         // next level
00365         mgc(smoother, matrix, pinfo, aggregates, lhs, update, rhs);
00366 
00367         if(matrix != matrices_->matrices().coarsest()){
00368           --smoother;
00369           --aggregates;
00370         }
00371         --level;
00372         //prolongate and add the correction (update is in coarse left hand side)
00373         --matrix;
00374 
00375         //typename Hierarchy<Domain,A>::Iterator coarseLhs = lhs--;
00376         --lhs;  
00377         --pinfo;
00378         
00379         
00380         *lhs=0;
00381         Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
00382           ::prolongate(*(*aggregates), *update, *lhs, 1.6, *pinfo);
00383         
00384         --update;
00385         --rhs;
00386 
00387         *update += *lhs;
00388         
00389 
00390 #endif  
00391         
00392         // update defect
00393         matrix->applyscaleadd(-1,static_cast<const Domain&>(*lhs), *rhs);
00394 
00395         // postsmoothing
00396         *lhs=0;
00397         for(std::size_t i=0; i < postSteps_; ++i)
00398           SmootherApplier<S>::postSmooth(*smoother, *lhs, *rhs);
00399         
00400         *update += *lhs;
00401       }
00402     }
00403 
00404     template<class M, class X, class S, class P, class A>
00405     void AMG<M,X,S,P,A>::additiveMgc(){
00406       
00407       // restrict residual to all levels
00408       typename ParallelInformationHierarchy::Iterator pinfo=matrices_->parallelInformation().finest();
00409       typename Hierarchy<Range,A>::Iterator rhs=rhs_->finest();      
00410       typename Hierarchy<Domain,A>::Iterator lhs = lhs_->finest();
00411       typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates=matrices_->aggregatesMaps().begin();
00412       
00413       for(typename Hierarchy<Range,A>::Iterator fineRhs=rhs++; fineRhs != rhs_->coarsest(); fineRhs=rhs++, ++aggregates){
00414         ++pinfo;
00415         Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
00416           ::restrict(*(*aggregates), *rhs, static_cast<const Range&>(*fineRhs), *pinfo);
00417       }
00418       
00419       // pinfo is invalid, set to coarsest level
00420       //pinfo = matrices_->parallelInformation().coarsest
00421       // calculate correction for all levels
00422       lhs = lhs_->finest();
00423       typename Hierarchy<Smoother,A>::Iterator smoother = smoothers_.finest();
00424       
00425       for(rhs=rhs_->finest(); rhs != rhs_->coarsest(); ++lhs, ++rhs, ++smoother){
00426         // presmoothing
00427         *lhs=0;
00428         smoother->apply(*lhs, *rhs);
00429       }
00430       
00431       // Coarse level solve
00432 #ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION 
00433       InverseOperatorResult res;
00434       pinfo->copyOwnerToAll(*rhs, *rhs);
00435       solver_->apply(*lhs, *rhs, res);
00436       
00437       if(!res.converged)
00438         DUNE_THROW(MathError, "Coarse solver did not converge");
00439 #else
00440       *lhs=0;
00441 #endif
00442       // Prologate and add up corrections from all levels
00443       --pinfo;
00444       --aggregates;
00445       
00446       for(typename Hierarchy<Domain,A>::Iterator coarseLhs = lhs--; coarseLhs != lhs_->finest(); coarseLhs = lhs--, --aggregates, --pinfo){
00447         Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
00448           ::prolongate(*(*aggregates), *coarseLhs, *lhs, 1, *pinfo);
00449       }
00450     }
00451 
00452     
00454     template<class M, class X, class S, class P, class A>
00455     void AMG<M,X,S,P,A>::post(Domain& x)
00456     {
00457       if(buildHierarchy_){
00458         delete solver_;
00459         ConstructionTraits<Smoother>::deconstruct(coarseSmoother_);
00460       }
00461       
00462       // Postprocess all smoothers
00463       typedef typename Hierarchy<Smoother,A>::Iterator Iterator;
00464       typedef typename Hierarchy<Range,A>::Iterator RIterator;
00465       typedef typename Hierarchy<Domain,A>::Iterator DIterator;
00466       Iterator coarsest = smoothers_.coarsest();
00467       Iterator smoother = smoothers_.finest();
00468       DIterator lhs = lhs_->finest();
00469       
00470       if(lhs!= lhs_->coarsest())
00471         smoother->post(*lhs);
00472       
00473       if(smoother != coarsest)
00474         for(++smoother; smoother != coarsest; ++smoother, ++lhs)
00475           smoother->post(*lhs);
00476 
00477       delete &(*lhs_->finest());
00478       delete lhs_;
00479       delete &(*update_->finest());
00480       delete update_;
00481       delete &(*rhs_->finest());
00482       delete rhs_;
00483     }
00484   } // end namespace Amg
00485 }// end namespace Dune
00486   
00487 #endif

Generated on Sun Nov 15 22:29:35 2009 for dune-istl by  doxygen 1.5.6