amg.hh

Go to the documentation of this file.
00001 // $Id: amg.hh 764 2007-05-03 12:39:19Z 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 namespace Dune
00013 {
00014   namespace Amg
00015   {
00036     template<class M, class X, class S, class PI=SequentialInformation,
00037              class A=std::allocator<X> >
00038     class AMG : public Preconditioner<X,X>
00039     {
00040     public:
00042       typedef M Operator;
00049       typedef PI ParallelInformation;
00051       typedef MatrixHierarchy<M, ParallelInformation, A> OperatorHierarchy;
00053       typedef typename OperatorHierarchy::ParallelInformationHierarchy ParallelInformationHierarchy;
00054 
00056       typedef X Domain;
00058       typedef X Range;
00060       typedef InverseOperator<X,X> CoarseSolver;
00066       typedef S Smoother;
00067   
00069       typedef typename SmootherTraits<Smoother>::Arguments SmootherArgs;
00070       
00071       enum {
00073         category = S::category
00074       };
00075 
00086       AMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver, 
00087           const SmootherArgs& smootherArgs, std::size_t gamma,
00088           std::size_t smoothingSteps);
00089 
00103       template<class C>
00104       AMG(const Operator& fineOperator, const C& criterion,
00105           const SmootherArgs& smootherArgs, std::size_t gamma=1,
00106           std::size_t smoothingSteps=2, const ParallelInformation& pinfo=ParallelInformation());
00107 
00108       ~AMG();
00109 
00111       void pre(Domain& x, Range& b);
00112 
00114       void apply(Domain& v, const Range& d);
00115       
00117       void post(Domain& x);
00118       
00119     private:
00121       void mgc(typename Hierarchy<Smoother,A>::Iterator& smoother,
00122                typename OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator& matrix,
00123                typename ParallelInformationHierarchy::Iterator& pinfo,
00124                typename OperatorHierarchy::AggregatesMapList::const_iterator& aggregates,
00125                typename Hierarchy<Domain,A>::Iterator& lhs,
00126                typename Hierarchy<Domain,A>::Iterator& update,
00127                typename Hierarchy<Range,A>::Iterator& rhs);
00128 
00129       //      void setupIndices(typename Matrix::ParallelIndexSet& indices, const Matrix& matrix);
00130       
00132       const OperatorHierarchy* matrices_;
00134       SmootherArgs smootherArgs_;
00136       Hierarchy<Smoother,A> smoothers_;
00138       CoarseSolver* solver_;
00140       Hierarchy<Range,A>* rhs_;
00142       Hierarchy<Domain,A>* lhs_;
00144       Hierarchy<Domain,A>* update_;
00146       typedef ScalarProductChooser<X,PI,M::category> ScalarProductChooser;
00148       typedef typename ScalarProductChooser::ScalarProduct ScalarProduct;
00150       ScalarProduct* scalarProduct_;
00152       std::size_t gamma_;
00154       std::size_t steps_;
00155       std::size_t level;
00156       bool buildHierarchy_;
00157       Smoother *coarseSmoother_;
00158     };
00159 
00160     template<class M, class X, class S, class P, class A>
00161     AMG<M,X,S,P,A>::AMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver, 
00162                          const SmootherArgs& smootherArgs,
00163                          std::size_t gamma, std::size_t smoothingSteps)
00164       : matrices_(&matrices), smootherArgs_(smootherArgs),
00165         smoothers_(), solver_(&coarseSolver), scalarProduct_(0),
00166         gamma_(gamma), steps_(smoothingSteps), buildHierarchy_(false)
00167     {
00168       assert(matrices_->isBuilt());
00169       //printMatrix(matrices_->finest());
00170       
00171       // build the necessary smoother hierarchies
00172       matrices_->coarsenSmoother(smoothers_, smootherArgs_);
00173       
00174     }
00175 
00176     template<class M, class X, class S, class P, class A>
00177     template<class C>
00178     AMG<M,X,S,P,A>::AMG(const Operator& matrix,
00179                         const C& criterion,
00180                         const SmootherArgs& smootherArgs,
00181                         std::size_t gamma, std::size_t smoothingSteps,
00182                         const P& pinfo)
00183       : smootherArgs_(smootherArgs),
00184         smoothers_(), scalarProduct_(0), gamma_(gamma),
00185         steps_(smoothingSteps), buildHierarchy_(true)
00186     {
00187       IsTrue<static_cast<int>(M::category)==static_cast<int>(S::category)>::yes();
00188       IsTrue<static_cast<int>(P::category)==static_cast<int>(S::category)>::yes();
00189       OperatorHierarchy* matrices = new OperatorHierarchy(const_cast<Operator&>(matrix), pinfo);
00190             
00191       matrices->template build<typename P::CopyFlags>(criterion);
00192       
00193       matrices_ = matrices;
00194       // build the necessary smoother hierarchies
00195       matrices_->coarsenSmoother(smoothers_, smootherArgs_);
00196             
00197     }
00198     
00199     template<class M, class X, class S, class P, class A>
00200     AMG<M,X,S,P,A>::~AMG()
00201     {
00202       if(buildHierarchy_){
00203         delete matrices_;
00204       }
00205       if(scalarProduct_)
00206         delete scalarProduct_;
00207     }
00208     
00209     template<class M, class X, class S, class P, class A>
00210     void AMG<M,X,S,P,A>::pre(Domain& x, Range& b)
00211     {
00212       Range* copy = new Range(b);
00213       rhs_ = new Hierarchy<Range,A>(*copy);
00214       Domain* dcopy = new Domain(x);
00215       lhs_ = new Hierarchy<Domain,A>(*dcopy);
00216       dcopy = new Domain(x);
00217       update_ = new Hierarchy<Domain,A>(*dcopy);
00218       matrices_->coarsenVector(*rhs_);
00219       matrices_->coarsenVector(*lhs_);
00220       matrices_->coarsenVector(*update_);
00221       
00222       // Preprocess all smoothers
00223       typedef typename Hierarchy<Smoother,A>::Iterator Iterator;
00224       typedef typename Hierarchy<Range,A>::Iterator RIterator;
00225       typedef typename Hierarchy<Domain,A>::Iterator DIterator;
00226       Iterator coarsest = smoothers_.coarsest();
00227       Iterator smoother = smoothers_.finest();
00228       RIterator rhs = rhs_->finest();
00229       DIterator lhs = lhs_->finest();
00230       
00231       if(rhs!=rhs_->coarsest()){
00232         for(; smoother != coarsest; ++smoother, ++lhs, ++rhs)
00233           smoother->pre(*lhs,*rhs);
00234         smoother->pre(*lhs,*rhs);
00235       }
00236       
00237       // The precondtioner might change x and b. So we have to 
00238       // copy the changes to the original vectors.
00239       x = *lhs_->finest();
00240       b = *rhs_->finest();
00241       
00242       if(buildHierarchy_){
00243         // Create the coarse Solver
00244         SmootherArgs sargs(smootherArgs_);
00245         sargs.iterations = 1;
00246         
00247         typename ConstructionTraits<Smoother>::Arguments cargs;
00248         cargs.setArgs(sargs);
00249         cargs.setMatrix(matrices_->matrices().coarsest()->getmat());
00250         cargs.setComm(*matrices_->parallelInformation().coarsest());
00251         
00252         coarseSmoother_ = ConstructionTraits<Smoother>::construct(cargs);
00253         scalarProduct_ = ScalarProductChooser::construct(*matrices_->parallelInformation().coarsest());
00254 
00255         solver_ = new BiCGSTABSolver<X>(const_cast<M&>(*matrices_->matrices().coarsest()), 
00256                                   *scalarProduct_, 
00257                                   *coarseSmoother_, 1E-2, 10000, 0);
00258       }
00259     }
00260 
00262     template<class M, class X, class S, class P, class A>
00263     void AMG<M,X,S,P,A>::apply(Domain& v, const Range& d)
00264     {
00265       typename Hierarchy<Smoother,A>::Iterator smoother = smoothers_.finest();
00266       typename OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator matrix = matrices_->matrices().finest();
00267       typename ParallelInformationHierarchy::Iterator pinfo = matrices_->parallelInformation().finest();
00268       typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates = matrices_->aggregatesMaps().begin();
00269       typename Hierarchy<Domain,A>::Iterator lhs = lhs_->finest();
00270       typename Hierarchy<Domain,A>::Iterator update = update_->finest();
00271       typename Hierarchy<Range,A>::Iterator rhs = rhs_->finest();
00272         
00273       *lhs = v;
00274       *rhs = d;
00275       *update=0;
00276       level=0;
00277       mgc(smoother, matrix, pinfo, aggregates, lhs, update, rhs);
00278       v=*update;
00279     }
00280     
00281     template<class M, class X, class S, class P, class A>
00282     void AMG<M,X,S,P,A>::mgc(typename Hierarchy<Smoother,A>::Iterator& smoother,
00283                              typename OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator& matrix,
00284                              typename ParallelInformationHierarchy::Iterator& pinfo,
00285                              typename OperatorHierarchy::AggregatesMapList::const_iterator& aggregates,
00286                              typename Hierarchy<Domain,A>::Iterator& lhs,
00287                              typename Hierarchy<Domain,A>::Iterator& update,
00288                              typename Hierarchy<Range,A>::Iterator& rhs){
00289       if(matrix == matrices_->matrices().coarsest()){
00290         // Solve directly
00291         InverseOperatorResult res;
00292         pinfo->copyOwnerToAll(*rhs, *rhs);
00293         solver_->apply(*update, *rhs, res);
00294         if(!res.converged)
00295           DUNE_THROW(MathError, "Coarse solver did not converge");
00296       }else{
00297         // presmoothing
00298         *lhs=0;
00299         for(std::size_t i=0; i < steps_; ++i)
00300           smoother->apply(*lhs, *rhs);
00301       
00302         // Accumulate update
00303         *update += *lhs;
00304         // update defect
00305         matrix->applyscaleadd(-1,static_cast<const Domain&>(*lhs), *rhs);
00306 
00307 #ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION                
00308         //restrict defect to coarse level right hand side.
00309         typename Hierarchy<Range,A>::Iterator fineRhs = rhs++;
00310         ++pinfo;
00311         
00312         Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
00313           ::restrict(*(*aggregates), *rhs, static_cast<const Range&>(*fineRhs), *pinfo);
00314         
00315         // prepare coarse system
00316         ++lhs;
00317         ++update;
00318         ++matrix;
00319         ++level;
00320 
00321         if(matrix != matrices_->matrices().coarsest()){
00322           ++smoother;
00323           ++aggregates;
00324         }
00325         
00326         // prepare the update on the next level
00327         *update=0;
00328         
00329         // next level
00330         mgc(smoother, matrix, pinfo, aggregates, lhs, update, rhs);
00331 
00332         if(matrix != matrices_->matrices().coarsest()){
00333           --smoother;
00334           --aggregates;
00335         }
00336         --level;
00337         //prolongate and add the correction (update is in coarse left hand side)
00338         --matrix;
00339 
00340         //typename Hierarchy<Domain,A>::Iterator coarseLhs = lhs--;
00341         --lhs;  
00342         --pinfo;
00343         
00344         
00345         *lhs=0;
00346         Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
00347           ::prolongate(*(*aggregates), *update, *lhs, 1.6);
00348         
00349         --update;
00350         --rhs;
00351 
00352         *update += *lhs;
00353         
00354         // update defect
00355         matrix->applyscaleadd(-1,static_cast<const Domain&>(*lhs), *rhs);
00356 
00357 #endif  
00358         // postsmoothing
00359         *lhs=0;
00360         for(std::size_t i=0; i < steps_; ++i)
00361           smoother->apply(*lhs, *rhs);
00362         *update += *lhs;
00363       }
00364     }
00365 
00367     template<class M, class X, class S, class P, class A>
00368     void AMG<M,X,S,P,A>::post(Domain& x)
00369     {
00370       if(buildHierarchy_){
00371         delete solver_;
00372         ConstructionTraits<Smoother>::deconstruct(coarseSmoother_);
00373       }
00374       
00375       // Postprocess all smoothers
00376       typedef typename Hierarchy<Smoother,A>::Iterator Iterator;
00377       typedef typename Hierarchy<Range,A>::Iterator RIterator;
00378       typedef typename Hierarchy<Domain,A>::Iterator DIterator;
00379       Iterator coarsest = smoothers_.coarsest();
00380       Iterator smoother = smoothers_.finest();
00381       DIterator lhs = lhs_->finest();
00382       
00383       if(lhs!= lhs_->coarsest())
00384         smoother->post(*lhs);
00385       
00386       if(smoother != coarsest)
00387         for(++smoother; smoother != coarsest; ++smoother, ++lhs)
00388           smoother->post(*lhs);
00389 
00390       delete &(*lhs_->finest());
00391       delete lhs_;
00392       delete &(*update_->finest());
00393       delete update_;
00394       delete &(*rhs_->finest());
00395       delete rhs_;
00396     }
00397   } // end namespace Amg
00398 }// end namespace Dune
00399   
00400 #endif

Generated on 12 Dec 2007 with Doxygen (ver 1.5.1)