dune-istl  2.1.1
hierarchy.hh
Go to the documentation of this file.
00001 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
00002 // vi: set et ts=4 sw=2 sts=2:
00003 // $Id: hierarchy.hh 1512 2012-01-03 14:21:12Z mblatt $
00004 #ifndef DUNE_AMGHIERARCHY_HH
00005 #define DUNE_AMGHIERARCHY_HH
00006 
00007 #include<list>
00008 #include<memory>
00009 #include<limits>
00010 #include<algorithm>
00011 #include"aggregates.hh"
00012 #include"graph.hh"
00013 #include"galerkin.hh"
00014 #include"renumberer.hh"
00015 #include"graphcreator.hh"
00016 #include<dune/common/stdstreams.hh>
00017 #include<dune/common/timer.hh>
00018 #include<dune/common/tuples.hh>
00019 #include<dune/common/bigunsignedint.hh>
00020 #include<dune/istl/bvector.hh>
00021 #include<dune/common/parallel/indexset.hh>
00022 #include<dune/istl/matrixutils.hh>
00023 #include<dune/istl/matrixredistribute.hh>
00024 #include<dune/istl/paamg/dependency.hh>
00025 #include<dune/istl/paamg/graph.hh>
00026 #include<dune/istl/paamg/indicescoarsener.hh>
00027 #include<dune/istl/paamg/globalaggregates.hh>
00028 #include<dune/istl/paamg/construction.hh>
00029 #include<dune/istl/paamg/smoother.hh>
00030 #include<dune/istl/paamg/transfer.hh>
00031 
00032 namespace Dune
00033 {
00034   namespace Amg
00035   {
00047     enum{ 
00055       MAX_PROCESSES = 72000};
00056 
00064     template<typename T, typename A=std::allocator<T> >
00065     class Hierarchy
00066     {
00067     public:
00071       typedef T MemberType;
00072 
00073       template<typename T1, typename T2>
00074       class LevelIterator;
00075 
00076     private:
00080       struct Element
00081       {
00082         friend class LevelIterator<Hierarchy<T,A>, T>;
00083         friend class LevelIterator<const Hierarchy<T,A>, const T>;
00084 
00086         Element* coarser_;
00087         
00089         Element* finer_;
00090         
00092         MemberType* element_;
00093         
00095         MemberType* redistributed_;
00096       };
00097     public:
00098 //       enum{
00099 //      /**
00100 //       * @brief If true only the method addCoarser will be usable
00101 //       * otherwise only the method addFiner will be usable.
00102 //       */
00103 //      coarsen = b
00104 //        };
00105       
00109       typedef typename A::template rebind<Element>::other Allocator;
00110 
00111       typedef typename ConstructionTraits<T>::Arguments Arguments;
00112       
00117       Hierarchy(MemberType& first);
00118       
00122       Hierarchy();
00123 
00128       void addCoarser(Arguments& args);
00129 
00130       void addRedistributedOnCoarsest(T* t);
00131       
00136       void addFiner(Arguments& args);
00137 
00144       template<class C, class T1>
00145       class LevelIterator 
00146         : public BidirectionalIteratorFacade<LevelIterator<C,T1>,T1,T1&>
00147       {
00148         friend class LevelIterator<typename remove_const<C>::type,
00149                                    typename remove_const<T1>::type >;
00150         friend class LevelIterator<const typename remove_const<C>::type, 
00151                                    const typename remove_const<T1>::type >;
00152 
00153       public:
00155         LevelIterator()
00156           : element_(0)
00157         {}
00158         
00159         LevelIterator(Element* element)
00160           : element_(element)
00161         {}
00162         
00164         LevelIterator(const LevelIterator<typename remove_const<C>::type,
00165                       typename remove_const<T1>::type>& other)
00166           : element_(other.element_)
00167         {}
00168         
00170         LevelIterator(const LevelIterator<const typename remove_const<C>::type,
00171                       const typename remove_const<T1>::type>& other)
00172           : element_(other.element_)
00173         {}
00174 
00178         bool equals(const LevelIterator<typename remove_const<C>::type, 
00179                     typename remove_const<T1>::type>& other) const
00180         {
00181           return element_ == other.element_;
00182         }
00183         
00187         bool equals(const LevelIterator<const typename remove_const<C>::type, 
00188                     const typename remove_const<T1>::type>& other) const
00189         {
00190           return element_ == other.element_;
00191         }
00192 
00194         T1& dereference() const
00195         {
00196           return *(element_->element_);
00197         }
00198         
00200         void increment()
00201         {
00202           element_ = element_->coarser_;
00203         }
00204         
00206         void decrement()
00207         {
00208           element_ = element_->finer_;
00209         }
00210         
00215         bool isRedistributed() const
00216         {
00217           return element_->redistributed_;
00218         }
00219         
00224         T1& getRedistributed() const
00225         {
00226           assert(element_->redistributed_);
00227           return *element_->redistributed_;
00228         }
00229         void addRedistributed(T1* t)
00230         {
00231           element_->redistributed_ = t;
00232         }
00233 
00234         void deleteRedistributed()
00235         {
00236           element_->redistributed_ = 0;
00237         }
00238         
00239       private:
00240         Element* element_;
00241       };
00242       
00244       typedef LevelIterator<Hierarchy<T,A>,T> Iterator;
00245 
00247       typedef LevelIterator<const Hierarchy<T,A>, const T> ConstIterator;
00248 
00253       Iterator finest();
00254       
00259       Iterator coarsest();
00260       
00261       
00266       ConstIterator finest() const;
00267       
00272       ConstIterator coarsest() const;
00273 
00278       std::size_t levels() const;
00279       
00281       ~Hierarchy();
00282       
00283     private:
00285       Element* finest_;
00287       Element* coarsest_;
00289       Element* nonAllocated_;
00291       Allocator allocator_;
00293       int levels_;
00294     };
00295     
00302     template<class M, class PI, class A=std::allocator<M> >
00303     class MatrixHierarchy
00304     {
00305     public:
00307       typedef M MatrixOperator;
00308       
00310       typedef typename MatrixOperator::matrix_type Matrix;
00311 
00313       typedef PI ParallelInformation;
00314       
00316       typedef A Allocator;
00317 
00319       typedef Dune::Amg::AggregatesMap<typename MatrixGraph<Matrix>::VertexDescriptor> AggregatesMap;
00320 
00322       typedef Dune::Amg::Hierarchy<MatrixOperator,Allocator> ParallelMatrixHierarchy;
00323 
00325       typedef Dune::Amg::Hierarchy<ParallelInformation,Allocator> ParallelInformationHierarchy;
00326 
00328       typedef typename Allocator::template rebind<AggregatesMap*>::other AAllocator;
00329       
00331       typedef std::list<AggregatesMap*,AAllocator> AggregatesMapList;
00332 
00334       typedef RedistributeInformation<ParallelInformation> RedistributeInfoType;
00335       
00337       typedef typename Allocator::template rebind<RedistributeInfoType>::other RILAllocator;
00338       
00340       typedef std::list<RedistributeInfoType,RILAllocator> RedistributeInfoList;
00341 
00347       MatrixHierarchy(const MatrixOperator& fineMatrix,
00348                       const ParallelInformation& pinfo=ParallelInformation());
00349       
00350       
00351       ~MatrixHierarchy();
00352       
00358       template<typename O, typename T>
00359       void build(const T& criterion);
00360 
00368       template<class F>
00369       void recalculateGalerkin(const F& copyFlags);
00370 
00375       template<class V, class TA>
00376       void coarsenVector(Hierarchy<BlockVector<V,TA> >& hierarchy) const;
00377 
00383       template<class S, class TA>
00384       void coarsenSmoother(Hierarchy<S,TA>& smoothers, 
00385                            const typename SmootherTraits<S>::Arguments& args) const;
00386       
00391       std::size_t levels() const;
00392       
00397       std::size_t maxlevels() const;
00398       
00399       bool hasCoarsest() const;
00400       
00405       bool isBuilt() const;
00406 
00411       const ParallelMatrixHierarchy& matrices() const;
00412       
00417       const ParallelInformationHierarchy& parallelInformation() const;
00418       
00423       const AggregatesMapList& aggregatesMaps() const;
00424       
00430       const RedistributeInfoList& redistributeInformation() const;
00431       
00432 
00433       typename MatrixOperator::field_type getProlongationDampingFactor() const
00434       {
00435         return prolongDamp_;
00436       }
00437       
00448       void getCoarsestAggregatesOnFinest(std::vector<std::size_t>& data) const;
00449       
00450     private:
00451       typedef typename ConstructionTraits<MatrixOperator>::Arguments MatrixArgs;
00452       typedef typename ConstructionTraits<ParallelInformation>::Arguments CommunicationArgs;
00454       AggregatesMapList aggregatesMaps_;
00456       RedistributeInfoList redistributes_;
00458       ParallelMatrixHierarchy matrices_;
00460       ParallelInformationHierarchy parallelInformation_;
00461 
00463       bool built_;
00464 
00466       int maxlevels_;
00467       
00468       typename MatrixOperator::field_type prolongDamp_;
00469       
00473       template<class Matrix, bool print>
00474       struct MatrixStats
00475       {
00476         
00480         static void stats(const Matrix& matrix)
00481         {}
00482       };
00483       
00484       template<class Matrix>
00485       struct MatrixStats<Matrix,true>
00486       {
00487         struct calc
00488         {
00489           typedef typename Matrix::size_type size_type;
00490           typedef typename Matrix::row_type matrix_row;
00491           
00492           calc()
00493           {
00494             min=std::numeric_limits<size_type>::max();
00495             max=0;
00496             sum=0;
00497           }
00498           
00499           void operator()(const matrix_row& row)
00500           {
00501             min=std::min(min, row.size());
00502             max=std::max(max, row.size());
00503             sum += row.size();
00504           }
00505           
00506           size_type min;
00507           size_type max;
00508           size_type sum;
00509         };
00513         static void stats(const Matrix& matrix)
00514         {
00515           calc c= for_each(matrix.begin(), matrix.end(), calc());
00516           dinfo<<"Matrix row: min="<<c.min<<" max="<<c.max
00517                <<" average="<<static_cast<double>(c.sum)/matrix.N()
00518                <<std::endl;
00519         }
00520       };
00521     };
00522 
00526     enum AccumulationMode{
00532       noAccu = 0,
00538       atOnceAccu=1,
00542       successiveAccu=2
00543     };
00544           
00545 
00549     template<class T>
00550     class CoarsenCriterion : public T
00551     {
00552     public:
00557       typedef T DependencyCriterion;
00558       
00562       void setMaxLevel(int l)
00563       {
00564         maxLevel_ = l;
00565       }
00569       int maxLevel() const
00570       {
00571         return maxLevel_;
00572       }
00573       
00577       void setCoarsenTarget(int nodes)
00578       {
00579         coarsenTarget_ = nodes;
00580       }
00581       
00585       int coarsenTarget() const
00586       {
00587         return coarsenTarget_;
00588       }
00589       
00595       void setMinCoarsenRate(double rate)
00596       {
00597         minCoarsenRate_ = rate;
00598       }
00599       
00603       double minCoarsenRate() const
00604       {
00605         return minCoarsenRate_;
00606       }
00607 
00611       AccumulationMode accumulate() const
00612       {
00613         return accumulate_;
00614       }
00618       void setAccumulate(AccumulationMode accu)
00619       {
00620         accumulate_=accu;
00621       }
00622 
00623       void setAccumulate(bool accu){
00624         accumulate_=accu?successiveAccu:noAccu;
00625       }
00631        void setProlongationDampingFactor(double d)
00632       {
00633         dampingFactor_ = d;
00634       }
00635 
00641       double getProlongationDampingFactor() const
00642       {
00643         return dampingFactor_;
00644       }
00655       CoarsenCriterion(int maxLevel=100, int coarsenTarget=1000, double minCoarsenRate=1.2,
00656                        double prolongDamp=1.6, AccumulationMode accumulate=successiveAccu)
00657         : T(), maxLevel_(maxLevel), coarsenTarget_(coarsenTarget), minCoarsenRate_(minCoarsenRate),
00658           dampingFactor_(prolongDamp), accumulate_( accumulate)
00659       {}
00660       
00661     private:
00665       int maxLevel_;
00669       int coarsenTarget_;
00673       double minCoarsenRate_;
00677       double dampingFactor_;
00682       AccumulationMode accumulate_;
00683     };
00684     
00685     template<typename M, typename C1>
00686     bool repartitionAndDistributeMatrix(const M& origMatrix, M& newMatrix, 
00687                                         SequentialInformation& origSequentialInformationomm, 
00688                                         SequentialInformation*& newComm, 
00689                                         RedistributeInformation<SequentialInformation>& ri,
00690                                         int nparts, C1& criterion)
00691     {
00692       DUNE_THROW(NotImplemented, "Redistribution does not make sense in sequential code!");
00693     }
00694     
00695     
00696     template<typename M, typename C, typename C1>
00697     bool repartitionAndDistributeMatrix(const M& origMatrix, M& newMatrix, C& origComm, C*& newComm, 
00698                                         RedistributeInformation<C>& ri,
00699                                         int nparts, C1& criterion)
00700   {
00701     Timer time;
00702 #ifdef AMG_REPART_ON_COMM_GRAPH
00703     // Done not repartition the matrix graph, but a graph of the communication scheme.
00704     bool existentOnRedist=Dune::commGraphRepartition(origMatrix, origComm, nparts, newComm, 
00705                                                      ri.getInterface(),
00706                                                      criterion.debugLevel()>1);
00707     
00708 #else
00709     typedef Dune::Amg::MatrixGraph<const M> MatrixGraph;
00710     typedef Dune::Amg::PropertiesGraph<MatrixGraph,
00711       VertexProperties,
00712       EdgeProperties,
00713       IdentityMap,
00714       IdentityMap> PropertiesGraph;
00715     MatrixGraph graph(origMatrix);
00716     PropertiesGraph pgraph(graph);
00717     buildDependency(pgraph, origMatrix, criterion);
00718     
00719 #ifdef DEBUG_REPART
00720     if(origComm.communicator().rank()==0)
00721       std::cout<<"Original matrix"<<std::endl;
00722     origComm.communicator().barrier();
00723   printGlobalSparseMatrix(origMatrix, origComm, std::cout);
00724 #endif
00725     bool existentOnRedist=Dune::graphRepartition(pgraph, origComm, nparts,
00726                                                  newComm, ri.getInterface(),
00727                                                  criterion.debugLevel()>1);
00728 #endif // if else AMG_REPART
00729     
00730     if(origComm.communicator().rank()==0  && criterion.debugLevel()>1)
00731       std::cout<<"Repartitioning took "<<time.elapsed()<<" seconds."<<std::endl;
00732 
00733     ri.setSetup();
00734 
00735 #ifdef DEBUG_REPART
00736     ri.checkInterface(origComm.indexSet(), newComm->indexSet(), origComm.communicator());
00737 #endif
00738 
00739     redistributeMatrix(const_cast<M&>(origMatrix), newMatrix, origComm, *newComm, ri);
00740 
00741 #ifdef DEBUG_REPART
00742     if(origComm.communicator().rank()==0)
00743       std::cout<<"Original matrix"<<std::endl;
00744     origComm.communicator().barrier();
00745     if(newComm->communicator().size()>0)
00746       printGlobalSparseMatrix(newMatrix, *newComm, std::cout);
00747     origComm.communicator().barrier();
00748 #endif
00749 
00750     if(origComm.communicator().rank()==0  && criterion.debugLevel()>1)
00751       std::cout<<"Redistributing matrix took "<<time.elapsed()<<" seconds."<<std::endl;
00752     return existentOnRedist;
00753     
00754   }
00755   
00756   template<typename M>
00757   bool repartitionAndDistributeMatrix(M& origMatrix, M& newMatrix, 
00758                                       SequentialInformation& origComm, 
00759                                       SequentialInformation& newComm, 
00760                                       RedistributeInformation<SequentialInformation>& ri)
00761   {
00762     return true;
00763   }
00764 
00765     template<class M, class IS, class A>
00766     MatrixHierarchy<M,IS,A>::MatrixHierarchy(const MatrixOperator& fineOperator,
00767                                              const ParallelInformation& pinfo)
00768       : matrices_(const_cast<MatrixOperator&>(fineOperator)),
00769         parallelInformation_(const_cast<ParallelInformation&>(pinfo))
00770     {
00771       dune_static_assert((static_cast<int>(MatrixOperator::category) == 
00772                           static_cast<int>(SolverCategory::sequential) ||
00773                           static_cast<int>(MatrixOperator::category) == 
00774                           static_cast<int>(SolverCategory::overlapping) ||
00775                           static_cast<int>(MatrixOperator::category) == 
00776                           static_cast<int>(SolverCategory::nonoverlapping)),
00777                          "MatrixOperator must be of category sequential or overlapping or nonoverlapping");
00778       if (static_cast<int>(MatrixOperator::category) != static_cast<int>(pinfo.getSolverCategory()))
00779         DUNE_THROW(ISTLError, "MatrixOperator and ParallelInformation must belong to the same category!");
00780       
00781     }
00782     
00783     template<class M, class IS, class A>
00784     template<typename O, typename T>
00785     void MatrixHierarchy<M,IS,A>::build(const T& criterion)
00786     {
00787 
00788       prolongDamp_ = criterion.getProlongationDampingFactor();
00789       typedef O OverlapFlags;
00790       typedef typename ParallelMatrixHierarchy::Iterator MatIterator;
00791       typedef typename ParallelInformationHierarchy::Iterator PInfoIterator;
00792 
00793       static const int noints=(Dune::Amg::MAX_PROCESSES/4096>0)?(Dune::Amg::MAX_PROCESSES/4096):1;
00794       
00795       typedef bigunsignedint<sizeof(int)*8*noints> BIGINT;      
00796       GalerkinProduct<ParallelInformation> productBuilder;
00797       MatIterator mlevel = matrices_.finest();
00798       MatrixStats<typename M::matrix_type,MINIMAL_DEBUG_LEVEL<=INFO_DEBUG_LEVEL>::stats(mlevel->getmat());
00799       
00800       PInfoIterator infoLevel = parallelInformation_.finest();
00801       BIGINT finenonzeros=countNonZeros(mlevel->getmat());
00802       finenonzeros = infoLevel->communicator().sum(finenonzeros);
00803       BIGINT allnonzeros = finenonzeros;
00804 
00805       
00806       int procs = infoLevel->communicator().size();
00807       int level = 0;
00808       int rank = 0;
00809 
00810       BIGINT unknowns = mlevel->getmat().N();
00811 
00812       unknowns = infoLevel->communicator().sum(unknowns);
00813       double dunknowns=unknowns.todouble();
00814       infoLevel->buildGlobalLookup(mlevel->getmat().N());
00815       redistributes_.push_back(RedistributeInfoType());
00816 
00817       for(; level < criterion.maxLevel(); ++level, ++mlevel){
00818         assert(matrices_.levels()==redistributes_.size());
00819         rank = infoLevel->communicator().rank();
00820         if(rank==0 && criterion.debugLevel()>1)
00821           std::cout<<"Level "<<level<<" has "<<dunknowns<<" unknowns, "<<dunknowns/infoLevel->communicator().size()
00822                    <<" unknowns per proc (procs="<<infoLevel->communicator().size()<<")"<<std::endl;
00823 
00824         MatrixOperator* matrix=&(*mlevel);
00825         ParallelInformation* info =&(*infoLevel);
00826 
00827         if((
00828 #if HAVE_PARMETIS
00829               criterion.accumulate()==successiveAccu
00830 #else
00831               false
00832 #endif
00833              || (criterion.accumulate()==atOnceAccu
00834                  && dunknowns < 30*infoLevel->communicator().size()))
00835            && infoLevel->communicator().size()>1 && 
00836            dunknowns/infoLevel->communicator().size() <= criterion.coarsenTarget())
00837           {
00838             // accumulate to fewer processors
00839             Matrix* redistMat= new Matrix();
00840             ParallelInformation* redistComm=0;
00841             std::size_t nodomains = (std::size_t)std::ceil(dunknowns/(criterion.minAggregateSize()
00842                                                                   *criterion.coarsenTarget()));
00843             if( nodomains<=criterion.minAggregateSize()/2 || 
00844                 dunknowns <= criterion.coarsenTarget() )
00845               nodomains=1;
00846 
00847             bool existentOnNextLevel = 
00848               repartitionAndDistributeMatrix(mlevel->getmat(), *redistMat, *infoLevel,
00849                                              redistComm, redistributes_.back(), nodomains,
00850                                              criterion);
00851             BIGINT unknowns = redistMat->N();
00852             unknowns = infoLevel->communicator().sum(unknowns);
00853             dunknowns= unknowns.todouble();
00854             if(redistComm->communicator().rank()==0 && criterion.debugLevel()>1)
00855               std::cout<<"Level "<<level<<" (redistributed) has "<<dunknowns<<" unknowns, "<<dunknowns/redistComm->communicator().size()
00856                        <<" unknowns per proc (procs="<<redistComm->communicator().size()<<")"<<std::endl;
00857             MatrixArgs args(*redistMat, *redistComm);
00858             mlevel.addRedistributed(ConstructionTraits<MatrixOperator>::construct(args));
00859             assert(mlevel.isRedistributed());
00860             infoLevel.addRedistributed(redistComm);
00861             infoLevel->freeGlobalLookup();
00862             
00863             if(!existentOnNextLevel)
00864               // We do not hold any data on the redistributed partitioning
00865               break;
00866             
00867             // Work on the redistributed Matrix from now on
00868             matrix = &(mlevel.getRedistributed());
00869             info = &(infoLevel.getRedistributed());
00870             info->buildGlobalLookup(matrix->getmat().N());
00871           }
00872         
00873         rank = info->communicator().rank();
00874         procs = info->communicator().size();
00875         if(dunknowns <= criterion.coarsenTarget())
00876           // No further coarsening needed
00877            break;
00878 
00879         typedef PropertiesGraphCreator<MatrixOperator> GraphCreator;
00880         typedef typename GraphCreator::PropertiesGraph PropertiesGraph;
00881         typedef typename GraphCreator::MatrixGraph MatrixGraph;
00882         typedef typename GraphCreator::GraphTuple  GraphTuple;
00883 
00884         typedef typename PropertiesGraph::VertexDescriptor Vertex;
00885         
00886         std::vector<bool> excluded(matrix->getmat().N(), false);
00887 
00888         GraphTuple graphs = GraphCreator::create(*matrix, excluded, *info, OverlapFlags());
00889         
00890         AggregatesMap* aggregatesMap=new AggregatesMap(get<1>(graphs)->maxVertex()+1);
00891 
00892         aggregatesMaps_.push_back(aggregatesMap);
00893 
00894         Timer watch;
00895         watch.reset();
00896         int noAggregates, isoAggregates, oneAggregates, skippedAggregates;
00897         
00898         tie(noAggregates, isoAggregates, oneAggregates, skippedAggregates) =    
00899           aggregatesMap->buildAggregates(matrix->getmat(), *(get<1>(graphs)), criterion);
00900 
00901         if(rank==0 && criterion.debugLevel()>2)
00902           std::cout<<" Have built "<<noAggregates<<" aggregates totally ("<<isoAggregates<<" isolated aggregates, "<<
00903             oneAggregates<<" aggregates of one vertex,  and skipped "<<
00904             skippedAggregates<<" aggregates)."<<std::endl;
00905 #ifdef TEST_AGGLO
00906         {
00907         // calculate size of local matrix in the distributed direction
00908         int start, end, overlapStart, overlapEnd;
00909         int n = UNKNOWNS/procs; // number of unknowns per process
00910         int bigger = UNKNOWNS%procs; // number of process with n+1 unknows
00911         int procs=info->communicator().rank();
00912 
00913         // Compute owner region
00914         if(rank<bigger){
00915           start = rank*(n+1);
00916           end   = (rank+1)*(n+1);
00917         }else{
00918           start = bigger + rank * n;
00919           end   = bigger + (rank + 1) * n;
00920         }
00921   
00922         // Compute overlap region
00923         if(start>0)
00924           overlapStart = start - 1;
00925         else
00926           overlapStart = start;
00927         
00928         if(end<UNKNOWNS)
00929           overlapEnd = end + 1;
00930         else
00931           overlapEnd = end;
00932   
00933         assert((UNKNOWNS)*(overlapEnd-overlapStart)==aggregatesMap->noVertices());
00934         for(int j=0; j< UNKNOWNS; ++j)
00935           for(int i=0; i < UNKNOWNS; ++i)
00936             {
00937               if(i>=overlapStart && i<overlapEnd)
00938                 {
00939                   int no = (j/2)*((UNKNOWNS)/2)+i/2;
00940                   (*aggregatesMap)[j*(overlapEnd-overlapStart)+i-overlapStart]=no;
00941                 }
00942             }
00943         }
00944 #endif
00945         if(info->communicator().rank()==0)
00946           std::cout<<"aggregating finished."<<std::endl;
00947 
00948         BIGINT gnoAggregates=noAggregates;
00949         gnoAggregates = info->communicator().sum(gnoAggregates);
00950         double dgnoAggregates = gnoAggregates.todouble();
00951 #ifdef TEST_AGGLO
00952         BIGINT gnoAggregates=((UNKNOWNS)/2)*((UNKNOWNS)/2);
00953 #endif
00954 
00955         if(criterion.debugLevel()>2 && rank==0)
00956           std::cout << "Building "<<dgnoAggregates<<" aggregates took "<<watch.elapsed()<<" seconds."<<std::endl;       
00957 
00958         if(dgnoAggregates==0 || dunknowns/dgnoAggregates<criterion.minCoarsenRate())
00959     {
00960             if(rank==0)
00961         {
00962               if(dgnoAggregates>0)
00963                 std::cerr << "Stopped coarsening because of rate breakdown "<<dunknowns<<"/"<<dgnoAggregates
00964                           <<"="<<dunknowns/dgnoAggregates<<"<"
00965                   <<criterion.minCoarsenRate()<<std::endl;
00966               else
00967             std::cerr<< "Could not build any aggregates. Probably no connected nodes."<<std::endl;
00968         }
00969             aggregatesMap->free();
00970             delete aggregatesMap;
00971             aggregatesMaps_.pop_back();
00972 
00973             if(criterion.accumulate() && mlevel.isRedistributed() && info->communicator().size()>1){
00974               // coarse level matrix was already redistributed, but to more than 1 process
00975               // Therefore need to delete the redistribution. Further down it will
00976               // then be redistributed to 1 process
00977               delete &(mlevel.getRedistributed().getmat());
00978               mlevel.deleteRedistributed();
00979               delete &(infoLevel.getRedistributed());
00980               infoLevel.deleteRedistributed();
00981               redistributes_.back().resetSetup();
00982             }
00983                 
00984             break;
00985     }
00986         unknowns =  noAggregates;
00987     dunknowns = dgnoAggregates;
00988 
00989     CommunicationArgs commargs(info->communicator(),info->getSolverCategory());
00990     parallelInformation_.addCoarser(commargs);
00991 
00992         ++infoLevel; // parallel information on coarse level
00993                 
00994         typename PropertyMapTypeSelector<VertexVisitedTag,PropertiesGraph>::Type visitedMap = 
00995           get(VertexVisitedTag(), *(get<1>(graphs)));
00996   
00997         watch.reset();
00998         int aggregates = IndicesCoarsener<ParallelInformation,OverlapFlags>
00999 	  ::coarsen(*info,
01000                     *(get<1>(graphs)),
01001                     visitedMap,
01002                     *aggregatesMap,
01003                     *infoLevel);
01004         
01005         GraphCreator::free(graphs);
01006         
01007         if(criterion.debugLevel()>2){
01008           if(rank==0)
01009             std::cout<<"Coarsening of index sets took "<<watch.elapsed()<<" seconds."<<std::endl;
01010         }
01011         
01012         watch.reset();
01013 
01014         infoLevel->buildGlobalLookup(aggregates);
01015         AggregatesPublisher<Vertex,OverlapFlags,ParallelInformation>::publish(*aggregatesMap,
01016                                                                               *info,
01017                                                                               infoLevel->globalLookup());
01018         
01019                 
01020         if(criterion.debugLevel()>2){
01021           if(rank==0)
01022             std::cout<<"Communicating global aggregate numbers took "<<watch.elapsed()<<" seconds."<<std::endl;
01023         }
01024 
01025         watch.reset();
01026         std::vector<bool>& visited=excluded;
01027         
01028         typedef std::vector<bool>::iterator Iterator;
01029         typedef IteratorPropertyMap<Iterator, IdentityMap> VisitedMap2;
01030         Iterator end = visited.end();
01031         for(Iterator iter= visited.begin(); iter != end; ++iter)
01032           *iter=false;
01033         
01034         VisitedMap2 visitedMap2(visited.begin(), Dune::IdentityMap());
01035 
01036         typename MatrixOperator::matrix_type* coarseMatrix;
01037 
01038         coarseMatrix = productBuilder.build(matrix->getmat(), *(get<0>(graphs)), visitedMap2, 
01039                                             *info, 
01040                                             *aggregatesMap,
01041                                             aggregates,
01042                                             OverlapFlags());
01043         
01044         info->freeGlobalLookup();
01045         
01046         delete get<0>(graphs);
01047         productBuilder.calculate(matrix->getmat(), *aggregatesMap, *coarseMatrix, *infoLevel, OverlapFlags());
01048         
01049         if(criterion.debugLevel()>2){
01050           if(rank==0)
01051             std::cout<<"Calculation of Galerkin product took "<<watch.elapsed()<<" seconds."<<std::endl;
01052         }
01053         
01054         BIGINT nonzeros = countNonZeros(*coarseMatrix);
01055         allnonzeros = allnonzeros + infoLevel->communicator().sum(nonzeros);
01056         MatrixArgs args(*coarseMatrix, *infoLevel);
01057         
01058         matrices_.addCoarser(args);
01059         redistributes_.push_back(RedistributeInfoType());
01060       } // end level loop
01061       
01062 
01063       infoLevel->freeGlobalLookup();
01064       
01065       built_=true;
01066       AggregatesMap* aggregatesMap=new AggregatesMap(0);
01067       aggregatesMaps_.push_back(aggregatesMap);
01068 
01069       if(criterion.debugLevel()>0){
01070         if(level==criterion.maxLevel()){
01071           BIGINT unknowns = mlevel->getmat().N();
01072           unknowns = infoLevel->communicator().sum(unknowns);
01073           double dunknowns = unknowns.todouble();
01074           if(rank==0 && criterion.debugLevel()>1){
01075             std::cout<<"Level "<<level<<" has "<<dunknowns<<" unknowns, "<<dunknowns/infoLevel->communicator().size()
01076                      <<" unknowns per proc (procs="<<infoLevel->communicator().size()<<")"<<std::endl;
01077           }
01078         }
01079       }
01080 
01081       if(criterion.accumulate() && !redistributes_.back().isSetup() && 
01082          infoLevel->communicator().size()>1){ 
01083 #if HAVE_MPI && !HAVE_PARMETIS
01084         if(criterion.accumulate()==successiveAccu &&
01085            infoLevel->communicator().rank()==0)
01086           std::cerr<<"Successive accumulation of data on coarse levels only works with ParMETIS installed."
01087                    <<"  Fell back to accumulation to one domain on coarsest level"<<std::endl;
01088 #endif
01089           
01090         // accumulate to fewer processors
01091         Matrix* redistMat= new Matrix();
01092         ParallelInformation* redistComm=0;
01093         int nodomains = 1;
01094             
01095         repartitionAndDistributeMatrix(mlevel->getmat(), *redistMat, *infoLevel,
01096                                        redistComm, redistributes_.back(), nodomains,criterion);
01097         MatrixArgs args(*redistMat, *redistComm);
01098         BIGINT unknowns = redistMat->N();
01099         unknowns = infoLevel->communicator().sum(unknowns);
01100         
01101         if(redistComm->communicator().rank()==0 && criterion.debugLevel()>1){
01102           double dunknowns= unknowns.todouble();
01103           std::cout<<"Level "<<level<<" redistributed has "<<dunknowns<<" unknowns, "<<dunknowns/redistComm->communicator().size()
01104                      <<" unknowns per proc (procs="<<redistComm->communicator().size()<<")"<<std::endl;
01105         }
01106         mlevel.addRedistributed(ConstructionTraits<MatrixOperator>::construct(args));
01107         infoLevel.addRedistributed(redistComm);
01108         infoLevel->freeGlobalLookup();
01109       }
01110 
01111         int levels = matrices_.levels();
01112         maxlevels_ = parallelInformation_.finest()->communicator().max(levels);
01113         assert(matrices_.levels()==redistributes_.size());
01114         if(hasCoarsest() && rank==0 && criterion.debugLevel()>1)
01115           std::cout<<"operator complexity: "<<allnonzeros.todouble()/finenonzeros.todouble()<<std::endl;
01116 
01117     }
01118     
01119     template<class M, class IS, class A>
01120     const typename MatrixHierarchy<M,IS,A>::ParallelMatrixHierarchy& 
01121     MatrixHierarchy<M,IS,A>::matrices() const
01122     {
01123       return matrices_;
01124     }
01125     
01126     template<class M, class IS, class A>
01127     const typename MatrixHierarchy<M,IS,A>::ParallelInformationHierarchy& 
01128     MatrixHierarchy<M,IS,A>::parallelInformation() const
01129     {
01130       return parallelInformation_;
01131     }
01132     
01133     template<class M, class IS, class A>
01134     void MatrixHierarchy<M,IS,A>::getCoarsestAggregatesOnFinest(std::vector<std::size_t>& data) const
01135     {
01136       int levels=aggregatesMaps().size();
01137       int maxlevels=parallelInformation_.finest()->communicator().max(levels);
01138       std::size_t size=(*(aggregatesMaps().begin()))->noVertices();
01139       // We need an auxiliary vector for the consecutive prolongation.
01140       std::vector<std::size_t> tmp;
01141       std::vector<std::size_t> *coarse, *fine;
01142 
01143       // make sure the allocated space suffices.
01144       tmp.reserve(size);
01145       data.reserve(size);
01146 
01147       // Correctly assign coarse and fine for the first prolongation such that
01148       // we end up in data in the end.
01149       if(levels%2==0){
01150         coarse=&tmp;
01151         fine=&data;
01152       }else{
01153         coarse=&data;
01154         fine=&tmp;
01155       }
01156 
01157       // Number the unknowns on the coarsest level consecutively for each process.
01158       if(levels==maxlevels){
01159         const AggregatesMap& map = *(*(++aggregatesMaps().rbegin()));
01160         std::size_t m=0;
01161         
01162         for(typename AggregatesMap::const_iterator iter = map.begin(); iter != map.end(); ++iter)
01163           if(*iter< AggregatesMap::ISOLATED)
01164             m=std::max(*iter,m);
01165 
01166         coarse->resize(m+1);
01167         std::size_t i=0;
01168         srand((unsigned)std::clock());
01169         std::set<size_t> used;
01170         for(typename std::vector<std::size_t>::iterator iter=coarse->begin(); iter != coarse->end();
01171             ++iter, ++i)
01172           {
01173             std::pair<std::set<std::size_t>::iterator,bool> ibpair
01174               = used.insert(static_cast<std::size_t>((((double)rand())/(RAND_MAX+1.0)))*coarse->size());
01175             
01176             while(!ibpair.second)
01177               ibpair = used.insert(static_cast<std::size_t>((((double)rand())/(RAND_MAX+1.0))*coarse->size()));
01178             *iter=*(ibpair.first);
01179           }
01180       }
01181 
01182       typename ParallelInformationHierarchy::Iterator pinfo = parallelInformation().coarsest();
01183       --pinfo;
01184             
01185       // Now consecutively project the numbers to the finest level.
01186       for(typename AggregatesMapList::const_reverse_iterator aggregates=++aggregatesMaps().rbegin();
01187           aggregates != aggregatesMaps().rend(); ++aggregates,--levels){
01188         
01189         fine->resize((*aggregates)->noVertices());
01190         fine->assign(fine->size(), 0);
01191         Transfer<typename AggregatesMap::AggregateDescriptor, std::vector<std::size_t>, ParallelInformation>
01192           ::prolongate(*(*aggregates), *coarse, *fine, static_cast<std::size_t>(1), *pinfo);
01193         --pinfo;
01194         std::swap(coarse, fine);
01195       }
01196 
01197       // Assertion to check that we really projected to data on the last step.
01198       assert(coarse==&data);
01199     }
01200     
01201     template<class M, class IS, class A>
01202     const typename MatrixHierarchy<M,IS,A>::AggregatesMapList& 
01203     MatrixHierarchy<M,IS,A>::aggregatesMaps() const
01204     {
01205       return aggregatesMaps_;
01206     }
01207     template<class M, class IS, class A>
01208     const typename MatrixHierarchy<M,IS,A>::RedistributeInfoList& 
01209     MatrixHierarchy<M,IS,A>::redistributeInformation() const 
01210     {
01211       return redistributes_;
01212     }
01213     
01214     template<class M, class IS, class A>
01215     MatrixHierarchy<M,IS,A>::~MatrixHierarchy()
01216     {
01217       typedef typename AggregatesMapList::reverse_iterator AggregatesMapIterator;
01218       typedef typename ParallelMatrixHierarchy::Iterator Iterator;
01219       typedef typename ParallelInformationHierarchy::Iterator InfoIterator;
01220       
01221       AggregatesMapIterator amap = aggregatesMaps_.rbegin();
01222       InfoIterator info = parallelInformation_.coarsest();
01223       for(Iterator level=matrices_.coarsest(), finest=matrices_.finest(); level != finest;  --level, --info, ++amap){
01224         (*amap)->free();
01225         delete *amap;
01226         delete &level->getmat();
01227         if(level.isRedistributed())
01228           delete &(level.getRedistributed().getmat());
01229       }
01230       delete *amap;
01231     }
01232 
01233     template<class M, class IS, class A>
01234     template<class V, class TA>
01235     void MatrixHierarchy<M,IS,A>::coarsenVector(Hierarchy<BlockVector<V,TA> >& hierarchy) const
01236     {
01237       assert(hierarchy.levels()==1);
01238       typedef typename ParallelMatrixHierarchy::ConstIterator Iterator;
01239       typedef typename RedistributeInfoList::const_iterator RIter;
01240       RIter redist = redistributes_.begin();
01241 
01242       Iterator matrix = matrices_.finest(), coarsest = matrices_.coarsest();
01243       int level=0;
01244       if(redist->isSetup())
01245           hierarchy.addRedistributedOnCoarsest(new BlockVector<V,TA>(matrix.getRedistributed().getmat().N()));
01246       Dune::dvverb<<"Level "<<level<<" has "<<matrices_.finest()->getmat().N()<<" unknowns!"<<std::endl;
01247 
01248       while(matrix != coarsest){
01249         ++matrix; ++level; ++redist;
01250         Dune::dvverb<<"Level "<<level<<" has "<<matrix->getmat().N()<<" unknowns!"<<std::endl;
01251           
01252         hierarchy.addCoarser(matrix->getmat().N());
01253         if(redist->isSetup())
01254           hierarchy.addRedistributedOnCoarsest(new BlockVector<V,TA>(matrix.getRedistributed().getmat().N()));
01255         
01256       }
01257           
01258     }
01259     
01260     template<class M, class IS, class A>
01261     template<class S, class TA>
01262     void MatrixHierarchy<M,IS,A>::coarsenSmoother(Hierarchy<S,TA>& smoothers, 
01263                                                     const typename SmootherTraits<S>::Arguments& sargs) const
01264     {
01265       assert(smoothers.levels()==0);
01266       typedef typename ParallelMatrixHierarchy::ConstIterator MatrixIterator;
01267       typedef typename ParallelInformationHierarchy::ConstIterator PinfoIterator;
01268       typedef typename AggregatesMapList::const_iterator AggregatesIterator;
01269       
01270       typename ConstructionTraits<S>::Arguments cargs;
01271       cargs.setArgs(sargs);
01272       PinfoIterator pinfo = parallelInformation_.finest();
01273       AggregatesIterator aggregates = aggregatesMaps_.begin();
01274       int level=0;
01275       for(MatrixIterator matrix = matrices_.finest(), coarsest = matrices_.coarsest(); 
01276           matrix != coarsest; ++matrix, ++pinfo, ++aggregates, ++level){
01277         cargs.setMatrix(matrix->getmat(), **aggregates);
01278         cargs.setComm(*pinfo);
01279         smoothers.addCoarser(cargs);
01280       }
01281       if(maxlevels()>levels()){
01282         // This is not the globally coarsest level and therefore smoothing is needed
01283         cargs.setMatrix(matrices_.coarsest()->getmat(), **aggregates);
01284         cargs.setComm(*pinfo);
01285         smoothers.addCoarser(cargs);
01286         ++level;
01287       }
01288     }
01289     
01290     template<class M, class IS, class A>
01291     template<class F>
01292     void MatrixHierarchy<M,IS,A>::recalculateGalerkin(const F& copyFlags)
01293     {
01294       typedef typename AggregatesMapList::iterator AggregatesMapIterator;
01295       typedef typename ParallelMatrixHierarchy::Iterator Iterator;
01296       typedef typename ParallelInformationHierarchy::Iterator InfoIterator;
01297       
01298       AggregatesMapIterator amap = aggregatesMaps_.begin();
01299       BaseGalerkinProduct productBuilder;
01300       InfoIterator info = parallelInformation_.finest();
01301       typename RedistributeInfoList::iterator riIter = redistributes_.begin();
01302       Iterator level = matrices_.finest(), coarsest=matrices_.coarsest();
01303       if(level.isRedistributed()){
01304         info->buildGlobalLookup(info->indexSet().size());
01305         redistributeMatrixEntries(const_cast<Matrix&>(level->getmat()), 
01306                                   const_cast<Matrix&>(level.getRedistributed().getmat()),
01307                                   *info,info.getRedistributed(), *riIter);
01308         info->freeGlobalLookup();
01309       }
01310       
01311       for(; level!=coarsest; ++amap){
01312         const Matrix& fine = (level.isRedistributed()?level.getRedistributed():*level).getmat();
01313         ++level;
01314         ++info;
01315         ++riIter;
01316         productBuilder.calculate(fine, *(*amap), const_cast<Matrix&>(level->getmat()), *info, copyFlags);
01317         if(level.isRedistributed()){
01318         info->buildGlobalLookup(info->indexSet().size());
01319           redistributeMatrixEntries(const_cast<Matrix&>(level->getmat()), 
01320                                     const_cast<Matrix&>(level.getRedistributed().getmat()), *info,
01321                                     info.getRedistributed(), *riIter);
01322           info->freeGlobalLookup();
01323         }
01324       }
01325     }
01326 
01327     template<class M, class IS, class A>
01328     std::size_t MatrixHierarchy<M,IS,A>::levels() const
01329     {
01330       return matrices_.levels();
01331     }
01332 
01333     template<class M, class IS, class A>
01334     std::size_t MatrixHierarchy<M,IS,A>::maxlevels() const
01335     {
01336       return maxlevels_;
01337     }
01338 
01339     template<class M, class IS, class A>
01340     bool MatrixHierarchy<M,IS,A>::hasCoarsest() const
01341     {
01342       return levels()==maxlevels() && 
01343         (!matrices_.coarsest().isRedistributed() ||matrices_.coarsest()->getmat().N()>0);
01344     }
01345     
01346     template<class M, class IS, class A>
01347     bool MatrixHierarchy<M,IS,A>::isBuilt() const
01348     {
01349       return built_;
01350     }
01351     
01352     template<class T, class A>
01353     Hierarchy<T,A>::Hierarchy()
01354       : finest_(0), coarsest_(0), nonAllocated_(0), allocator_(), levels_(0)
01355     {}
01356 
01357     template<class T, class A>
01358     Hierarchy<T,A>::Hierarchy(MemberType& first)
01359       : allocator_()
01360     {
01361       finest_ = allocator_.allocate(1,0);
01362       finest_->element_ = &first;
01363       finest_->redistributed_ = 0;
01364       nonAllocated_ = finest_;
01365       coarsest_ = finest_;
01366       coarsest_->coarser_ = coarsest_->finer_ = 0;
01367       levels_ = 1;
01368     }
01369 
01370     template<class T, class A>
01371     Hierarchy<T,A>::~Hierarchy()
01372     {
01373       while(coarsest_){
01374         Element* current = coarsest_;
01375         coarsest_ = coarsest_->finer_;
01376         if(current != nonAllocated_){
01377           if(current->redistributed_)
01378             ConstructionTraits<T>::deconstruct(current->redistributed_);
01379           ConstructionTraits<T>::deconstruct(current->element_);
01380         }
01381         allocator_.deallocate(current, 1);
01382         //coarsest_->coarser_ = 0;
01383       }
01384     }
01385 
01386     template<class T, class A>
01387     std::size_t Hierarchy<T,A>::levels() const
01388     {
01389       return levels_;
01390     }
01391 
01392     template<class T, class A>
01393     void Hierarchy<T,A>::addRedistributedOnCoarsest(T* t)
01394     {
01395       coarsest_->redistributed_ = t;
01396     }
01397     
01398     template<class T, class A>
01399     void Hierarchy<T,A>::addCoarser(Arguments& args)
01400     {
01401       if(!coarsest_){
01402         assert(!finest_);
01403         coarsest_ = allocator_.allocate(1,0);
01404         coarsest_->element_ = ConstructionTraits<MemberType>::construct(args);
01405         finest_ = coarsest_;
01406         coarsest_->finer_ = 0;
01407       }else{
01408         coarsest_->coarser_ = allocator_.allocate(1,0);
01409         coarsest_->coarser_->finer_ = coarsest_;
01410         coarsest_ = coarsest_->coarser_;
01411         coarsest_->element_ = ConstructionTraits<MemberType>::construct(args);
01412       }
01413       coarsest_->redistributed_ = 0;
01414       coarsest_->coarser_=0;
01415       ++levels_;
01416     }
01417    
01418     
01419     template<class T, class A>
01420     void Hierarchy<T,A>::addFiner(Arguments& args)
01421     {
01422       if(!finest_){
01423         assert(!coarsest_);
01424         finest_ = allocator_.allocate(1,0);
01425         finest_->element = ConstructionTraits<T>::construct(args);
01426         coarsest_ = finest_;
01427         coarsest_->coarser_ = coarsest_->finer_ = 0;
01428       }else{
01429         finest_->finer_ = allocator_.allocate(1,0);
01430         finest_->finer_->coarser_ = finest_;
01431         finest_ = finest_->finer_;
01432         finest_->finer = 0;
01433         finest_->element = ConstructionTraits<T>::construct(args);
01434       }
01435       ++levels_;
01436     }
01437 
01438     template<class T, class A>
01439     typename Hierarchy<T,A>::Iterator Hierarchy<T,A>::finest()
01440     {
01441       return Iterator(finest_);
01442     }
01443 
01444     template<class T, class A>
01445     typename Hierarchy<T,A>::Iterator Hierarchy<T,A>::coarsest()
01446     {
01447       return Iterator(coarsest_);
01448     }
01449 
01450     template<class T, class A>
01451     typename Hierarchy<T,A>::ConstIterator Hierarchy<T,A>::finest() const
01452     {
01453       return ConstIterator(finest_);
01454     }
01455 
01456     template<class T, class A>
01457     typename Hierarchy<T,A>::ConstIterator Hierarchy<T,A>::coarsest() const
01458     {
01459       return ConstIterator(coarsest_);
01460     }
01462   }// namespace Amg
01463 }// namespace Dune
01464 
01465 #endif