hierarchy.hh

Go to the documentation of this file.
00001 // $Id: hierarchy.hh 1151 2010-01-20 14:43:24Z mblatt $
00002 #ifndef DUNE_AMGHIERARCHY_HH
00003 #define DUNE_AMGHIERARCHY_HH
00004 
00005 #include<list>
00006 #include<memory>
00007 #include<limits>
00008 #include<algorithm>
00009 #include"aggregates.hh"
00010 #include"graph.hh"
00011 #include"galerkin.hh"
00012 #include"renumberer.hh"
00013 #include"graphcreator.hh"
00014 #include<dune/common/stdstreams.hh>
00015 #include<dune/common/timer.hh>
00016 #include<dune/common/tuples.hh>
00017 #include<dune/istl/bvector.hh>
00018 #include<dune/istl/indexset.hh>
00019 #include<dune/istl/matrixutils.hh>
00020 #include<dune/istl/matrixredistribute.hh>
00021 #include<dune/istl/paamg/dependency.hh>
00022 #include<dune/istl/paamg/graph.hh>
00023 #include<dune/istl/paamg/indicescoarsener.hh>
00024 #include<dune/istl/paamg/globalaggregates.hh>
00025 #include<dune/istl/paamg/construction.hh>
00026 #include<dune/istl/paamg/smoother.hh>
00027 #include<dune/istl/paamg/transfer.hh>
00028 
00029 namespace Dune
00030 {
00031   namespace Amg
00032   {
00050     template<typename T, typename A=std::allocator<T> >
00051     class Hierarchy
00052     {
00053     public:
00057       typedef T MemberType;
00058 
00059       template<typename T1, typename T2>
00060       class LevelIterator;
00061 
00062     private:
00066       struct Element
00067       {
00068         friend class LevelIterator<Hierarchy<T,A>, T>;
00069         friend class LevelIterator<const Hierarchy<T,A>, const T>;
00070 
00072         Element* coarser_;
00073         
00075         Element* finer_;
00076         
00078         MemberType* element_;
00079         
00081         MemberType* redistributed_;
00082       };
00083     public:
00084 //       enum{
00085 //      /**
00086 //       * @brief If true only the method addCoarser will be usable
00087 //       * otherwise only the method addFiner will be usable.
00088 //       */
00089 //      coarsen = b
00090 //        };
00091       
00095       typedef typename A::template rebind<Element>::other Allocator;
00096 
00097       typedef typename ConstructionTraits<T>::Arguments Arguments;
00098       
00103       Hierarchy(MemberType& first);
00104       
00108       Hierarchy();
00109 
00114       void addCoarser(Arguments& args);
00115 
00116       void addRedistributedOnCoarsest(T* t);
00117       
00122       void addFiner(Arguments& args);
00123 
00130       template<class C, class T1>
00131       class LevelIterator 
00132         : public BidirectionalIteratorFacade<LevelIterator<C,T1>,T1,T1&>
00133       {
00134         friend class LevelIterator<typename remove_const<C>::type,
00135                                    typename remove_const<T1>::type >;
00136         friend class LevelIterator<const typename remove_const<C>::type, 
00137                                    const typename remove_const<T1>::type >;
00138 
00139       public:
00141         LevelIterator()
00142           : element_(0)
00143         {}
00144         
00145         LevelIterator(Element* element)
00146           : element_(element)
00147         {}
00148         
00150         LevelIterator(const LevelIterator<typename remove_const<C>::type,
00151                       typename remove_const<T1>::type>& other)
00152           : element_(other.element_)
00153         {}
00154         
00156         LevelIterator(const LevelIterator<const typename remove_const<C>::type,
00157                       const typename remove_const<T1>::type>& other)
00158           : element_(other.element_)
00159         {}
00160 
00164         bool equals(const LevelIterator<typename remove_const<C>::type, 
00165                     typename remove_const<T1>::type>& other) const
00166         {
00167           return element_ == other.element_;
00168         }
00169         
00173         bool equals(const LevelIterator<const typename remove_const<C>::type, 
00174                     const typename remove_const<T1>::type>& other) const
00175         {
00176           return element_ == other.element_;
00177         }
00178 
00180         T1& dereference() const
00181         {
00182           return *(element_->element_);
00183         }
00184         
00186         void increment()
00187         {
00188           element_ = element_->coarser_;
00189         }
00190         
00192         void decrement()
00193         {
00194           element_ = element_->finer_;
00195         }
00196         
00201         bool isRedistributed() const
00202         {
00203           return element_->redistributed_;
00204         }
00205         
00210         T1& getRedistributed() const
00211         {
00212           assert(element_->redistributed_);
00213           return *element_->redistributed_;
00214         }
00215         void addRedistributed(T1* t)
00216         {
00217           element_->redistributed_ = t;
00218         }
00219 
00220       private:
00221         Element* element_;
00222       };
00223       
00225       typedef LevelIterator<Hierarchy<T,A>,T> Iterator;
00226 
00228       typedef LevelIterator<const Hierarchy<T,A>, const T> ConstIterator;
00229 
00234       Iterator finest();
00235       
00240       Iterator coarsest();
00241       
00242       
00247       ConstIterator finest() const;
00248       
00253       ConstIterator coarsest() const;
00254 
00259       std::size_t levels() const;
00260       
00262       ~Hierarchy();
00263       
00264     private:
00266       Element* finest_;
00268       Element* coarsest_;
00270       Element* nonAllocated_;
00272       Allocator allocator_;
00274       int levels_;
00275     };
00276     
00283     template<class M, class PI, class A=std::allocator<M> >
00284     class MatrixHierarchy
00285     {
00286     public:
00288       typedef M MatrixOperator;
00289       
00291       typedef typename MatrixOperator::matrix_type Matrix;
00292 
00294       typedef PI ParallelInformation;
00295       
00297       typedef A Allocator;
00298 
00300       typedef Dune::Amg::AggregatesMap<typename MatrixGraph<Matrix>::VertexDescriptor> AggregatesMap;
00301 
00303       typedef Dune::Amg::Hierarchy<MatrixOperator,Allocator> ParallelMatrixHierarchy;
00304 
00306       typedef Dune::Amg::Hierarchy<ParallelInformation,Allocator> ParallelInformationHierarchy;
00307 
00309       typedef typename Allocator::template rebind<AggregatesMap*>::other AAllocator;
00310       
00312       typedef std::list<AggregatesMap*,AAllocator> AggregatesMapList;
00313 
00315       typedef RedistributeInformation<ParallelInformation> RedistributeInfoType;
00316       
00318       typedef std::list<RedistributeInfoType,AAllocator> RedistributeInfoList;
00319 
00325       MatrixHierarchy(const MatrixOperator& fineMatrix,
00326                       const ParallelInformation& pinfo=ParallelInformation());
00327       
00328       
00329       ~MatrixHierarchy();
00330       
00336       template<typename O, typename T>
00337       void build(const T& criterion);
00338 
00346       template<class F>
00347       void recalculateGalerkin(const F& copyFlags);
00348 
00353       template<class V, class TA>
00354       void coarsenVector(Hierarchy<BlockVector<V,TA> >& hierarchy) const;
00355 
00361       template<class S, class TA>
00362       void coarsenSmoother(Hierarchy<S,TA>& smoothers, 
00363                            const typename SmootherTraits<S>::Arguments& args) const;
00364       
00369       std::size_t levels() const;
00370       
00375       std::size_t maxlevels() const;
00376       
00377       bool hasCoarsest() const;
00378       
00383       bool isBuilt() const;
00384 
00389       const ParallelMatrixHierarchy& matrices() const;
00390       
00395       const ParallelInformationHierarchy& parallelInformation() const;
00396       
00401       const AggregatesMapList& aggregatesMaps() const;
00402       
00408       const RedistributeInfoList& redistributeInformation() const;
00409       
00410 
00411       typename MatrixOperator::field_type getProlongationDampingFactor() const
00412       {
00413         return prolongDamp_;
00414       }
00415       
00426       void getCoarsestAggregatesOnFinest(std::vector<std::size_t>& data) const;
00427       
00428     private:
00429       typedef typename ConstructionTraits<MatrixOperator>::Arguments MatrixArgs;
00431       AggregatesMapList aggregatesMaps_;
00433       RedistributeInfoList redistributes_;
00435       ParallelMatrixHierarchy matrices_;
00437       ParallelInformationHierarchy parallelInformation_;
00438 
00440       bool built_;
00441 
00443       int maxlevels_;
00444       
00445       typename MatrixOperator::field_type prolongDamp_;
00446       
00450       template<class Matrix, bool print>
00451       struct MatrixStats
00452       {
00453         
00457         static void stats(const Matrix& matrix)
00458         {}
00459       };
00460       
00461       template<class Matrix>
00462       struct MatrixStats<Matrix,true>
00463       {
00464         struct calc
00465         {
00466           typedef typename Matrix::size_type size_type;
00467           typedef typename Matrix::row_type matrix_row;
00468           
00469           calc()
00470           {
00471             min=std::numeric_limits<size_type>::max();
00472             max=0;
00473             sum=0;
00474           }
00475           
00476           void operator()(const matrix_row& row)
00477           {
00478             min=std::min(min, row.size());
00479             max=std::max(max, row.size());
00480             sum += row.size();
00481           }
00482           
00483           size_type min;
00484           size_type max;
00485           size_type sum;
00486         };
00490         static void stats(const Matrix& matrix)
00491         {
00492           calc c= for_each(matrix.begin(), matrix.end(), calc());
00493           dinfo<<"Matrix row: min="<<c.min<<" max="<<c.max
00494                <<" average="<<static_cast<double>(c.sum)/matrix.N()
00495                <<std::endl;
00496         }
00497       };
00498     };
00499 
00503     enum AccumulationMode{
00509       noAccu = 0,
00515       atOnceAccu=1,
00519       successiveAccu=2
00520     };
00521           
00522 
00526     template<class T>
00527     class CoarsenCriterion : public T
00528     {
00529     public:
00534       typedef T DependencyCriterion;
00535       
00539       void setMaxLevel(int l)
00540       {
00541         maxLevel_ = l;
00542       }
00546       int maxLevel() const
00547       {
00548         return maxLevel_;
00549       }
00550       
00554       void setCoarsenTarget(int nodes)
00555       {
00556         coarsenTarget_ = nodes;
00557       }
00558       
00562       int coarsenTarget() const
00563       {
00564         return coarsenTarget_;
00565       }
00566       
00572       void setMinCoarsenRate(double rate)
00573       {
00574         minCoarsenRate_ = rate;
00575       }
00576       
00580       double minCoarsenRate() const
00581       {
00582         return minCoarsenRate_;
00583       }
00584 
00588       AccumulationMode accumulate() const
00589       {
00590         return accumulate_;
00591       }
00595       void setAccumulate(AccumulationMode accu)
00596       {
00597         accumulate_=accu;
00598       }
00599 
00600       void setAccumulate(bool accu){
00601         accumulate_=accu?successiveAccu:noAccu;
00602       }
00608        void setProlongationDampingFactor(double d)
00609       {
00610         dampingFactor_ = d;
00611       }
00612 
00618       double getProlongationDampingFactor() const
00619       {
00620         return dampingFactor_;
00621       }
00632       CoarsenCriterion(int maxLevel=100, int coarsenTarget=1000, double minCoarsenRate=1.2,
00633                        double prolongDamp=1.6, AccumulationMode accumulate=successiveAccu)
00634         : T(), maxLevel_(maxLevel), coarsenTarget_(coarsenTarget), minCoarsenRate_(minCoarsenRate),
00635           dampingFactor_(prolongDamp), accumulate_( accumulate)
00636       {}
00637       
00638     private:
00642       int maxLevel_;
00646       int coarsenTarget_;
00650       double minCoarsenRate_;
00654       double dampingFactor_;
00659       AccumulationMode accumulate_;
00660     };
00661     
00662     template<typename M, typename C1>
00663     bool repartitionAndDistributeMatrix(const M& origMatrix, M& newMatrix, 
00664                                         SequentialInformation& origSequentialInformationomm, 
00665                                         SequentialInformation*& newComm, 
00666                                         RedistributeInformation<SequentialInformation>& ri,
00667                                         int nparts, C1& criterion)
00668     {
00669       DUNE_THROW(NotImplemented, "Redistribution does not make sense in sequential code!");
00670     }
00671     
00672     
00673     template<typename M, typename C, typename C1>
00674     bool repartitionAndDistributeMatrix(const M& origMatrix, M& newMatrix, C& origComm, C*& newComm, 
00675                                         RedistributeInformation<C>& ri,
00676                                         int nparts, C1& criterion)
00677   {
00678     typedef Dune::Amg::MatrixGraph<const M> MatrixGraph;
00679     typedef Dune::Amg::PropertiesGraph<MatrixGraph,
00680       VertexProperties,
00681       EdgeProperties,
00682       IdentityMap,
00683       IdentityMap> PropertiesGraph;
00684     MatrixGraph graph(origMatrix);
00685     PropertiesGraph pgraph(graph);
00686     buildDependency(pgraph, origMatrix, criterion);
00687     
00688 #ifdef DEBUG_REPART
00689     if(origComm.communicator().rank()==0)
00690       std::cout<<"Original matrix"<<std::endl;
00691     origComm.communicator().barrier();
00692   printGlobalSparseMatrix(origMatrix, origComm, std::cout);
00693 #endif
00694     bool existentOnRedist=Dune::graphRepartition(pgraph, origComm, nparts,
00695                                                  newComm, ri.getInterface());
00696     ri.setSetup();
00697 #ifdef DEBUG_REPART
00698     ri.checkInterface(origComm.indexSet(), newComm->indexSet(), origComm.communicator());
00699 #endif
00700     redistributeMatrix(const_cast<M&>(origMatrix), newMatrix, origComm, *newComm, ri);
00701 #if !HAVE_PARMETIS && HAVE_MPI
00702     #warning Parmetis is not installed or used. Did you use the parmetis flags? It is stronly recommend to use parallel AMG with parmetis.
00703 #endif
00704 
00705 #ifdef DEBUG_REPART
00706     if(origComm.communicator().rank()==0)
00707       std::cout<<"Original matrix"<<std::endl;
00708     origComm.communicator().barrier();
00709     if(newComm->communicator().size()>0)
00710       printGlobalSparseMatrix(newMatrix, *newComm, std::cout);
00711     origComm.communicator().barrier();
00712 #endif
00713 
00714     return existentOnRedist;
00715     
00716   }
00717   
00718   template<typename M>
00719   bool repartitionAndDistributeMatrix(M& origMatrix, M& newMatrix, 
00720                                       SequentialInformation& origComm, 
00721                                       SequentialInformation& newComm, 
00722                                       RedistributeInformation<SequentialInformation>& ri)
00723   {
00724     return true;
00725     
00726   }
00727 
00728     template<class M, class IS, class A>
00729     MatrixHierarchy<M,IS,A>::MatrixHierarchy(const MatrixOperator& fineOperator,
00730                                              const ParallelInformation& pinfo)
00731       : matrices_(const_cast<MatrixOperator&>(fineOperator)),
00732         parallelInformation_(const_cast<ParallelInformation&>(pinfo))
00733     {
00734       dune_static_assert((static_cast<int>(MatrixOperator::category) == static_cast<int>(SolverCategory::sequential) ||
00735         static_cast<int>(MatrixOperator::category) == static_cast<int>(SolverCategory::overlapping)),
00736                          "MatrixOperator must be of category sequential or overlapping");
00737       dune_static_assert((static_cast<int>(MatrixOperator::category) == static_cast<int>(ParallelInformation::category)),
00738                          "MatrixOperator and ParallelInformation must belong to the same category!");
00739     }
00740     
00741     template<class M, class IS, class A>
00742     template<typename O, typename T>
00743     void MatrixHierarchy<M,IS,A>::build(const T& criterion)
00744     {
00745 
00746       prolongDamp_ = criterion.getProlongationDampingFactor();
00747       typedef O OverlapFlags;
00748       typedef typename ParallelMatrixHierarchy::Iterator MatIterator;
00749       typedef typename ParallelInformationHierarchy::Iterator PInfoIterator;
00750       
00751       GalerkinProduct<ParallelInformation> productBuilder;
00752       MatIterator mlevel = matrices_.finest();
00753       MatrixStats<typename M::matrix_type,MINIMAL_DEBUG_LEVEL<=INFO_DEBUG_LEVEL>::stats(mlevel->getmat());
00754       
00755       PInfoIterator infoLevel = parallelInformation_.finest();
00756       std::size_t finenonzeros=countNonZeros(mlevel->getmat());
00757       finenonzeros = infoLevel->communicator().sum(finenonzeros);
00758       std::size_t allnonzeros = finenonzeros;
00759 
00760       
00761       int procs = infoLevel->communicator().size();
00762       int level = 0;
00763       int rank = 0;
00764       int unknowns = mlevel->getmat().N();
00765 
00766       unknowns = infoLevel->communicator().sum(unknowns);
00767       infoLevel->buildGlobalLookup(mlevel->getmat().N());
00768       redistributes_.push_back(RedistributeInfoType());
00769 
00770       for(; level < criterion.maxLevel(); ++level, ++mlevel){
00771         assert(matrices_.levels()==redistributes_.size());
00772         rank = infoLevel->communicator().rank();
00773         if(rank==0 && criterion.debugLevel()>1)
00774           std::cout<<"Level "<<level<<" has "<<unknowns<<" unknowns, "<<unknowns/infoLevel->communicator().size()
00775                    <<" unknowns per proc (procs="<<infoLevel->communicator().size()<<")"<<std::endl;
00776 
00777         MatrixOperator* matrix=&(*mlevel);
00778         ParallelInformation* info =&(*infoLevel);
00779 
00780         if((
00781 #ifdef ENABLE_PARMETIS 
00782               (HAVE_PARMETIS && criterion.accumulate()==successiveAccu)
00783 #else
00784               false
00785 #endif
00786              || (criterion.accumulate()==atOnceAccu
00787                  && unknowns < 30*infoLevel->communicator().size()))
00788            && infoLevel->communicator().size()>1 && 
00789            unknowns/infoLevel->communicator().size() <= criterion.coarsenTarget())
00790           {
00791             // accumulate to fewer processors
00792             Matrix* redistMat= new Matrix();
00793             ParallelInformation* redistComm=0;
00794             std::size_t nodomains = unknowns/(criterion.minAggregateSize()
00795                                      *criterion.coarsenTarget());
00796             if(nodomains<=criterion.minAggregateSize()/2)
00797               nodomains=1;
00798             
00799             bool existentOnNextLevel = 
00800               repartitionAndDistributeMatrix(mlevel->getmat(), *redistMat, *infoLevel,
00801                                              redistComm, redistributes_.back(), nodomains,
00802                                              criterion);
00803             int unknowns = redistMat->N();
00804             unknowns = infoLevel->communicator().sum(unknowns);
00805             if(redistComm->communicator().rank()==0 && criterion.debugLevel()>1)
00806               std::cout<<"Level "<<level<<" (redistributed) has "<<unknowns<<" unknowns, "<<unknowns/redistComm->communicator().size()
00807                        <<" unknowns per proc (procs="<<redistComm->communicator().size()<<")"<<std::endl;
00808             MatrixArgs args(*redistMat, *redistComm);
00809             mlevel.addRedistributed(ConstructionTraits<MatrixOperator>::construct(args));
00810             assert(mlevel.isRedistributed());
00811             infoLevel.addRedistributed(redistComm);
00812             infoLevel->freeGlobalLookup();
00813             
00814             if(!existentOnNextLevel)
00815               // We do not hold any data on the redistributed partitioning
00816               break;
00817             
00818             // Work on the redistributed Matrix from now on
00819             matrix = &(mlevel.getRedistributed());
00820             info = &(infoLevel.getRedistributed());
00821             info->buildGlobalLookup(matrix->getmat().N());
00822           }
00823         
00824         rank = info->communicator().rank();
00825         procs = info->communicator().size();
00826         if(unknowns <= criterion.coarsenTarget())
00827            // No further coarsening needed
00828            break;
00829 
00830         typedef PropertiesGraphCreator<MatrixOperator> GraphCreator;
00831         typedef typename GraphCreator::PropertiesGraph PropertiesGraph;
00832         typedef typename GraphCreator::MatrixGraph MatrixGraph;
00833         typedef typename GraphCreator::GraphTuple  GraphTuple;
00834 
00835         typedef typename PropertiesGraph::VertexDescriptor Vertex;
00836         
00837         std::vector<bool> excluded(matrix->getmat().N(), false);
00838 
00839         GraphTuple graphs = GraphCreator::create(*matrix, excluded, *info, OverlapFlags());
00840         
00841         AggregatesMap* aggregatesMap=new AggregatesMap(get<1>(graphs)->maxVertex()+1);
00842 
00843         aggregatesMaps_.push_back(aggregatesMap);
00844 
00845         Timer watch;
00846         watch.reset();
00847         int noAggregates, isoAggregates, oneAggregates, skippedAggregates;
00848         
00849         tie(noAggregates, isoAggregates, oneAggregates, skippedAggregates) =    
00850           aggregatesMap->buildAggregates(matrix->getmat(), *(get<1>(graphs)), criterion);
00851 
00852 #ifdef TEST_AGGLO
00853         {
00854         // calculate size of local matrix in the distributed direction
00855         int start, end, overlapStart, overlapEnd;
00856         int n = UNKNOWNS/procs; // number of unknowns per process
00857         int bigger = UNKNOWNS%procs; // number of process with n+1 unknows
00858         int procs=info->communicator().rank();
00859 
00860         // Compute owner region
00861         if(rank<bigger){
00862           start = rank*(n+1);
00863           end   = (rank+1)*(n+1);
00864         }else{
00865           start = bigger + rank * n;
00866           end   = bigger + (rank + 1) * n;
00867         }
00868   
00869         // Compute overlap region
00870         if(start>0)
00871           overlapStart = start - 1;
00872         else
00873           overlapStart = start;
00874         
00875         if(end<UNKNOWNS)
00876           overlapEnd = end + 1;
00877         else
00878           overlapEnd = end;
00879   
00880         assert((UNKNOWNS)*(overlapEnd-overlapStart)==aggregatesMap->noVertices());
00881         for(int j=0; j< UNKNOWNS; ++j)
00882           for(int i=0; i < UNKNOWNS; ++i)
00883             {
00884               if(i>=overlapStart && i<overlapEnd)
00885                 {
00886                   int no = (j/2)*((UNKNOWNS)/2)+i/2;
00887                   (*aggregatesMap)[j*(overlapEnd-overlapStart)+i-overlapStart]=no;
00888                 }
00889             }
00890         }
00891 #endif
00892         noAggregates = info->communicator().sum(noAggregates);
00893 
00894 #ifdef TEST_AGGLO
00895         noAggregates=((UNKNOWNS)/2)*((UNKNOWNS)/2);
00896 #endif
00897 
00898         if(criterion.debugLevel()>2 && rank==0)
00899           std::cout << "Building "<<noAggregates<<" aggregates took "<<watch.elapsed()<<" seconds."<<std::endl; 
00900 
00901         if(!noAggregates || unknowns/noAggregates<criterion.minCoarsenRate())
00902     {
00903             if(rank==0)
00904         {
00905               if(noAggregates)
00906             std::cerr << "Stopped coarsening because of rate breakdown "<<unknowns/noAggregates<<"<"
00907                   <<criterion.minCoarsenRate()<<std::endl;
00908               else
00909             std::cerr<< "Could not build any aggregates. Probably no connected nodes."<<std::endl;
00910         }
00911             aggregatesMap->free();
00912             delete aggregatesMap;
00913             aggregatesMaps_.pop_back();
00914             break;
00915     }
00916         unknowns =  noAggregates;
00917                         
00918         parallelInformation_.addCoarser(info->communicator());
00919 
00920         ++infoLevel; // parallel information on coarse level
00921                 
00922         typename PropertyMapTypeSelector<VertexVisitedTag,PropertiesGraph>::Type visitedMap = 
00923           get(VertexVisitedTag(), *(get<1>(graphs)));
00924   
00925         watch.reset();
00926         int aggregates = IndicesCoarsener<ParallelInformation,OverlapFlags>
00927 	  ::coarsen(*info,
00928                     *(get<1>(graphs)),
00929                     visitedMap,
00930                     *aggregatesMap,
00931                     *infoLevel);
00932         
00933         GraphCreator::free(graphs);
00934         
00935         if(criterion.debugLevel()>2){
00936           if(rank==0)
00937             std::cout<<"Coarsening of index sets took "<<watch.elapsed()<<" seconds."<<std::endl;
00938         }
00939         
00940         watch.reset();
00941 
00942         infoLevel->buildGlobalLookup(aggregates);
00943         AggregatesPublisher<Vertex,OverlapFlags,ParallelInformation>::publish(*aggregatesMap,
00944                                                                               *info,
00945                                                                               infoLevel->globalLookup());
00946         
00947                 
00948         if(criterion.debugLevel()>2){
00949           if(rank==0)
00950             std::cout<<"Communicating global aggregate numbers took "<<watch.elapsed()<<" seconds."<<std::endl;
00951         }
00952 
00953         watch.reset();
00954         std::vector<bool>& visited=excluded;
00955         
00956         typedef std::vector<bool>::iterator Iterator;
00957         typedef IteratorPropertyMap<Iterator, IdentityMap> VisitedMap2;
00958         Iterator end = visited.end();
00959         for(Iterator iter= visited.begin(); iter != end; ++iter)
00960           *iter=false;
00961         
00962         VisitedMap2 visitedMap2(visited.begin(), Dune::IdentityMap());
00963 
00964         typename MatrixOperator::matrix_type* coarseMatrix;
00965 
00966         coarseMatrix = productBuilder.build(matrix->getmat(), *(get<0>(graphs)), visitedMap2, 
00967                                             *info, 
00968                                             *aggregatesMap,
00969                                             aggregates,
00970                                             OverlapFlags());
00971         
00972         info->freeGlobalLookup();
00973         
00974         delete get<0>(graphs);
00975         productBuilder.calculate(matrix->getmat(), *aggregatesMap, *coarseMatrix, *infoLevel, OverlapFlags());
00976         
00977         if(criterion.debugLevel()>2){
00978           if(rank==0)
00979             std::cout<<"Calculation of Galerkin product took "<<watch.elapsed()<<" seconds."<<std::endl;
00980         }
00981         
00982         std::size_t nonzeros = countNonZeros(*coarseMatrix);
00983         allnonzeros += infoLevel->communicator().sum(nonzeros);
00984         MatrixArgs args(*coarseMatrix, *infoLevel);
00985         
00986         matrices_.addCoarser(args);
00987         redistributes_.push_back(RedistributeInfoType());
00988       }
00989       
00990 
00991       infoLevel->freeGlobalLookup();
00992       
00993       built_=true;
00994       AggregatesMap* aggregatesMap=new AggregatesMap(0);
00995       aggregatesMaps_.push_back(aggregatesMap);
00996 
00997       if(criterion.debugLevel()>0){
00998         if(level==criterion.maxLevel()){
00999           int unknowns = mlevel->getmat().N();
01000           unknowns = infoLevel->communicator().sum(unknowns);
01001           if(rank==0 && criterion.debugLevel()>1)
01002             std::cout<<"Level "<<level<<" has "<<unknowns<<" unknowns, "<<unknowns/infoLevel->communicator().size()
01003                      <<" unknowns per proc (procs="<<infoLevel->communicator().size()<<")"<<std::endl;
01004         }
01005       }
01006 
01007       if(criterion.accumulate() && !redistributes_.back().isSetup() && 
01008          infoLevel->communicator().size()>1){ 
01009 #if HAVE_MPI && !HAVE_PARMETIS
01010         if(criterion.accumulate()==successiveAccu &&
01011            infoLevel->communicator().rank()==0)
01012           std::cerr<<"Successive accumulation of data on coarse levels only works with ParMETIS installed."
01013                    <<"  Fell back to accumulation to one domain on coarsest level"<<std::endl;
01014 #endif
01015           
01016         // accumulate to fewer processors
01017         Matrix* redistMat= new Matrix();
01018         ParallelInformation* redistComm=0;
01019         int nodomains = 1;
01020             
01021         repartitionAndDistributeMatrix(mlevel->getmat(), *redistMat, *infoLevel,
01022                                        redistComm, redistributes_.back(), nodomains,criterion);
01023         MatrixArgs args(*redistMat, *redistComm);
01024         int unknowns = redistMat->N();
01025         unknowns = infoLevel->communicator().sum(unknowns);
01026         if(redistComm->communicator().rank()==0 && criterion.debugLevel()>1)
01027           std::cout<<"Level "<<level<<" redistributed has "<<unknowns<<" unknowns, "<<unknowns/redistComm->communicator().size()
01028                      <<" unknowns per proc (procs="<<redistComm->communicator().size()<<")"<<std::endl;
01029         mlevel.addRedistributed(ConstructionTraits<MatrixOperator>::construct(args));
01030         infoLevel.addRedistributed(redistComm);
01031         infoLevel->freeGlobalLookup();
01032       }
01033 
01034         int levels = matrices_.levels();
01035         maxlevels_ = parallelInformation_.finest()->communicator().max(levels);
01036         assert(matrices_.levels()==redistributes_.size());
01037         if(hasCoarsest() && rank==0 && criterion.debugLevel()>1)
01038           std::cout<<"operator complexity: "<<((double)allnonzeros)/finenonzeros<<std::endl;
01039 
01040     }
01041     
01042     template<class M, class IS, class A>
01043     const typename MatrixHierarchy<M,IS,A>::ParallelMatrixHierarchy& 
01044     MatrixHierarchy<M,IS,A>::matrices() const
01045     {
01046       return matrices_;
01047     }
01048     
01049     template<class M, class IS, class A>
01050     const typename MatrixHierarchy<M,IS,A>::ParallelInformationHierarchy& 
01051     MatrixHierarchy<M,IS,A>::parallelInformation() const
01052     {
01053       return parallelInformation_;
01054     }
01055     
01056     template<class M, class IS, class A>
01057     void MatrixHierarchy<M,IS,A>::getCoarsestAggregatesOnFinest(std::vector<std::size_t>& data) const
01058     {
01059       int levels=aggregatesMaps().size();
01060       int maxlevels=parallelInformation_.finest()->communicator().max(levels);
01061       std::size_t size=(*(aggregatesMaps().begin()))->noVertices();
01062       // We need an auxiliary vector for the consecutive prolongation.
01063       std::vector<std::size_t> tmp;
01064       std::vector<std::size_t> *coarse, *fine;
01065 
01066       // make sure the allocated space suffices.
01067       tmp.reserve(size);
01068       data.reserve(size);
01069 
01070       // Correctly assign coarse and fine for the first prolongation such that
01071       // we end up in data in the end.
01072       if(levels%2==0){
01073         coarse=&tmp;
01074         fine=&data;
01075       }else{
01076         coarse=&data;
01077         fine=&tmp;
01078       }
01079 
01080       // Number the unknowns on the coarsest level consecutively for each process.
01081       if(levels==maxlevels){
01082         const AggregatesMap& map = *(*(++aggregatesMaps().rbegin()));
01083         std::size_t m=0;
01084         
01085         for(typename AggregatesMap::const_iterator iter = map.begin(); iter != map.end(); ++iter)
01086           if(*iter< AggregatesMap::ISOLATED)
01087             m=std::max(*iter,m);
01088 
01089         coarse->resize(m+1);
01090         std::size_t i=0;
01091         srand((unsigned)std::clock());
01092         std::set<size_t> used;
01093         for(typename std::vector<std::size_t>::iterator iter=coarse->begin(); iter != coarse->end();
01094             ++iter, ++i)
01095           {
01096             std::pair<std::set<std::size_t>::iterator,bool> ibpair
01097               = used.insert(static_cast<std::size_t>((((double)rand())/(RAND_MAX+1.0)))*coarse->size());
01098             
01099             while(!ibpair.second)
01100               ibpair = used.insert(static_cast<std::size_t>((((double)rand())/(RAND_MAX+1.0))*coarse->size()));
01101             *iter=*(ibpair.first);
01102           }
01103       }
01104 
01105       typename ParallelInformationHierarchy::Iterator pinfo = parallelInformation().coarsest();
01106       --pinfo;
01107             
01108       // Now consecutively project the numbers to the finest level.
01109       for(typename AggregatesMapList::const_reverse_iterator aggregates=++aggregatesMaps().rbegin();
01110           aggregates != aggregatesMaps().rend(); ++aggregates,--levels){
01111         
01112         fine->resize((*aggregates)->noVertices());
01113         fine->assign(fine->size(), 0);
01114         Transfer<typename AggregatesMap::AggregateDescriptor, std::vector<std::size_t>, ParallelInformation>
01115           ::prolongate(*(*aggregates), *coarse, *fine, static_cast<std::size_t>(1), *pinfo);
01116         --pinfo;
01117         std::swap(coarse, fine);
01118       }
01119 
01120       // Assertion to check that we really projected to data on the last step.
01121       assert(coarse==&data);
01122     }
01123     
01124     template<class M, class IS, class A>
01125     const typename MatrixHierarchy<M,IS,A>::AggregatesMapList& 
01126     MatrixHierarchy<M,IS,A>::aggregatesMaps() const
01127     {
01128       return aggregatesMaps_;
01129     }
01130     template<class M, class IS, class A>
01131     const typename MatrixHierarchy<M,IS,A>::RedistributeInfoList& 
01132     MatrixHierarchy<M,IS,A>::redistributeInformation() const 
01133     {
01134       return redistributes_;
01135     }
01136     
01137     template<class M, class IS, class A>
01138     MatrixHierarchy<M,IS,A>::~MatrixHierarchy()
01139     {
01140       typedef typename AggregatesMapList::reverse_iterator AggregatesMapIterator;
01141       typedef typename ParallelMatrixHierarchy::Iterator Iterator;
01142       typedef typename ParallelInformationHierarchy::Iterator InfoIterator;
01143       
01144       AggregatesMapIterator amap = aggregatesMaps_.rbegin();
01145       InfoIterator info = parallelInformation_.coarsest();
01146       for(Iterator level=matrices_.coarsest(), finest=matrices_.finest(); level != finest;  --level, --info, ++amap){
01147         (*amap)->free();
01148         delete *amap;
01149         delete &level->getmat();
01150         if(level.isRedistributed())
01151           delete &(level.getRedistributed().getmat());
01152       }
01153       delete *amap;
01154     }
01155 
01156     template<class M, class IS, class A>
01157     template<class V, class TA>
01158     void MatrixHierarchy<M,IS,A>::coarsenVector(Hierarchy<BlockVector<V,TA> >& hierarchy) const
01159     {
01160       assert(hierarchy.levels()==1);
01161       typedef typename ParallelMatrixHierarchy::ConstIterator Iterator;
01162       typedef typename RedistributeInfoList::const_iterator RIter;
01163       RIter redist = redistributes_.begin();
01164 
01165       Iterator matrix = matrices_.finest(), coarsest = matrices_.coarsest();
01166       int level=0;
01167       if(redist->isSetup())
01168           hierarchy.addRedistributedOnCoarsest(new BlockVector<V,TA>(matrix.getRedistributed().getmat().N()));
01169       Dune::dvverb<<"Level "<<level<<" has "<<matrices_.finest()->getmat().N()<<" unknowns!"<<std::endl;
01170 
01171       while(matrix != coarsest){
01172         ++matrix; ++level; ++redist;
01173         Dune::dvverb<<"Level "<<level<<" has "<<matrix->getmat().N()<<" unknowns!"<<std::endl;
01174           
01175         hierarchy.addCoarser(matrix->getmat().N());
01176         if(redist->isSetup())
01177           hierarchy.addRedistributedOnCoarsest(new BlockVector<V,TA>(matrix.getRedistributed().getmat().N()));
01178         
01179       }
01180           
01181     }
01182     
01183     template<class M, class IS, class A>
01184     template<class S, class TA>
01185     void MatrixHierarchy<M,IS,A>::coarsenSmoother(Hierarchy<S,TA>& smoothers, 
01186                                                     const typename SmootherTraits<S>::Arguments& sargs) const
01187     {
01188       assert(smoothers.levels()==0);
01189       typedef typename ParallelMatrixHierarchy::ConstIterator MatrixIterator;
01190       typedef typename ParallelInformationHierarchy::ConstIterator PinfoIterator;
01191       typedef typename AggregatesMapList::const_iterator AggregatesIterator;
01192       
01193       typename ConstructionTraits<S>::Arguments cargs;
01194       cargs.setArgs(sargs);
01195       PinfoIterator pinfo = parallelInformation_.finest();
01196       AggregatesIterator aggregates = aggregatesMaps_.begin();
01197       int level=0;
01198       for(MatrixIterator matrix = matrices_.finest(), coarsest = matrices_.coarsest(); 
01199           matrix != coarsest; ++matrix, ++pinfo, ++aggregates, ++level){
01200         cargs.setMatrix(matrix->getmat(), **aggregates);
01201         cargs.setComm(*pinfo);
01202         smoothers.addCoarser(cargs);
01203       }
01204       if(maxlevels()>levels()){
01205         // This is not the globally coarsest level and therefore smoothing is needed
01206         cargs.setMatrix(matrices_.coarsest()->getmat(), **aggregates);
01207         cargs.setComm(*pinfo);
01208         smoothers.addCoarser(cargs);
01209         ++level;
01210       }
01211     }
01212     
01213     template<class M, class IS, class A>
01214     template<class F>
01215     void MatrixHierarchy<M,IS,A>::recalculateGalerkin(const F& copyFlags)
01216     {
01217       typedef typename AggregatesMapList::iterator AggregatesMapIterator;
01218       typedef typename ParallelMatrixHierarchy::Iterator Iterator;
01219       typedef typename ParallelInformationHierarchy::Iterator InfoIterator;
01220       
01221       AggregatesMapIterator amap = aggregatesMaps_.begin();
01222       BaseGalerkinProduct productBuilder;
01223       InfoIterator info = parallelInformation_.finest();
01224       
01225       for(Iterator level = matrices_.finest(), coarsest=matrices_.coarsest(); level!=coarsest; ++amap){
01226         const Matrix& fine = level->getmat();
01227         ++level;
01228         ++info;
01229         productBuilder.calculate(fine, *(*amap), const_cast<Matrix&>(level->getmat()), *info, copyFlags);
01230         
01231       }
01232     }
01233 
01234     template<class M, class IS, class A>
01235     std::size_t MatrixHierarchy<M,IS,A>::levels() const
01236     {
01237       return matrices_.levels();
01238     }
01239 
01240     template<class M, class IS, class A>
01241     std::size_t MatrixHierarchy<M,IS,A>::maxlevels() const
01242     {
01243       return maxlevels_;
01244     }
01245 
01246     template<class M, class IS, class A>
01247     bool MatrixHierarchy<M,IS,A>::hasCoarsest() const
01248     {
01249       return levels()==maxlevels() && 
01250         (!matrices_.coarsest().isRedistributed() ||matrices_.coarsest()->getmat().N()>0);
01251     }
01252     
01253     template<class M, class IS, class A>
01254     bool MatrixHierarchy<M,IS,A>::isBuilt() const
01255     {
01256       return built_;
01257     }
01258     
01259     template<class T, class A>
01260     Hierarchy<T,A>::Hierarchy()
01261       : finest_(0), coarsest_(0), nonAllocated_(0), allocator_(), levels_(0)
01262     {}
01263 
01264     template<class T, class A>
01265     Hierarchy<T,A>::Hierarchy(MemberType& first)
01266       : allocator_()
01267     {
01268       finest_ = allocator_.allocate(1,0);
01269       finest_->element_ = &first;
01270       finest_->redistributed_ = 0;
01271       nonAllocated_ = finest_;
01272       coarsest_ = finest_;
01273       coarsest_->coarser_ = coarsest_->finer_ = 0;
01274       levels_ = 1;
01275     }
01276 
01277     template<class T, class A>
01278     Hierarchy<T,A>::~Hierarchy()
01279     {
01280       while(coarsest_){
01281         Element* current = coarsest_;
01282         coarsest_ = coarsest_->finer_;
01283         if(current != nonAllocated_){
01284           if(current->redistributed_)
01285             ConstructionTraits<T>::deconstruct(current->redistributed_);
01286           ConstructionTraits<T>::deconstruct(current->element_);
01287         }
01288         allocator_.deallocate(current, 1);
01289         //coarsest_->coarser_ = 0;
01290       }
01291     }
01292 
01293     template<class T, class A>
01294     std::size_t Hierarchy<T,A>::levels() const
01295     {
01296       return levels_;
01297     }
01298 
01299     template<class T, class A>
01300     void Hierarchy<T,A>::addRedistributedOnCoarsest(T* t)
01301     {
01302       coarsest_->redistributed_ = t;
01303     }
01304     
01305     template<class T, class A>
01306     void Hierarchy<T,A>::addCoarser(Arguments& args)
01307     {
01308       if(!coarsest_){
01309         assert(!finest_);
01310         coarsest_ = allocator_.allocate(1,0);
01311         coarsest_->element_ = ConstructionTraits<MemberType>::construct(args);
01312         finest_ = coarsest_;
01313         coarsest_->finer_ = 0;
01314       }else{
01315         coarsest_->coarser_ = allocator_.allocate(1,0);
01316         coarsest_->coarser_->finer_ = coarsest_;
01317         coarsest_ = coarsest_->coarser_;
01318         coarsest_->element_ = ConstructionTraits<MemberType>::construct(args);
01319       }
01320       coarsest_->redistributed_ = 0;
01321       coarsest_->coarser_=0;
01322       ++levels_;
01323     }
01324    
01325     
01326     template<class T, class A>
01327     void Hierarchy<T,A>::addFiner(Arguments& args)
01328     {
01329       if(!finest_){
01330         assert(!coarsest_);
01331         finest_ = allocator_.allocate(1,0);
01332         finest_->element = ConstructionTraits<T>::construct(args);
01333         coarsest_ = finest_;
01334         coarsest_->coarser_ = coarsest_->finer_ = 0;
01335       }else{
01336         finest_->finer_ = allocator_.allocate(1,0);
01337         finest_->finer_->coarser_ = finest_;
01338         finest_ = finest_->finer_;
01339         finest_->finer = 0;
01340         finest_->element = ConstructionTraits<T>::construct(args);
01341       }
01342       ++levels_;
01343     }
01344 
01345     template<class T, class A>
01346     typename Hierarchy<T,A>::Iterator Hierarchy<T,A>::finest()
01347     {
01348       return Iterator(finest_);
01349     }
01350 
01351     template<class T, class A>
01352     typename Hierarchy<T,A>::Iterator Hierarchy<T,A>::coarsest()
01353     {
01354       return Iterator(coarsest_);
01355     }
01356 
01357     template<class T, class A>
01358     typename Hierarchy<T,A>::ConstIterator Hierarchy<T,A>::finest() const
01359     {
01360       return ConstIterator(finest_);
01361     }
01362 
01363     template<class T, class A>
01364     typename Hierarchy<T,A>::ConstIterator Hierarchy<T,A>::coarsest() const
01365     {
01366       return ConstIterator(coarsest_);
01367     }
01369   }// namespace Amg
01370 }// namespace Dune
01371 
01372 #endif
Generated on Sat Apr 24 11:13:46 2010 for dune-istl by  doxygen 1.6.3