indicessyncer.hh

Go to the documentation of this file.
00001 // $Id: indicessyncer.hh 1075 2009-09-18 11:13:20Z mblatt $
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 
00016 #if HAVE_MPI
00017 namespace Dune
00018 {
00035   template<typename T>
00036   class IndicesSyncer
00037   {
00038   public:
00039 
00041     typedef T ParallelIndexSet;
00042 
00044     typedef typename ParallelIndexSet::IndexPair IndexPair;
00045     
00047     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00048     
00050     typedef typename ParallelIndexSet::LocalIndex::Attribute Attribute;
00051 
00055     typedef Dune::RemoteIndices<ParallelIndexSet> RemoteIndices;
00056     
00066     IndicesSyncer(ParallelIndexSet& indexSet, 
00067                   RemoteIndices& remoteIndices);
00068     
00077     void sync();
00078 
00088     template<typename T1>
00089     void sync(T1& numberer);
00090     
00091   private:    
00092     
00094     ParallelIndexSet& indexSet_;
00095 
00097     RemoteIndices& remoteIndices_;
00098     
00100     char** sendBuffers_;
00101 
00103     char* receiveBuffer_;
00104     
00106     std::size_t* sendBufferSizes_;
00107     
00109     int receiveBufferSize_; // int because of MPI
00110     
00114     struct MessageInformation
00115     {
00116       MessageInformation()
00117         : publish(), pairs()
00118       {}
00120       int publish;
00125       int pairs;
00126     };
00127 
00131     class DefaultNumberer
00132     {
00133     public:
00139       std::size_t operator()(const GlobalIndex& global)
00140       {
00141         return std::numeric_limits<size_t>::max();
00142       }
00143     };
00144     
00146     MPI_Datatype datatype_;
00147 
00149     int rank_;
00150     
00155     typedef SLList<GlobalIndex, typename RemoteIndices::Allocator> GlobalIndexList;
00156     
00158     typedef typename GlobalIndexList::ModifyIterator GlobalIndexModifier;
00159     
00163     typedef typename SLList<GlobalIndex, typename RemoteIndices::Allocator>::iterator 
00164     GlobalIndexIterator;
00165 
00167     typedef std::map<int, GlobalIndexList> GlobalIndicesMap;
00168     
00177     GlobalIndicesMap globalMap_;
00178 
00182     typedef SLList<bool, typename RemoteIndices::Allocator> BoolList;
00183 
00187     typedef typename BoolList::iterator BoolIterator;
00188 
00190     typedef typename BoolList::ModifyIterator BoolListModifier;
00191     
00193     typedef std::map<int,BoolList> BoolMap;
00194 
00199     BoolMap oldMap_;
00200     
00202     std::map<int,MessageInformation> infoSend_;
00203     
00205     typedef typename RemoteIndices::RemoteIndexList RemoteIndexList;
00206 
00208     typedef typename RemoteIndexList::ModifyIterator RemoteIndexModifier;
00209     
00211     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00212     
00214     typedef typename RemoteIndexList::iterator RemoteIndexIterator;
00215 
00217     typedef typename RemoteIndexList::const_iterator ConstRemoteIndexIterator;
00218 
00220     typedef Dune::tuple<RemoteIndexModifier,GlobalIndexModifier,BoolListModifier,
00221                   const ConstRemoteIndexIterator> IteratorTuple;
00222 
00230     class Iterators
00231     {
00232       friend class IndicesSyncer<T>;
00233     public:
00243       Iterators(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
00244                 BoolList& booleans);
00245       
00249       Iterators();
00250       
00254       Iterators& operator++();
00255       
00261       void insert(const RemoteIndex& index, const GlobalIndex& global);
00262       
00267       RemoteIndex& remoteIndex() const;
00268       
00273       GlobalIndex& globalIndex() const;
00274       
00280       bool isOld() const;
00281       
00291       void reset(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
00292                  BoolList& booleans);
00293 
00299       bool isNotAtEnd() const;
00300       
00306       bool isAtEnd() const;
00307       
00308     private:
00318       IteratorTuple iterators_;
00319     };
00320     
00322     typedef std::map<int,Iterators> IteratorsMap;
00323     
00335     IteratorsMap iteratorsMap_;
00336         
00338     void calculateMessageSizes();
00339 
00347     void packAndSend(int destination, char* buffer, std::size_t bufferSize, MPI_Request& req);
00348     
00353     template<typename T1>
00354     void recvAndUnpack(T1& numberer);
00355 
00359     void registerMessageDatatype();
00360     
00364     void insertIntoRemoteIndexList(int process, const GlobalIndex& global, char attribute);
00365 
00369     void resetIteratorsMap();
00370     
00375     bool checkReset();
00376     
00385     bool checkReset(const Iterators& iterators, RemoteIndexList& rlist, GlobalIndexList& gList,
00386                     BoolList& bList);
00387   };
00388   
00405   template<typename T, typename A, typename A1>
00406   void storeGlobalIndicesOfRemoteIndices(std::map<int,SLList<typename T::GlobalIndex,A> >& globalMap,
00407                                        const RemoteIndices<T,A1>& remoteIndices,
00408                                        const T& indexSet)
00409   {
00410     typedef typename RemoteIndices<T,A1>::const_iterator RemoteIterator;
00411     
00412     for(RemoteIterator remote = remoteIndices.begin(), end =remoteIndices.end(); remote != end; ++remote){
00413       typedef typename RemoteIndices<T,A1>::RemoteIndexList RemoteIndexList;
00414       typedef typename RemoteIndexList::const_iterator RemoteIndexIterator;
00415       typedef SLList<typename T::GlobalIndex,A> GlobalIndexList;
00416       GlobalIndexList& global = globalMap[remote->first];
00417       RemoteIndexList& rList = *(remote->second.first);
00418       
00419       for(RemoteIndexIterator index = rList.begin(), riEnd = rList.end();
00420           index != riEnd; ++index){
00421         global.push_back(index->localIndexPair().global());
00422       }
00423     }
00424   }
00425   
00434   template<typename T, typename A, typename A1>
00435   inline void repairLocalIndexPointers(std::map<int,SLList<typename T::GlobalIndex,A> >& globalMap,
00436                                        RemoteIndices<T,A1>& remoteIndices,
00437                                        const T& indexSet)
00438   {
00439     typedef typename RemoteIndices<T>::RemoteIndexMap::iterator RemoteIterator;
00440     typedef typename RemoteIndices<T>::RemoteIndexList::iterator RemoteIndexIterator;
00441     typedef typename T::GlobalIndex GlobalIndex;
00442     typedef typename SLList<GlobalIndex,A>::iterator GlobalIndexIterator;
00443 
00444     assert(globalMap.size()==remoteIndices.remoteIndices_.size());
00445     // Repair pointers to index set in remote indices.
00446     typename std::map<int,SLList<GlobalIndex,A> >::iterator global = globalMap.begin();
00447     RemoteIterator end = remoteIndices.remoteIndices_.end();
00448     
00449     for(RemoteIterator remote = remoteIndices.remoteIndices_.begin(); remote != end; ++remote, ++global){
00450       typedef typename T::const_iterator IndexIterator;
00451 
00452       assert(remote->first==global->first);
00453       assert(remote->second.first->size() == global->second.size());
00454       
00455       RemoteIndexIterator riEnd  = remote->second.first->end();
00456       RemoteIndexIterator rIndex = remote->second.first->begin();
00457       GlobalIndexIterator gIndex = global->second.begin();
00458       GlobalIndexIterator gEnd   = global->second.end();
00459       IndexIterator       index  = indexSet.begin();
00460 
00461       assert(rIndex==riEnd || gIndex != global->second.end());
00462       while(rIndex != riEnd){
00463         // Search for the index in the set.
00464         assert(gIndex != global->second.end());
00465         
00466         while(index->global() < *gIndex)
00467           ++index;
00468         
00469         assert(index != indexSet.end() && index->global() == *gIndex);
00470         
00471         rIndex->localIndex_ = &(*index);
00472         
00473         ++rIndex;
00474         ++gIndex;
00475       }
00476     }
00477     remoteIndices.sourceSeqNo_ = remoteIndices.source_->seqNo();
00478     remoteIndices.destSeqNo_ = remoteIndices.target_->seqNo();
00479   }
00480   
00481   template<typename T>
00482   IndicesSyncer<T>::IndicesSyncer(ParallelIndexSet& indexSet, 
00483                                       RemoteIndices& remoteIndices)
00484     : indexSet_(indexSet), remoteIndices_(remoteIndices)
00485   {
00486     // index sets must match.
00487     assert(remoteIndices.source_ == remoteIndices.target_);
00488     assert(remoteIndices.source_ == &indexSet);
00489     MPI_Comm_rank(remoteIndices_.communicator(), &rank_);
00490   }
00491     
00492   template<typename T>
00493   IndicesSyncer<T>::Iterators::Iterators(RemoteIndexList& remoteIndices, 
00494                                                GlobalIndexList& globalIndices,
00495                                                BoolList& booleans)
00496     : iterators_(remoteIndices.beginModify(), globalIndices.beginModify(),
00497                     booleans.beginModify(), remoteIndices.end())
00498   { }
00499 
00500   template<typename T>
00501   IndicesSyncer<T>::Iterators::Iterators()
00502     : iterators_()
00503   {}
00504   
00505   template<typename T>
00506   inline typename IndicesSyncer<T>::Iterators& IndicesSyncer<T>::Iterators::operator++()
00507   {
00508       ++(get<0>(iterators_));
00509       ++(get<1>(iterators_));
00510       ++(get<2>(iterators_));
00511       return *this;
00512   }
00513 
00514   template<typename T>
00515   inline void IndicesSyncer<T>::Iterators::insert(const RemoteIndex& index, 
00516                                                  const GlobalIndex& global)
00517   {
00518     get<0>(iterators_).insert(index);
00519     get<1>(iterators_).insert(global);
00520     get<2>(iterators_).insert(false);
00521   }
00522   
00523   template<typename T>
00524   inline typename IndicesSyncer<T>::RemoteIndex& 
00525   IndicesSyncer<T>::Iterators::remoteIndex() const
00526   {
00527     return *(get<0>(iterators_));
00528   }
00529   
00530   template<typename T>
00531   inline typename IndicesSyncer<T>::GlobalIndex&  IndicesSyncer<T>::Iterators::globalIndex() const
00532   {
00533     return *(get<1>(iterators_));
00534   }
00535   
00536   template<typename T>
00537   inline bool IndicesSyncer<T>::Iterators::isOld() const
00538   {
00539     return *(get<2>(iterators_));
00540   }
00541   
00542   template<typename T>
00543   inline void IndicesSyncer<T>::Iterators::reset(RemoteIndexList& remoteIndices, 
00544                                                 GlobalIndexList& globalIndices,
00545                                                 BoolList& booleans)
00546   {
00547     get<0>(iterators_) = remoteIndices.beginModify();
00548     get<1>(iterators_) = globalIndices.beginModify();
00549     get<2>(iterators_) = booleans.beginModify();
00550   }
00551   
00552   template<typename T>
00553   inline bool IndicesSyncer<T>::Iterators::isNotAtEnd() const
00554   {
00555     return get<0>(iterators_)!=get<3>(iterators_);
00556   }
00557   
00558   template<typename T>
00559   inline bool IndicesSyncer<T>::Iterators::isAtEnd() const
00560   {
00561     return get<0>(iterators_)==get<3>(iterators_);
00562   }
00563 
00564   template<typename T>
00565   void IndicesSyncer<T>::registerMessageDatatype()
00566   {
00567     MPI_Datatype type[2] = {MPI_INT, MPI_INT};
00568     int blocklength[2] = {1,1};
00569     MPI_Aint displacement[2];
00570     MPI_Aint base;
00571     
00572     // Compute displacement
00573     MessageInformation message;
00574     
00575     MPI_Address( &(message.publish), displacement);
00576     MPI_Address( &(message.pairs), displacement+1);
00577     
00578     // Make the displacement relative
00579     MPI_Address(&message, &base);
00580     displacement[0] -= base;
00581     displacement[1] -= base;
00582     
00583     MPI_Type_struct( 2, blocklength, displacement, type, &datatype_); 
00584     MPI_Type_commit(&datatype_);
00585   }
00586   
00587   template<typename T>
00588   void IndicesSyncer<T>::calculateMessageSizes()
00589   {    
00590     typedef typename ParallelIndexSet::const_iterator IndexIterator;
00591     typedef CollectiveIterator<T,typename RemoteIndices::Allocator> CollectiveIterator;
00592     
00593     IndexIterator iEnd = indexSet_.end();
00594     CollectiveIterator collIter = remoteIndices_.template iterator<true>();
00595     
00596     for(IndexIterator index = indexSet_.begin(); index != iEnd; ++index){
00597       collIter.advance(index->global());
00598       if(collIter.empty())
00599         break;
00600       int knownRemote=0;
00601 
00602       typedef typename CollectiveIterator::iterator ValidIterator;
00603       ValidIterator end = collIter.end();
00604 
00605       // Count the remote indices we know.
00606       for(ValidIterator valid = collIter.begin(); valid != end; ++valid)
00607         ++knownRemote;
00608       
00609       if(knownRemote>0){
00610         Dune::dverb<<rank_<<": publishing "<<knownRemote<<" for index "<<index->global()<< " for processes ";
00611         
00612         // Update MessageInformation
00613         for(ValidIterator valid = collIter.begin(); valid != end; ++valid){
00614           ++(infoSend_[valid.process()].publish);
00615           (infoSend_[valid.process()].pairs) += knownRemote;
00616           Dune::dverb<<valid.process()<<" ";
00617           Dune::dverb<<"(publish="<<infoSend_[valid.process()].publish<<", pairs="<<infoSend_[valid.process()].pairs
00618                <<") ";
00619         }
00620         Dune::dverb<<std::endl;
00621       }
00622     }
00623     
00624     typedef typename std::map<int,MessageInformation>::const_iterator 
00625       MessageIterator;
00626     
00627     const MessageIterator end = infoSend_.end();
00628 
00629     // registerMessageDatatype();
00630     
00631     // Now determine the buffersizes needed for each neighbour using MPI_Pack_size
00632     MessageInformation dummy;
00633     
00634     MessageIterator messageIter= infoSend_.begin();
00635 
00636     typedef typename RemoteIndices::RemoteIndexMap::const_iterator RemoteIterator;
00637     const RemoteIterator rend = remoteIndices_.end();
00638     int neighbour=0;
00639     
00640     for(RemoteIterator remote = remoteIndices_.begin(); remote != rend; ++remote, ++neighbour){
00641       MessageInformation* message;
00642       MessageInformation recv;
00643       
00644       if(messageIter != end && messageIter->first==remote->first){
00645         // We want to send message information to that process
00646         message = const_cast<MessageInformation*>(&(messageIter->second));
00647         ++messageIter;
00648       }else
00649         // We do not want to send information but the other process might.
00650         message = &dummy;
00651 
00652       sendBufferSizes_[neighbour]=0;
00653       int tsize;
00654       // The number of indices published
00655       MPI_Pack_size(1, MPI_INT,remoteIndices_.communicator(), &tsize);
00656       sendBufferSizes_[neighbour] += tsize;
00657 
00658       for(int i=0; i < message->publish; ++i){
00659         // The global index
00660         MPI_Pack_size(1, MPITraits<GlobalIndex>::getType(), remoteIndices_.communicator(), &tsize);
00661         sendBufferSizes_[neighbour] += tsize;
00662         // The attribute in the local index
00663         MPI_Pack_size(1, MPI_CHAR, remoteIndices_.communicator(), &tsize);
00664         sendBufferSizes_[neighbour] += tsize;
00665         // The number of corresponding remote indices
00666         MPI_Pack_size(1, MPI_INT, remoteIndices_.communicator(), &tsize);
00667         sendBufferSizes_[neighbour] += tsize;
00668       }
00669       for(int i=0; i < message->pairs; ++i){
00670         // The process of the remote index
00671         MPI_Pack_size(1, MPI_INT, remoteIndices_.communicator(), &tsize);
00672         sendBufferSizes_[neighbour] += tsize;
00673         // The attribute of the remote index
00674         MPI_Pack_size(1, MPI_CHAR, remoteIndices_.communicator(), &tsize);
00675         sendBufferSizes_[neighbour] += tsize;
00676       }
00677           
00678       Dune::dverb<<rank_<<": Buffer (neighbour="<<remote->first<<") size is "<< sendBufferSizes_[neighbour]<<" for publish="<<message->publish<<" pairs="<<message->pairs<<std::endl;
00679     }
00680 
00681   }
00682   
00683   template<typename T>
00684   inline void IndicesSyncer<T>::sync()
00685   {
00686     DefaultNumberer numberer;
00687     sync(numberer);
00688   }
00689   
00690   template<typename T>
00691   template<typename T1>
00692   void IndicesSyncer<T>::sync(T1& numberer)
00693   {
00694     
00695     // The pointers to the local indices in the remote indices
00696     // will become invalid due to the resorting of the index set.
00697     // Therefore store the corresponding global indices.
00698     // Mark all indices as not added
00699     
00700     typedef typename RemoteIndices::RemoteIndexMap::const_iterator 
00701       RemoteIterator;
00702 
00703     const RemoteIterator end = remoteIndices_.end();
00704     
00705     // Number of neighbours might change during the syncing.
00706     // save the old neighbours
00707     std::size_t noOldNeighbours = remoteIndices_.neighbours();
00708     int* oldNeighbours = new int[noOldNeighbours];
00709     sendBufferSizes_ = new std::size_t[noOldNeighbours];
00710     std::size_t i=0;
00711     
00712     for(RemoteIterator remote = remoteIndices_.begin(); remote != end; ++remote, ++i){
00713       typedef typename RemoteIndices::RemoteIndexList::const_iterator
00714         RemoteIndexIterator;
00715 
00716       oldNeighbours[i]=remote->first;
00717 
00718       // Make sure we only have one remote index list.
00719       assert(remote->second.first==remote->second.second);
00720 
00721       RemoteIndexList& rList = *(remote->second.first);
00722             
00723       // Store the corresponding global indices.
00724       GlobalIndexList& global = globalMap_[remote->first];
00725       BoolList& added = oldMap_[remote->first];
00726       RemoteIndexIterator riEnd = rList.end();
00727       
00728       for(RemoteIndexIterator index = rList.begin();
00729           index != riEnd; ++index){
00730         global.push_back(index->localIndexPair().global());
00731         added.push_back(true);
00732       }
00733 
00734       Iterators iterators(rList, global, added);
00735       iteratorsMap_.insert(std::make_pair(remote->first, iterators));
00736       assert(checkReset(iteratorsMap_[remote->first], rList,global,added));
00737     }
00738     
00739     // Exchange indices with each neighbour
00740     const RemoteIterator rend = remoteIndices_.end();
00741 
00742     calculateMessageSizes();
00743     
00744     // Allocate the buffers
00745     receiveBufferSize_=1;
00746     sendBuffers_ = new char*[noOldNeighbours];
00747     
00748     for(std::size_t i=0; i<noOldNeighbours; ++i){
00749       sendBuffers_[i] = new char[sendBufferSizes_[i]];
00750       receiveBufferSize_ = std::max(receiveBufferSize_, static_cast<int>(sendBufferSizes_[i]));
00751     }
00752     
00753     receiveBuffer_=new char[receiveBufferSize_];
00754     
00755     indexSet_.beginResize();
00756     
00757     Dune::dverb<<rank_<<": Neighbours: ";
00758     
00759     for(i = 0; i<noOldNeighbours; ++i)
00760       Dune::dverb<<oldNeighbours[i]<<" ";
00761     
00762     Dune::dverb<<std::endl;
00763 
00764     MPI_Request* requests = new MPI_Request[noOldNeighbours];
00765     MPI_Status* statuses = new MPI_Status[noOldNeighbours];
00766 
00767     // Pack Message data and start the sends
00768     for(i = 0; i<noOldNeighbours; ++i)
00769         packAndSend(oldNeighbours[i], sendBuffers_[i], sendBufferSizes_[i], requests[i]);
00770 
00771     // Probe for incoming messages, receive and unpack them
00772     for(i = 0; i<noOldNeighbours; ++i)
00773       recvAndUnpack(numberer);
00774 //       }else{
00775 //      recvAndUnpack(oldNeighbours[i], numberer);
00776 //      packAndSend(oldNeighbours[i]);
00777 //       }
00778 //     }
00779 
00780     delete[] receiveBuffer_;
00781     
00782     // Wait for the completion of the sends
00783     // Wait for completion of sends
00784     if(MPI_SUCCESS!=MPI_Waitall(noOldNeighbours, requests, statuses)){
00785       std::cerr<<": MPI_Error occurred while sending message"<<std::endl;
00786       for(i=0;i< noOldNeighbours;i++)
00787         if(MPI_SUCCESS!=statuses[i].MPI_ERROR)
00788           std::cerr<<"Destination "<<statuses[i].MPI_SOURCE<<" error code: "<<statuses[i].MPI_ERROR<<std::endl;
00789     }
00790   
00791     delete[] statuses;
00792     delete[] requests;
00793     
00794     for(std::size_t i=0; i<noOldNeighbours; ++i)
00795       delete[] sendBuffers_[i];
00796 
00797     delete[] sendBuffers_;
00798     delete[] sendBufferSizes_;
00799     
00800     // No need for the iterator tuples any more
00801     iteratorsMap_.clear();
00802     
00803     indexSet_.endResize();
00804     
00805     delete[] oldNeighbours;
00806     
00807     repairLocalIndexPointers<T,typename RemoteIndices::Allocator>(globalMap_, remoteIndices_, indexSet_);
00808     
00809     oldMap_.clear();    
00810     globalMap_.clear();
00811         
00812     // update the sequence number
00813     remoteIndices_.sourceSeqNo_ = remoteIndices_.destSeqNo_ = indexSet_.seqNo();    
00814   }
00815   
00816   template<typename T>
00817   void IndicesSyncer<T>::packAndSend(int destination, char* buffer, std::size_t bufferSize, MPI_Request& request)
00818   {
00819     typedef typename ParallelIndexSet::const_iterator IndexIterator;
00820     typedef typename RemoteIndexList::const_iterator RemoteIndexIterator;
00821     typedef typename GlobalIndexList::const_iterator GlobalIterator;
00822     typedef typename BoolList::const_iterator BoolIterator;
00823         
00824     IndexIterator      iEnd       = indexSet_.end();
00825     int                bpos       = 0;
00826     int                published  = 0;
00827     int                pairs      = 0;
00828     
00829     assert(checkReset());
00830     
00831     // Pack the number of indices we publish
00832     MPI_Pack(&(infoSend_[destination].publish), 1, MPI_INT, buffer, bufferSize, &bpos, 
00833              remoteIndices_.communicator());
00834 
00835     for(IndexIterator index = indexSet_.begin(); index != iEnd; ++index){
00836       // Search for corresponding remote indices in all iterator tuples
00837       typedef typename IteratorsMap::iterator Iterator;
00838       Iterator iteratorsEnd = iteratorsMap_.end();
00839       
00840       // advance all iterators to a position with global index >= index->global()
00841       for(Iterator iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
00842         while(iterators->second.isNotAtEnd() && iterators->second.globalIndex() < index->global())
00843           ++(iterators->second);
00844           
00845       // Add all remote indices positioned at global which were already present before calling sync
00846       // to the message.
00847       // Count how many remote indices we will send
00848       int indices = 0;
00849       bool knownRemote = false; // Is the remote process supposed to know this index?
00850       
00851       for(Iterator iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
00852         if(iterators->second.isNotAtEnd() && iterators->second.isOld() 
00853            && iterators->second.globalIndex() == index->global()){
00854           indices++;
00855           if(destination == iterators->first)
00856             knownRemote = true;
00857         }
00858       
00859       if(!knownRemote || indices==0)
00860         // We do not need to send any indices
00861         continue;
00862       
00863       Dune::dverb<<rank_<<": sending "<<indices<<" for index "<<index->global()<<" to "<<destination<<std::endl;
00864 
00865       pairs+=indices;
00866       assert(pairs <= infoSend_[destination].pairs);
00867 
00868       // Pack the global index, the attribute and the number
00869       MPI_Pack(const_cast<GlobalIndex*>(&(index->global())), 1, MPITraits<GlobalIndex>::getType(), buffer, bufferSize, &bpos, 
00870                remoteIndices_.communicator());
00871       
00872       char attr = index->local().attribute();
00873       MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos, 
00874                remoteIndices_.communicator());
00875 
00876       // Pack the number of remote indices we send.
00877       MPI_Pack(&indices, 1, MPI_INT, buffer, bufferSize, &bpos, 
00878                remoteIndices_.communicator());
00879         
00880       // Pack the information about the remote indices
00881       for(Iterator iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
00882         if(iterators->second.isNotAtEnd() && iterators->second.isOld() 
00883            && iterators->second.globalIndex() == index->global()){
00884           int process = iterators->first;
00885           MPI_Pack(&process, 1, MPI_INT, buffer, bufferSize, &bpos, 
00886                    remoteIndices_.communicator());
00887           char attr = iterators->second.remoteIndex().attribute();
00888           
00889           MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos, 
00890                    remoteIndices_.communicator());
00891         }
00892       ++published;
00893       assert(published <= infoSend_[destination].publish);
00894     }
00895 
00896     // Make sure we send all expected entries
00897     assert(published == infoSend_[destination].publish);
00898 
00899     resetIteratorsMap();
00900 
00901     Dune::dverb << rank_<<": Sending message of "<<bpos<<" bytes to "<<destination<<std::endl;
00902     
00903     MPI_Issend(buffer, bpos, MPI_PACKED, destination, 345, remoteIndices_.communicator(),&request);
00904   }
00905 
00906   template<typename T>
00907   inline void IndicesSyncer<T>::insertIntoRemoteIndexList(int process, const GlobalIndex& global, 
00908                                                                 char attribute)
00909   {
00910     Dune::dvverb<<"Inserting from "<<process<<" "<<global<<" "<<attribute<<std::endl;
00911     
00912     // There might be cases where there no remote indices for that process yet
00913     typename IteratorsMap::iterator found = iteratorsMap_.find(process);
00914     
00915     if( found == iteratorsMap_.end() ){
00916       Dune::dverb<<"Discovered new neighbour "<<process<<std::endl;
00917       RemoteIndexList* rlist = new RemoteIndexList();
00918       remoteIndices_.remoteIndices_.insert(std::make_pair(process,std::make_pair(rlist,rlist)));
00919       Iterators iterators = Iterators(*rlist, globalMap_[process], oldMap_[process]);
00920       found = iteratorsMap_.insert(std::make_pair(process, iterators)).first;
00921     }
00922     
00923     Iterators& iterators = found->second;
00924     
00925     // Search for the remote index
00926     while(iterators.isNotAtEnd() && iterators.globalIndex() < global){
00927       // Increment all iterators
00928       ++iterators;
00929       
00930     }
00931     
00932     if(iterators.isAtEnd() || iterators.globalIndex() != global){
00933       // The entry is not yet known
00934       // Insert in the the list and do not change the first iterator.
00935       iterators.insert(RemoteIndex(Attribute(attribute)),global);
00936     }else if(iterators.isNotAtEnd()){
00937       // Assert that the attributes match 
00938       assert(iterators.remoteIndex().attribute() == attribute);
00939     }
00940   }
00941   
00942   template<typename T>
00943   template<typename T1>
00944   void IndicesSyncer<T>::recvAndUnpack(T1& numberer)
00945   {
00946     typedef typename ParallelIndexSet::const_iterator IndexIterator;
00947     typedef typename RemoteIndexList::iterator RemoteIndexIterator;
00948     typedef typename GlobalIndexList::iterator GlobalIndexIterator;
00949 
00950     IndexIterator    iEnd   = indexSet_.end();
00951     IndexIterator    index  = indexSet_.begin();
00952     int              bpos   = 0;
00953     int              publish;
00954     
00955     assert(checkReset());
00956 
00957     MPI_Status status;
00958 
00959     // We have to determine the message size and source before the receive
00960     MPI_Probe(MPI_ANY_SOURCE, 345, remoteIndices_.communicator(), &status);
00961     
00962     int source=status.MPI_SOURCE;
00963     int count;
00964     MPI_Get_count(&status, MPI_PACKED, &count);
00965 
00966     Dune::dvverb<<rank_<<": Receiving message from "<< source<<" with "<<count<<" bytes"<<std::endl;
00967 
00968     if(count>receiveBufferSize_){
00969       receiveBufferSize_=count;
00970       delete[] receiveBuffer_;
00971       receiveBuffer_ = new char[receiveBufferSize_];
00972     }
00973 
00974     MPI_Recv(receiveBuffer_, count, MPI_PACKED, source, 345, remoteIndices_.communicator(), &status);
00975     
00976     // How many global entries were published?
00977     MPI_Unpack(receiveBuffer_,  count, &bpos, &publish, 1, MPI_INT, remoteIndices_.communicator());
00978 
00979     // Now unpack the remote indices and add them.
00980     while(publish>0){
00981       
00982       // Unpack information about the local index on the source process
00983       GlobalIndex  global;          // global index of the current entry
00984       char sourceAttribute; // Attribute on the source process
00985       int pairs;
00986       
00987       MPI_Unpack(receiveBuffer_,  count, &bpos, &global, 1, MPITraits<GlobalIndex>::getType(), 
00988                  remoteIndices_.communicator());
00989       MPI_Unpack(receiveBuffer_,  count, &bpos, &sourceAttribute, 1, MPI_CHAR, 
00990                  remoteIndices_.communicator());
00991       MPI_Unpack(receiveBuffer_,  count, &bpos, &pairs, 1, MPI_INT, 
00992                  remoteIndices_.communicator());
00993     
00994       // Insert the entry on the remote process to our
00995       // remote index list
00996       insertIntoRemoteIndexList(source, global, sourceAttribute);
00997 
00998       // Unpack the remote indices
00999       while(pairs>0){
01000         // Unpack the process id that knows the index
01001         int process;
01002         char attribute;
01003         MPI_Unpack(receiveBuffer_,  count, &bpos, &process, 1, MPI_INT, 
01004                  remoteIndices_.communicator());
01005         // Unpack the attribute
01006         MPI_Unpack(receiveBuffer_,  count, &bpos, &attribute, 1, MPI_CHAR, 
01007                    remoteIndices_.communicator());
01008         
01009         if(process==rank_){
01010           // Now we know the local attribute of the global index
01011           // Do we know that global index already?
01012           IndexIterator pos = std::lower_bound(index, iEnd, IndexPair(global));
01013                   
01014           if(pos == iEnd || pos->global() != global){
01015             // No, we do not. Add it!
01016             indexSet_.add(global,ParallelLocalIndex<Attribute>(numberer(global),
01017                                                         Attribute(attribute), true));
01018           }else{
01019             // Attributes have to match!
01020             assert(attribute==pos->local().attribute());
01021           }
01022           index=pos;
01023         }else{
01024           insertIntoRemoteIndexList(process, global, attribute);
01025         }
01026         
01027         --pairs;
01028       }
01029       --publish;
01030     }
01031     
01032     resetIteratorsMap();
01033   }
01034   
01035   template<typename T>
01036   void IndicesSyncer<T>::resetIteratorsMap(){
01037     
01038     // Reset iterators in all tuples.
01039     typedef typename IteratorsMap::iterator Iterator;
01040     typedef typename RemoteIndices::RemoteIndexMap::iterator 
01041       RemoteIterator;
01042     typedef typename GlobalIndicesMap::iterator GlobalIterator;
01043     typedef typename BoolMap::iterator BoolIterator;
01044     
01045     const RemoteIterator remoteEnd = remoteIndices_.remoteIndices_.end();
01046     Iterator iterators = iteratorsMap_.begin();
01047     GlobalIterator global = globalMap_.begin();
01048     BoolIterator added = oldMap_.begin();
01049     
01050     for(RemoteIterator remote = remoteIndices_.remoteIndices_.begin();
01051         remote != remoteEnd; ++remote, ++global, ++added, ++iterators){
01052       iterators->second.reset(*(remote->second.first), global->second, added->second);
01053     }
01054   }
01055  
01056   template<typename T>
01057   bool IndicesSyncer<T>::checkReset(const Iterators& iterators, RemoteIndexList& rList, GlobalIndexList& gList,
01058                                           BoolList& bList){
01059     
01060     if(get<0>(iterators.iterators_)!=rList.begin())
01061       return false;
01062     if(get<1>(iterators.iterators_)!=gList.begin())
01063       return false;
01064     if(get<2>(iterators.iterators_)!=bList.begin())
01065       return false;
01066     return true;
01067   }
01068   
01069 
01070   template<typename T>
01071   bool IndicesSyncer<T>::checkReset(){
01072     
01073     // Reset iterators in all tuples.
01074     typedef typename IteratorsMap::iterator Iterator;
01075     typedef typename RemoteIndices::RemoteIndexMap::iterator 
01076       RemoteIterator;
01077     typedef typename GlobalIndicesMap::iterator GlobalIterator;
01078     typedef typename BoolMap::iterator BoolIterator;
01079     
01080     const RemoteIterator remoteEnd = remoteIndices_.remoteIndices_.end();
01081     Iterator iterators = iteratorsMap_.begin();
01082     GlobalIterator global = globalMap_.begin();
01083     BoolIterator added = oldMap_.begin();
01084     bool ret = true;
01085     
01086     for(RemoteIterator remote = remoteIndices_.remoteIndices_.begin();
01087         remote != remoteEnd; ++remote, ++global, ++added, ++iterators){
01088       if(!checkReset(iterators->second, *(remote->second.first), global->second,
01089                      added->second))
01090         ret=false;
01091     }
01092     return ret;
01093   } 
01094 }
01095 
01096 #endif
01097 #endif
Generated on Sat Apr 24 11:13:46 2010 for dune-istl by  doxygen 1.6.3