indicessyncer.hh

Go to the documentation of this file.
00001 // $Id$
00002 #ifndef DUNE_INDICESSYNCER_HH
00003 #define DUNE_INDICESSYNCER_HH
00004 
00005 #include"indexset.hh"
00006 #include"remoteindices.hh"
00007 #include<dune/common/stdstreams.hh>
00008 #include<dune/common/tuples.hh>
00009 #include<dune/common/sllist.hh>
00010 #include<cassert>
00011 #include<cmath>
00012 #include<limits>
00013 #include<algorithm>
00014 #include<functional>
00015 #include<map>
00016 
00017 #if HAVE_MPI
00018 namespace Dune
00019 {
00036   template<typename T>
00037   class IndicesSyncer
00038   {
00039   public:
00040 
00042     typedef T ParallelIndexSet;
00043 
00045     typedef typename ParallelIndexSet::IndexPair IndexPair;
00046     
00048     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00049     
00051     typedef typename ParallelIndexSet::LocalIndex::Attribute Attribute;
00052 
00056     typedef Dune::RemoteIndices<ParallelIndexSet> RemoteIndices;
00057     
00067     IndicesSyncer(ParallelIndexSet& indexSet, 
00068                   RemoteIndices& remoteIndices);
00069     
00079     void sync();
00080 
00091     template<typename T1>
00092     void sync(T1& numberer);
00093     
00094   private:    
00095     
00097     ParallelIndexSet& indexSet_;
00098 
00100     RemoteIndices& remoteIndices_;
00101     
00103     char** sendBuffers_;
00104 
00106     char* receiveBuffer_;
00107     
00109     std::size_t* sendBufferSizes_;
00110     
00112     int receiveBufferSize_; // int because of MPI
00113     
00117     struct MessageInformation
00118     {
00119       MessageInformation()
00120         : publish(), pairs()
00121       {}
00123       int publish;
00128       int pairs;
00129     };
00130 
00134     class DefaultNumberer
00135     {
00136     public:
00142       std::size_t operator()(const GlobalIndex& global)
00143       {
00144         return std::numeric_limits<size_t>::max();
00145       }
00146     };
00147     
00149     MPI_Datatype datatype_;
00150 
00152     int rank_;
00153     
00158     typedef SLList<std::pair<GlobalIndex,Attribute>, typename RemoteIndices::Allocator> GlobalIndexList;
00159     
00161     typedef typename GlobalIndexList::ModifyIterator GlobalIndexModifier;
00162     
00166     typedef typename SLList<GlobalIndex, typename RemoteIndices::Allocator>::iterator 
00167     GlobalIndexIterator;
00168 
00170     typedef std::map<int, GlobalIndexList> GlobalIndicesMap;
00171     
00180     GlobalIndicesMap globalMap_;
00181 
00185     typedef SLList<bool, typename RemoteIndices::Allocator> BoolList;
00186 
00190     typedef typename BoolList::iterator BoolIterator;
00191 
00193     typedef typename BoolList::ModifyIterator BoolListModifier;
00194     
00196     typedef std::map<int,BoolList> BoolMap;
00197 
00202     BoolMap oldMap_;
00203     
00205     std::map<int,MessageInformation> infoSend_;
00206     
00208     typedef typename RemoteIndices::RemoteIndexList RemoteIndexList;
00209 
00211     typedef typename RemoteIndexList::ModifyIterator RemoteIndexModifier;
00212     
00214     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00215     
00217     typedef typename RemoteIndexList::iterator RemoteIndexIterator;
00218 
00220     typedef typename RemoteIndexList::const_iterator ConstRemoteIndexIterator;
00221 
00223     typedef Dune::tuple<RemoteIndexModifier,GlobalIndexModifier,BoolListModifier,
00224                   const ConstRemoteIndexIterator> IteratorTuple;
00225 
00233     class Iterators
00234     {
00235       friend class IndicesSyncer<T>;
00236     public:
00246       Iterators(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
00247                 BoolList& booleans);
00248       
00252       Iterators();
00253       
00257       Iterators& operator++();
00258       
00264       void insert(const RemoteIndex& index, 
00265                   const std::pair<GlobalIndex,Attribute>& global);
00266       
00271       RemoteIndex& remoteIndex() const;
00272       
00277       std::pair<GlobalIndex,Attribute>& globalIndexPair() const;
00278 
00279       Attribute& attribute() const;
00280       
00286       bool isOld() const;
00287       
00297       void reset(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
00298                  BoolList& booleans);
00299 
00305       bool isNotAtEnd() const;
00306       
00312       bool isAtEnd() const;
00313       
00314     private:
00324       IteratorTuple iterators_;
00325     };
00326     
00328     typedef std::map<int,Iterators> IteratorsMap;
00329     
00341     IteratorsMap iteratorsMap_;
00342         
00344     void calculateMessageSizes();
00345 
00353     void packAndSend(int destination, char* buffer, std::size_t bufferSize, MPI_Request& req);
00354     
00359     template<typename T1>
00360     void recvAndUnpack(T1& numberer);
00361 
00365     void registerMessageDatatype();
00366     
00370     void insertIntoRemoteIndexList(int process, 
00371                                    const std::pair<GlobalIndex,Attribute>& global, 
00372                                    char attribute);
00373 
00377     void resetIteratorsMap();
00378     
00383     bool checkReset();
00384     
00393     bool checkReset(const Iterators& iterators, RemoteIndexList& rlist, GlobalIndexList& gList,
00394                     BoolList& bList);
00395   };
00396 
00397   template<typename TG, typename TA>
00398   bool operator<(const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
00399                  const std::pair<TG,TA>& i2)
00400   {
00401     return i1.global() < i2.first || 
00402       (i1.global() == i2.first && i1.local().attribute()<i2.second);
00403   }
00404 
00405   template<typename TG, typename TA>
00406   bool operator<(const std::pair<TG,TA>& i1,
00407                   const IndexPair<TG,ParallelLocalIndex<TA> >& i2)
00408   {
00409     return i1.first < i2.global() || 
00410       (i1.first == i2.global() && i1.second<i2.local().attribute());
00411   }
00412 
00413   template<typename TG, typename TA>
00414   bool operator==(const IndexPair<TG,ParallelLocalIndex<TA> >& i1, 
00415                   const std::pair<TG,TA>& i2)
00416   {
00417     return (i1.global() == i2.first && i1.local().attribute()==i2.second);
00418   }
00419   
00420   template<typename TG, typename TA>
00421   bool operator!=(const IndexPair<TG,ParallelLocalIndex<TA> >& i1, 
00422                   const std::pair<TG,TA>& i2)
00423   {
00424     return (i1.global() != i2.first || i1.local().attribute()!=i2.second);
00425   }
00426 
00427   template<typename TG, typename TA>
00428   bool operator==(const std::pair<TG,TA>& i2,
00429                   const IndexPair<TG,ParallelLocalIndex<TA> >& i1)
00430   {
00431     return (i1.global() == i2.first && i1.local().attribute()==i2.second);
00432   }
00433   
00434   template<typename TG, typename TA>
00435   bool operator!=(const std::pair<TG,TA>& i2,
00436                   const IndexPair<TG,ParallelLocalIndex<TA> >& i1)
00437   {
00438     return (i1.global() != i2.first || i1.local().attribute()!=i2.second);
00439   }
00440 
00457   template<typename T, typename A, typename A1>
00458   void storeGlobalIndicesOfRemoteIndices(std::map<int,SLList<std::pair<typename T::GlobalIndex, typename T::LocalIndex::Attribute>,A> >& globalMap,
00459                                        const RemoteIndices<T,A1>& remoteIndices,
00460                                        const T& indexSet)
00461   {
00462     typedef typename RemoteIndices<T,A1>::const_iterator RemoteIterator;
00463     
00464     for(RemoteIterator remote = remoteIndices.begin(), end =remoteIndices.end(); remote != end; ++remote){
00465       typedef typename RemoteIndices<T,A1>::RemoteIndexList RemoteIndexList;
00466       typedef typename RemoteIndexList::const_iterator RemoteIndexIterator;
00467       typedef SLList<std::pair<typename T::GlobalIndex, typename T::LocalIndex::Attribute>,A> GlobalIndexList;
00468       GlobalIndexList& global = globalMap[remote->first];
00469       RemoteIndexList& rList = *(remote->second.first);
00470       
00471       for(RemoteIndexIterator index = rList.begin(), riEnd = rList.end();
00472           index != riEnd; ++index){
00473         global.push_back(std::make_pair(index->localIndexPair().global(),
00474                                         index->localIndexPair().local().attribute()));
00475       }
00476     }
00477   }
00478   
00487   template<typename T, typename A, typename A1>
00488   inline void repairLocalIndexPointers(std::map<int,
00489                                        SLList<std::pair<typename T::GlobalIndex,
00490                                        typename T::LocalIndex::Attribute>,A> >& globalMap,
00491                                        RemoteIndices<T,A1>& remoteIndices,
00492                                        const T& indexSet)
00493   {
00494     typedef typename RemoteIndices<T,A1>::RemoteIndexMap::iterator RemoteIterator;
00495     typedef typename RemoteIndices<T,A1>::RemoteIndexList::iterator RemoteIndexIterator;
00496     typedef typename T::GlobalIndex GlobalIndex;
00497     typedef typename T::LocalIndex::Attribute Attribute;
00498     typedef std::pair<GlobalIndex,Attribute> GlobalIndexPair;
00499     typedef SLList<GlobalIndexPair,A> GlobalIndexPairList;
00500     typedef typename GlobalIndexPairList::iterator GlobalIndexIterator;
00501 
00502     assert(globalMap.size()==static_cast<std::size_t>(remoteIndices.neighbours()));
00503     // Repair pointers to index set in remote indices.
00504     typename std::map<int,GlobalIndexPairList>::iterator global = globalMap.begin();
00505     RemoteIterator end = remoteIndices.remoteIndices_.end();
00506     
00507     for(RemoteIterator remote = remoteIndices.remoteIndices_.begin(); remote != end; ++remote, ++global){
00508       typedef typename T::const_iterator IndexIterator;
00509 
00510       assert(remote->first==global->first);
00511       assert(remote->second.first->size() == global->second.size());
00512       
00513       RemoteIndexIterator riEnd  = remote->second.first->end();
00514       RemoteIndexIterator rIndex = remote->second.first->begin();
00515       GlobalIndexIterator gIndex = global->second.begin();
00516       GlobalIndexIterator gEnd   = global->second.end();
00517       IndexIterator       index  = indexSet.begin();
00518 
00519       assert(rIndex==riEnd || gIndex != global->second.end());
00520       while(rIndex != riEnd){
00521         // Search for the index in the set.
00522         assert(gIndex != global->second.end());
00523         
00524         while(*index < *gIndex)
00525           ++index;
00526         
00527         assert(index != indexSet.end() && *index == *gIndex);
00528         
00529         rIndex->localIndex_ = &(*index);
00530         
00531         ++rIndex;
00532         ++gIndex;
00533       }
00534     }
00535     remoteIndices.sourceSeqNo_ = remoteIndices.source_->seqNo();
00536     remoteIndices.destSeqNo_ = remoteIndices.target_->seqNo();
00537   }
00538   
00539   template<typename T>
00540   IndicesSyncer<T>::IndicesSyncer(ParallelIndexSet& indexSet, 
00541                                       RemoteIndices& remoteIndices)
00542     : indexSet_(indexSet), remoteIndices_(remoteIndices)
00543   {
00544     // index sets must match.
00545     assert(remoteIndices.source_ == remoteIndices.target_);
00546     assert(remoteIndices.source_ == &indexSet);
00547     MPI_Comm_rank(remoteIndices_.communicator(), &rank_);
00548   }
00549     
00550   template<typename T>
00551   IndicesSyncer<T>::Iterators::Iterators(RemoteIndexList& remoteIndices, 
00552                                                GlobalIndexList& globalIndices,
00553                                                BoolList& booleans)
00554     : iterators_(remoteIndices.beginModify(), globalIndices.beginModify(),
00555                     booleans.beginModify(), remoteIndices.end())
00556   { }
00557 
00558   template<typename T>
00559   IndicesSyncer<T>::Iterators::Iterators()
00560     : iterators_()
00561   {}
00562   
00563   template<typename T>
00564   inline typename IndicesSyncer<T>::Iterators& IndicesSyncer<T>::Iterators::operator++()
00565   {
00566       ++(get<0>(iterators_));
00567       ++(get<1>(iterators_));
00568       ++(get<2>(iterators_));
00569       return *this;
00570   }
00571 
00572   template<typename T>
00573   inline void IndicesSyncer<T>::Iterators::insert(const RemoteIndex& index, 
00574                                                  const std::pair<GlobalIndex,Attribute>& global)
00575   {
00576     get<0>(iterators_).insert(index);
00577     get<1>(iterators_).insert(global);
00578     get<2>(iterators_).insert(false);
00579   }
00580   
00581   template<typename T>
00582   inline typename IndicesSyncer<T>::RemoteIndex& 
00583   IndicesSyncer<T>::Iterators::remoteIndex() const
00584   {
00585     return *(get<0>(iterators_));
00586   }
00587   
00588   template<typename T>
00589     inline std::pair<typename IndicesSyncer<T>::GlobalIndex,typename IndicesSyncer<T>::Attribute>&  
00590     IndicesSyncer<T>::Iterators::globalIndexPair() const
00591   {
00592     return *(get<1>(iterators_));
00593   }
00594   
00595   template<typename T>
00596   inline bool IndicesSyncer<T>::Iterators::isOld() const
00597   {
00598     return *(get<2>(iterators_));
00599   }
00600   
00601   template<typename T>
00602   inline void IndicesSyncer<T>::Iterators::reset(RemoteIndexList& remoteIndices, 
00603                                                 GlobalIndexList& globalIndices,
00604                                                 BoolList& booleans)
00605   {
00606     get<0>(iterators_) = remoteIndices.beginModify();
00607     get<1>(iterators_) = globalIndices.beginModify();
00608     get<2>(iterators_) = booleans.beginModify();
00609   }
00610   
00611   template<typename T>
00612   inline bool IndicesSyncer<T>::Iterators::isNotAtEnd() const
00613   {
00614     return get<0>(iterators_)!=get<3>(iterators_);
00615   }
00616   
00617   template<typename T>
00618   inline bool IndicesSyncer<T>::Iterators::isAtEnd() const
00619   {
00620     return get<0>(iterators_)==get<3>(iterators_);
00621   }
00622 
00623   template<typename T>
00624   void IndicesSyncer<T>::registerMessageDatatype()
00625   {
00626     MPI_Datatype type[2] = {MPI_INT, MPI_INT};
00627     int blocklength[2] = {1,1};
00628     MPI_Aint displacement[2];
00629     MPI_Aint base;
00630     
00631     // Compute displacement
00632     MessageInformation message;
00633     
00634     MPI_Address( &(message.publish), displacement);
00635     MPI_Address( &(message.pairs), displacement+1);
00636     
00637     // Make the displacement relative
00638     MPI_Address(&message, &base);
00639     displacement[0] -= base;
00640     displacement[1] -= base;
00641     
00642     MPI_Type_struct( 2, blocklength, displacement, type, &datatype_); 
00643     MPI_Type_commit(&datatype_);
00644   }
00645   
00646   template<typename T>
00647   void IndicesSyncer<T>::calculateMessageSizes()
00648   {    
00649     typedef typename ParallelIndexSet::const_iterator IndexIterator;
00650     typedef CollectiveIterator<T,typename RemoteIndices::Allocator> CollectiveIterator;
00651     
00652     IndexIterator iEnd = indexSet_.end();
00653     CollectiveIterator collIter = remoteIndices_.template iterator<true>();
00654     
00655     for(IndexIterator index = indexSet_.begin(); index != iEnd; ++index){
00656       collIter.advance(index->global(), index->local().attribute());
00657       if(collIter.empty())
00658         break;
00659       int knownRemote=0;
00660 
00661       typedef typename CollectiveIterator::iterator ValidIterator;
00662       ValidIterator end = collIter.end();
00663 
00664       // Count the remote indices we know.
00665       for(ValidIterator valid = collIter.begin(); valid != end; ++valid){
00666         ++knownRemote;
00667       }
00668       
00669       if(knownRemote>0){
00670         Dune::dverb<<rank_<<": publishing "<<knownRemote<<" for index "<<index->global()<< " for processes ";
00671         
00672         // Update MessageInformation
00673         for(ValidIterator valid = collIter.begin(); valid != end; ++valid){
00674           ++(infoSend_[valid.process()].publish);
00675           (infoSend_[valid.process()].pairs) += knownRemote;
00676           Dune::dverb<<valid.process()<<" ";
00677           Dune::dverb<<"(publish="<<infoSend_[valid.process()].publish<<", pairs="<<infoSend_[valid.process()].pairs
00678                <<") ";
00679         }
00680         Dune::dverb<<std::endl;
00681       }
00682     }
00683     
00684     typedef typename std::map<int,MessageInformation>::const_iterator 
00685       MessageIterator;
00686     
00687     const MessageIterator end = infoSend_.end();
00688 
00689     // registerMessageDatatype();
00690     
00691     // Now determine the buffersizes needed for each neighbour using MPI_Pack_size
00692     MessageInformation dummy;
00693     
00694     MessageIterator messageIter= infoSend_.begin();
00695 
00696     typedef typename RemoteIndices::RemoteIndexMap::const_iterator RemoteIterator;
00697     const RemoteIterator rend = remoteIndices_.end();
00698     int neighbour=0;
00699     
00700     for(RemoteIterator remote = remoteIndices_.begin(); remote != rend; ++remote, ++neighbour){
00701       MessageInformation* message;
00702       MessageInformation recv;
00703       
00704       if(messageIter != end && messageIter->first==remote->first){
00705         // We want to send message information to that process
00706         message = const_cast<MessageInformation*>(&(messageIter->second));
00707         ++messageIter;
00708       }else
00709         // We do not want to send information but the other process might.
00710         message = &dummy;
00711 
00712       sendBufferSizes_[neighbour]=0;
00713       int tsize;
00714       // The number of indices published
00715       MPI_Pack_size(1, MPI_INT,remoteIndices_.communicator(), &tsize);
00716       sendBufferSizes_[neighbour] += tsize;
00717 
00718       for(int i=0; i < message->publish; ++i){
00719         // The global index
00720         MPI_Pack_size(1, MPITraits<GlobalIndex>::getType(), remoteIndices_.communicator(), &tsize);
00721         sendBufferSizes_[neighbour] += tsize;
00722         // The attribute in the local index
00723         MPI_Pack_size(1, MPI_CHAR, remoteIndices_.communicator(), &tsize);
00724         sendBufferSizes_[neighbour] += tsize;
00725         // The number of corresponding remote indices
00726         MPI_Pack_size(1, MPI_INT, remoteIndices_.communicator(), &tsize);
00727         sendBufferSizes_[neighbour] += tsize;
00728       }
00729       for(int i=0; i < message->pairs; ++i){
00730         // The process of the remote index
00731         MPI_Pack_size(1, MPI_INT, remoteIndices_.communicator(), &tsize);
00732         sendBufferSizes_[neighbour] += tsize;
00733         // The attribute of the remote index
00734         MPI_Pack_size(1, MPI_CHAR, remoteIndices_.communicator(), &tsize);
00735         sendBufferSizes_[neighbour] += tsize;
00736       }
00737           
00738       Dune::dverb<<rank_<<": Buffer (neighbour="<<remote->first<<") size is "<< sendBufferSizes_[neighbour]<<" for publish="<<message->publish<<" pairs="<<message->pairs<<std::endl;
00739     }
00740 
00741   }
00742   
00743   template<typename T>
00744   inline void IndicesSyncer<T>::sync()
00745   {
00746     DefaultNumberer numberer;
00747     sync(numberer);
00748   }
00749   
00750   template<typename T>
00751   template<typename T1>
00752   void IndicesSyncer<T>::sync(T1& numberer)
00753   {
00754     
00755     // The pointers to the local indices in the remote indices
00756     // will become invalid due to the resorting of the index set.
00757     // Therefore store the corresponding global indices.
00758     // Mark all indices as not added
00759     
00760     typedef typename RemoteIndices::RemoteIndexMap::const_iterator 
00761       RemoteIterator;
00762 
00763     const RemoteIterator end = remoteIndices_.end();
00764     
00765     // Number of neighbours might change during the syncing.
00766     // save the old neighbours
00767     std::size_t noOldNeighbours = remoteIndices_.neighbours();
00768     int* oldNeighbours = new int[noOldNeighbours];
00769     sendBufferSizes_ = new std::size_t[noOldNeighbours];
00770     std::size_t i=0;
00771     
00772     for(RemoteIterator remote = remoteIndices_.begin(); remote != end; ++remote, ++i){
00773       typedef typename RemoteIndices::RemoteIndexList::const_iterator
00774         RemoteIndexIterator;
00775 
00776       oldNeighbours[i]=remote->first;
00777 
00778       // Make sure we only have one remote index list.
00779       assert(remote->second.first==remote->second.second);
00780 
00781       RemoteIndexList& rList = *(remote->second.first);
00782             
00783       // Store the corresponding global indices.
00784       GlobalIndexList& global = globalMap_[remote->first];
00785       BoolList& added = oldMap_[remote->first];
00786       RemoteIndexIterator riEnd = rList.end();
00787       
00788       for(RemoteIndexIterator index = rList.begin();
00789           index != riEnd; ++index){
00790         global.push_back(std::make_pair(index->localIndexPair().global(),
00791                                         index->localIndexPair().local().attribute()));
00792         added.push_back(true);
00793       }
00794 
00795       Iterators iterators(rList, global, added);
00796       iteratorsMap_.insert(std::make_pair(remote->first, iterators));
00797       assert(checkReset(iteratorsMap_[remote->first], rList,global,added));
00798     }
00799     
00800     // Exchange indices with each neighbour
00801     const RemoteIterator rend = remoteIndices_.end();
00802 
00803     calculateMessageSizes();
00804     
00805     // Allocate the buffers
00806     receiveBufferSize_=1;
00807     sendBuffers_ = new char*[noOldNeighbours];
00808     
00809     for(std::size_t i=0; i<noOldNeighbours; ++i){
00810       sendBuffers_[i] = new char[sendBufferSizes_[i]];
00811       receiveBufferSize_ = std::max(receiveBufferSize_, static_cast<int>(sendBufferSizes_[i]));
00812     }
00813     
00814     receiveBuffer_=new char[receiveBufferSize_];
00815     
00816     indexSet_.beginResize();
00817     
00818     Dune::dverb<<rank_<<": Neighbours: ";
00819     
00820     for(i = 0; i<noOldNeighbours; ++i)
00821       Dune::dverb<<oldNeighbours[i]<<" ";
00822     
00823     Dune::dverb<<std::endl;
00824 
00825     MPI_Request* requests = new MPI_Request[noOldNeighbours];
00826     MPI_Status* statuses = new MPI_Status[noOldNeighbours];
00827 
00828     // Pack Message data and start the sends
00829     for(i = 0; i<noOldNeighbours; ++i)
00830         packAndSend(oldNeighbours[i], sendBuffers_[i], sendBufferSizes_[i], requests[i]);
00831 
00832     // Probe for incoming messages, receive and unpack them
00833     for(i = 0; i<noOldNeighbours; ++i)
00834       recvAndUnpack(numberer);
00835 //       }else{
00836 //      recvAndUnpack(oldNeighbours[i], numberer);
00837 //      packAndSend(oldNeighbours[i]);
00838 //       }
00839 //     }
00840 
00841     delete[] receiveBuffer_;
00842     
00843     // Wait for the completion of the sends
00844     // Wait for completion of sends
00845     if(MPI_SUCCESS!=MPI_Waitall(noOldNeighbours, requests, statuses)){
00846       std::cerr<<": MPI_Error occurred while sending message"<<std::endl;
00847       for(i=0;i< noOldNeighbours;i++)
00848         if(MPI_SUCCESS!=statuses[i].MPI_ERROR)
00849           std::cerr<<"Destination "<<statuses[i].MPI_SOURCE<<" error code: "<<statuses[i].MPI_ERROR<<std::endl;
00850     }
00851   
00852     delete[] statuses;
00853     delete[] requests;
00854     
00855     for(std::size_t i=0; i<noOldNeighbours; ++i)
00856       delete[] sendBuffers_[i];
00857 
00858     delete[] sendBuffers_;
00859     delete[] sendBufferSizes_;
00860     
00861     // No need for the iterator tuples any more
00862     iteratorsMap_.clear();
00863     
00864     indexSet_.endResize();
00865     
00866     delete[] oldNeighbours;
00867     
00868     repairLocalIndexPointers(globalMap_, remoteIndices_, indexSet_);
00869     
00870     oldMap_.clear();    
00871     globalMap_.clear();
00872         
00873     // update the sequence number
00874     remoteIndices_.sourceSeqNo_ = remoteIndices_.destSeqNo_ = indexSet_.seqNo();    
00875   }
00876   
00877   template<typename T>
00878   void IndicesSyncer<T>::packAndSend(int destination, char* buffer, std::size_t bufferSize, MPI_Request& request)
00879   {
00880     typedef typename ParallelIndexSet::const_iterator IndexIterator;
00881     typedef typename RemoteIndexList::const_iterator RemoteIndexIterator;
00882     typedef typename GlobalIndexList::const_iterator GlobalIterator;
00883     typedef typename BoolList::const_iterator BoolIterator;
00884         
00885     IndexIterator      iEnd       = indexSet_.end();
00886     int                bpos       = 0;
00887     int                published  = 0;
00888     int                pairs      = 0;
00889     
00890     assert(checkReset());
00891     
00892     // Pack the number of indices we publish
00893     MPI_Pack(&(infoSend_[destination].publish), 1, MPI_INT, buffer, bufferSize, &bpos, 
00894              remoteIndices_.communicator());
00895 
00896     for(IndexIterator index = indexSet_.begin(); index != iEnd; ++index){
00897       // Search for corresponding remote indices in all iterator tuples
00898       typedef typename IteratorsMap::iterator Iterator;
00899       Iterator iteratorsEnd = iteratorsMap_.end();
00900 
00901       // advance all iterators to a position with global index >= index->global()
00902       for(Iterator iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators){
00903         while(iterators->second.isNotAtEnd() && 
00904               iterators->second.globalIndexPair().first < index->global())
00905           ++(iterators->second);
00906         assert(!iterators->second.isNotAtEnd() || iterators->second.globalIndexPair().first >= index->global());
00907       }
00908       
00909       // Add all remote indices positioned at global which were already present before calling sync
00910       // to the message.
00911       // Count how many remote indices we will send
00912       int indices = 0;
00913       bool knownRemote = false; // Is the remote process supposed to know this index?
00914             
00915       for(Iterator iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
00916         {
00917           std::pair<GlobalIndex,Attribute> p;
00918           bool old, valid=iterators->second.isNotAtEnd();
00919           if(valid)
00920             {
00921               p =iterators->second.globalIndexPair();
00922               old=iterators->second.isOld();
00923             }
00924           
00925         if(iterators->second.isNotAtEnd() && iterators->second.isOld() 
00926            && iterators->second.globalIndexPair().first == index->global()){
00927           indices++;
00928           if(destination == iterators->first)
00929             knownRemote = true;
00930         }
00931         }
00932       
00933       if(!knownRemote)
00934         // We do not need to send any indices
00935         continue;
00936       
00937       Dune::dverb<<rank_<<": sending "<<indices<<" for index "<<index->global()<<" to "<<destination<<std::endl;
00938 
00939 
00940       // Pack the global index, the attribute and the number
00941       MPI_Pack(const_cast<GlobalIndex*>(&(index->global())), 1, MPITraits<GlobalIndex>::getType(), buffer, bufferSize, &bpos, 
00942                remoteIndices_.communicator());
00943       
00944       char attr = index->local().attribute();
00945       MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos, 
00946                remoteIndices_.communicator());
00947 
00948       // Pack the number of remote indices we send.
00949       MPI_Pack(&indices, 1, MPI_INT, buffer, bufferSize, &bpos, 
00950                remoteIndices_.communicator());
00951         
00952       // Pack the information about the remote indices
00953       for(Iterator iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
00954         if(iterators->second.isNotAtEnd() && iterators->second.isOld() 
00955            && iterators->second.globalIndexPair().first == index->global()){
00956           int process = iterators->first;
00957 
00958           ++pairs;
00959           assert(pairs <= infoSend_[destination].pairs);
00960           MPI_Pack(&process, 1, MPI_INT, buffer, bufferSize, &bpos, 
00961                    remoteIndices_.communicator());
00962           char attr = iterators->second.remoteIndex().attribute();
00963           
00964           MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos, 
00965                    remoteIndices_.communicator());
00966           --indices;
00967         }
00968       assert(indices==0);
00969       ++published;
00970       Dune::dvverb<<" (publish="<<published<<", pairs="<<pairs<<")"<<std::endl;
00971       assert(published <= infoSend_[destination].publish);
00972     }
00973 
00974     // Make sure we send all expected entries
00975     assert(published == infoSend_[destination].publish);
00976     assert(pairs == infoSend_[destination].pairs);
00977     resetIteratorsMap();
00978 
00979     Dune::dverb << rank_<<": Sending message of "<<bpos<<" bytes to "<<destination<<std::endl;
00980     
00981     MPI_Issend(buffer, bpos, MPI_PACKED, destination, 345, remoteIndices_.communicator(),&request);
00982   }
00983 
00984   template<typename T>
00985     inline void IndicesSyncer<T>::insertIntoRemoteIndexList(int process, 
00986                                                             const std::pair<GlobalIndex,Attribute>& globalPair,
00987                                                             char attribute)
00988   {
00989     Dune::dverb<<"Inserting from "<<process<<" "<<globalPair.first<<", "<<
00990                 globalPair.second<<" "<<attribute<<std::endl;
00991 
00992     resetIteratorsMap();
00993 
00994     // There might be cases where there no remote indices for that process yet
00995     typename IteratorsMap::iterator found = iteratorsMap_.find(process);
00996     
00997     if( found == iteratorsMap_.end() ){
00998       Dune::dverb<<"Discovered new neighbour "<<process<<std::endl;
00999       RemoteIndexList* rlist = new RemoteIndexList();
01000       remoteIndices_.remoteIndices_.insert(std::make_pair(process,std::make_pair(rlist,rlist)));
01001       Iterators iterators = Iterators(*rlist, globalMap_[process], oldMap_[process]);
01002       found = iteratorsMap_.insert(std::make_pair(process, iterators)).first;
01003     }
01004     
01005     Iterators& iterators = found->second;
01006     
01007     // Search for the remote index
01008     while(iterators.isNotAtEnd() && iterators.globalIndexPair() < globalPair){
01009       // Increment all iterators
01010       ++iterators;
01011       
01012     }
01013     
01014     if(iterators.isAtEnd() || iterators.globalIndexPair() != globalPair){
01015       // The entry is not yet known
01016       // Insert in the the list and do not change the first iterator.
01017       iterators.insert(RemoteIndex(Attribute(attribute)),globalPair);
01018       return;
01019     }
01020     
01021     // Global indices match
01022     bool indexIsThere=false;
01023     for(Iterators tmpIterators = iterators; 
01024         !tmpIterators.isAtEnd() && tmpIterators.globalIndexPair() == globalPair;
01025         ++tmpIterators)
01026       //entry already exists with the same attribute
01027       if(tmpIterators.globalIndexPair().second == attribute){
01028         indexIsThere=true;
01029         break;
01030       }
01031     
01032     if(!indexIsThere)
01033       // The entry is not yet known
01034       // Insert in the the list and do not change the first iterator.
01035       iterators.insert(RemoteIndex(Attribute(attribute)),globalPair);
01036   }
01037   
01038   template<typename T>
01039   template<typename T1>
01040   void IndicesSyncer<T>::recvAndUnpack(T1& numberer)
01041   {
01042     typedef typename ParallelIndexSet::const_iterator IndexIterator;
01043     typedef typename RemoteIndexList::iterator RemoteIndexIterator;
01044     typedef typename GlobalIndexList::iterator GlobalIndexIterator;
01045 
01046     IndexIterator    iEnd   = indexSet_.end();
01047     IndexIterator    index  = indexSet_.begin();
01048     int              bpos   = 0;
01049     int              publish;
01050     
01051     assert(checkReset());
01052 
01053     MPI_Status status;
01054 
01055     // We have to determine the message size and source before the receive
01056     MPI_Probe(MPI_ANY_SOURCE, 345, remoteIndices_.communicator(), &status);
01057     
01058     int source=status.MPI_SOURCE;
01059     int count;
01060     MPI_Get_count(&status, MPI_PACKED, &count);
01061 
01062     Dune::dvverb<<rank_<<": Receiving message from "<< source<<" with "<<count<<" bytes"<<std::endl;
01063 
01064     if(count>receiveBufferSize_){
01065       receiveBufferSize_=count;
01066       delete[] receiveBuffer_;
01067       receiveBuffer_ = new char[receiveBufferSize_];
01068     }
01069 
01070     MPI_Recv(receiveBuffer_, count, MPI_PACKED, source, 345, remoteIndices_.communicator(), &status);
01071     
01072     // How many global entries were published?
01073     MPI_Unpack(receiveBuffer_,  count, &bpos, &publish, 1, MPI_INT, remoteIndices_.communicator());
01074     
01075     // Now unpack the remote indices and add them.
01076     while(publish>0){
01077       
01078       // Unpack information about the local index on the source process
01079       GlobalIndex  global;          // global index of the current entry
01080       char sourceAttribute; // Attribute on the source process
01081       int pairs;
01082       
01083       MPI_Unpack(receiveBuffer_,  count, &bpos, &global, 1, MPITraits<GlobalIndex>::getType(), 
01084                  remoteIndices_.communicator());
01085       MPI_Unpack(receiveBuffer_,  count, &bpos, &sourceAttribute, 1, MPI_CHAR, 
01086                  remoteIndices_.communicator());
01087       MPI_Unpack(receiveBuffer_,  count, &bpos, &pairs, 1, MPI_INT, 
01088                  remoteIndices_.communicator());
01089     
01090       // Insert the entry on the remote process to our
01091       // remote index list
01092       SLList<std::pair<int,Attribute> > sourceAttributeList;
01093       sourceAttributeList.push_back(std::make_pair(source,Attribute(sourceAttribute)));
01094       bool foundSelf=false;
01095       Attribute myAttribute=Attribute();
01096       
01097       // Unpack the remote indices
01098       for(;pairs>0;--pairs){
01099         // Unpack the process id that knows the index
01100         int process;
01101         char attribute;
01102         MPI_Unpack(receiveBuffer_,  count, &bpos, &process, 1, MPI_INT, 
01103                  remoteIndices_.communicator());
01104         // Unpack the attribute
01105         MPI_Unpack(receiveBuffer_,  count, &bpos, &attribute, 1, MPI_CHAR, 
01106                    remoteIndices_.communicator());
01107 
01108         if(process==rank_){
01109           foundSelf=true;
01110           myAttribute=Attribute(attribute);
01111           // Now we know the local attribute of the global index
01112           //Only add the index if it is unknown.
01113           // Do we know that global index already?
01114           IndexIterator pos = std::lower_bound(index, iEnd, IndexPair(global));
01115 
01116           if(pos == iEnd || pos->global() != global){
01117             // no entry with this global index
01118             indexSet_.add(global,
01119                           ParallelLocalIndex<Attribute>(numberer(global),
01120                                                         myAttribute, true));
01121             Dune::dvverb << "Adding "<<global<<" "<<myAttribute<<std::endl;
01122             continue;
01123           }
01124           
01125           // because of above the global indices match. Add only if the attribute is different
01126           bool indexIsThere = false;
01127           index=pos;
01128 
01129           for(;pos->global()==global;++pos)
01130             if(pos->local().attribute() == myAttribute){
01131               Dune::dvverb<<"found "<<global<<" "<<myAttribute<<std::endl;
01132               indexIsThere = true;
01133               break;
01134             }
01135           
01136           if(!indexIsThere){
01137             indexSet_.add(global,
01138                           ParallelLocalIndex<Attribute>(numberer(global),
01139                                                         myAttribute, true));  
01140             Dune::dvverb << "Adding "<<global<<" "<<myAttribute<<std::endl;
01141           }
01142           
01143         }else{
01144           sourceAttributeList.push_back(std::make_pair(source,Attribute(attribute)));
01145         }
01146       }
01147       assert(foundSelf);
01148       // Insert remote indices
01149       typedef typename SLList<std::pair<int,Attribute> >::const_iterator Iter;
01150       for(Iter i=sourceAttributeList.begin(), end=sourceAttributeList.end();
01151           i!=end; ++i)
01152         insertIntoRemoteIndexList(i->first, std::make_pair(global, myAttribute),
01153                                   i->second);
01154       --publish;
01155     }
01156     
01157     resetIteratorsMap();
01158   }
01159   
01160   template<typename T>
01161   void IndicesSyncer<T>::resetIteratorsMap(){
01162     
01163     // Reset iterators in all tuples.
01164     typedef typename IteratorsMap::iterator Iterator;
01165     typedef typename RemoteIndices::RemoteIndexMap::iterator 
01166       RemoteIterator;
01167     typedef typename GlobalIndicesMap::iterator GlobalIterator;
01168     typedef typename BoolMap::iterator BoolIterator;
01169     
01170     const RemoteIterator remoteEnd = remoteIndices_.remoteIndices_.end();
01171     Iterator iterators = iteratorsMap_.begin();
01172     GlobalIterator global = globalMap_.begin();
01173     BoolIterator added = oldMap_.begin();
01174     
01175     for(RemoteIterator remote = remoteIndices_.remoteIndices_.begin();
01176         remote != remoteEnd; ++remote, ++global, ++added, ++iterators){
01177       iterators->second.reset(*(remote->second.first), global->second, added->second);
01178     }
01179   }
01180  
01181   template<typename T>
01182   bool IndicesSyncer<T>::checkReset(const Iterators& iterators, RemoteIndexList& rList, GlobalIndexList& gList,
01183                                           BoolList& bList){
01184     
01185     if(get<0>(iterators.iterators_)!=rList.begin())
01186       return false;
01187     if(get<1>(iterators.iterators_)!=gList.begin())
01188       return false;
01189     if(get<2>(iterators.iterators_)!=bList.begin())
01190       return false;
01191     return true;
01192   }
01193   
01194 
01195   template<typename T>
01196   bool IndicesSyncer<T>::checkReset(){
01197     
01198     // Reset iterators in all tuples.
01199     typedef typename IteratorsMap::iterator Iterator;
01200     typedef typename RemoteIndices::RemoteIndexMap::iterator 
01201       RemoteIterator;
01202     typedef typename GlobalIndicesMap::iterator GlobalIterator;
01203     typedef typename BoolMap::iterator BoolIterator;
01204     
01205     const RemoteIterator remoteEnd = remoteIndices_.remoteIndices_.end();
01206     Iterator iterators = iteratorsMap_.begin();
01207     GlobalIterator global = globalMap_.begin();
01208     BoolIterator added = oldMap_.begin();
01209     bool ret = true;
01210     
01211     for(RemoteIterator remote = remoteIndices_.remoteIndices_.begin();
01212         remote != remoteEnd; ++remote, ++global, ++added, ++iterators){
01213       if(!checkReset(iterators->second, *(remote->second.first), global->second,
01214                      added->second))
01215         ret=false;
01216     }
01217     return ret;
01218   } 
01219 }
01220 
01221 #endif
01222 #endif

Generated on Fri Apr 29 2011 with Doxygen (ver 1.7.1) [doxygen-log,error-log].