remoteindices.hh

Go to the documentation of this file.
00001 // $Id: remoteindices.hh 983 2009-01-28 17:20:31Z mblatt $
00002 #ifndef DUNE_REMOTEINDICES_HH
00003 #define DUNE_REMOTEINDICES_HH
00004 
00005 #include"indexset.hh"
00006 #include"istlexception.hh"
00007 #include"plocalindex.hh"
00008 #include<dune/common/poolallocator.hh>
00009 #include<dune/common/sllist.hh>
00010 #include<dune/common/static_assert.hh>
00011 #include<map>
00012 #include<utility>
00013 #include<iostream>
00014 
00015 #if HAVE_MPI
00016 #include"mpitraits.hh"
00017 #include"mpi.h"
00018 
00019 namespace Dune{
00030 
00031   template<typename TG, typename TA>
00032   class MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >
00033   {
00034   public:
00035     inline static MPI_Datatype getType();
00036   private:
00037     static MPI_Datatype type;
00038   };
00039   
00040   
00041   template<typename T>
00042   class RemoteIndices;
00043   
00044   template<typename T1, typename T2>
00045   class RemoteIndex;
00046   
00047   template<typename T>
00048   class IndicesSyncer;
00049   
00050   template<typename T1, typename T2>
00051   std::ostream& operator<<(std::ostream& os, const RemoteIndex<T1,T2>& index);
00052   
00053 
00054   template<typename T, bool mode>
00055   class RemoteIndexListModifier;
00056 
00060   template<typename T1, typename T2>
00061   class RemoteIndex
00062   {
00063     template<typename T>
00064     friend class IndicesSyncer;
00065     
00066     template<typename T, typename A>
00067     friend void repairLocalIndexPointers(std::map<int,SLList<typename T::GlobalIndex,A> >&, RemoteIndices<T>&,
00068                                          const T&); 
00069 
00070     template<typename T, bool mode>
00071     friend class RemoteIndexListModifier;
00072 
00073   public:    
00078     typedef T1 GlobalIndex;
00087     typedef T2 Attribute;
00088 
00092     typedef IndexPair<GlobalIndex,ParallelLocalIndex<Attribute> >
00093     PairType;
00094     
00099     const Attribute attribute() const;
00100 
00106     const PairType& localIndexPair() const;
00107 
00111     RemoteIndex();
00112     
00113     
00119     RemoteIndex(const T2& attribute, 
00120                 const PairType* local);
00121 
00122     
00128     RemoteIndex(const T2& attribute);
00129 
00130     bool operator==(const RemoteIndex& ri) const;
00131     
00132     bool operator!=(const RemoteIndex& ri) const;
00133   private:
00135     const PairType* localIndex_;
00136 
00138     char attribute_;
00139   };
00140   
00141   template<class T>
00142   std::ostream& operator<<(std::ostream& os, const RemoteIndices<T>& indices);
00143   
00144   template<typename T>
00145   class InterfaceBuilder;
00146 
00147   template<class T>
00148   class CollectiveIterator;
00149 
00150   template<class T>
00151     class IndicesSyncer;
00152 
00166   template<class T>
00167   class RemoteIndices
00168   {
00169     friend class InterfaceBuilder<T>;
00170     friend class IndicesSyncer<T>;
00171     template<typename T1, typename A>
00172     friend void repairLocalIndexPointers(std::map<int,SLList<typename T1::GlobalIndex,A> >&, 
00173                                          RemoteIndices<T1>&,
00174                                          const T1&);
00175     friend std::ostream& operator<<<>(std::ostream&, const RemoteIndices<T>&);
00176     
00177   public:
00178    
00182     typedef T ParallelIndexSet;
00183  
00187     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00188     
00189     
00193     typedef typename ParallelIndexSet::LocalIndex LocalIndex;
00194 
00198     typedef typename LocalIndex::Attribute Attribute;
00199  
00203     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00204 
00205     
00209     typedef Dune::PoolAllocator<RemoteIndex,1> Allocator;
00210 
00212     typedef Dune::SLList<RemoteIndex,Allocator>
00213     RemoteIndexList;
00214     
00216     typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
00217     RemoteIndexMap;    
00218     
00219     typedef typename RemoteIndexMap::const_iterator const_iterator;
00220     
00234     inline RemoteIndices(const ParallelIndexSet& source, const ParallelIndexSet& destination, 
00235                          const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>());
00236 
00237     RemoteIndices();
00238     
00255     void setIndexSets(const ParallelIndexSet& source, const ParallelIndexSet& destination, 
00256                       const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>());
00257 
00258     void setNeighbours(const std::vector<int>& neighbours)
00259     {
00260       neighbourIds=neighbours;
00261     }
00262     
00266     ~RemoteIndices();
00267     
00277     template<bool ignorePublic>
00278     void rebuild();
00279 
00280     bool operator==(const RemoteIndices& ri);
00281     
00289     inline bool isSynced() const;
00290 
00294     inline MPI_Comm communicator() const;
00295 
00305     template<bool mode, bool send>
00306     inline RemoteIndexListModifier<T,mode> getModifier(int process);
00307     
00312     inline const_iterator begin() const;
00313     
00318     inline const_iterator end() const;
00319     
00323     template<bool send>
00324     inline CollectiveIterator<T> iterator() const;
00325 
00329     inline void free();
00330     
00335     inline int neighbours() const;
00336     
00337   private:    
00339     const ParallelIndexSet* source_;
00340 
00342     const ParallelIndexSet* target_;
00343 
00345     MPI_Comm comm_;
00346 
00349     std::vector<int> neighbourIds;
00350     
00352     const static int commTag_=333;
00353     
00358     int sourceSeqNo_;
00359     
00364     int destSeqNo_;
00365 
00369     bool publicIgnored;
00370 
00374     bool firstBuild;
00375         
00377     typedef IndexPair<GlobalIndex, LocalIndex> 
00378     PairType;
00379     
00386     RemoteIndexMap remoteIndices_;
00387     
00394     template<bool ignorePublic>
00395     inline void buildRemote();
00396 
00402     inline int noPublic(const ParallelIndexSet& indexSet);
00403     
00415     template<bool ignorePublic>
00416     inline void packEntries(PairType** myPairs, const ParallelIndexSet& indexSet,
00417                             char* p_out, MPI_Datatype type, int bufferSize, 
00418                             int* position, int n);
00419     
00433     inline void unpackIndices(RemoteIndexList& remote, int remoteEntries,
00434                               PairType** local, int localEntries, char* p_in, 
00435                               MPI_Datatype type, int* positon, int bufferSize);
00436 
00437     inline void unpackIndices(RemoteIndexList& send, RemoteIndexList& receive,
00438                               int remoteEntries, PairType** localSource,
00439                               int localSourceEntries, PairType** localDest,
00440                               int localDestEntries, char* p_in,
00441                               MPI_Datatype type, int* position, int bufferSize);
00442 
00443     void unpackCreateRemote(char* p_in, PairType** sourcePairs, PairType** DestPairs, 
00444                             int remoteProc,  int sourcePublish, int destPublish, 
00445                             int bufferSize, bool sendTwo);
00446   };
00447   
00465   template<class T, bool mode>
00466   class RemoteIndexListModifier
00467   {
00468 
00469     template<typename T1>
00470     friend class RemoteIndices;
00471   
00472   public:
00473     class InvalidPosition : public ISTLError
00474     {};
00475 
00476     enum {
00485       MODIFYINDEXSET=mode
00486     };
00487   
00491     typedef T ParallelIndexSet;
00492 
00496     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00497 
00501     typedef typename ParallelIndexSet::LocalIndex LocalIndex;
00502     
00506     typedef typename LocalIndex::Attribute Attribute;
00507     
00511     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00512     
00516     typedef Dune::PoolAllocator<RemoteIndex,1> Allocator;
00517 
00519     typedef Dune::SLList<RemoteIndex,Allocator>
00520     RemoteIndexList;
00521 
00525     typedef SLListModifyIterator<RemoteIndex,Allocator> ModifyIterator;
00526     
00530     typedef typename RemoteIndexList::const_iterator ConstIterator;
00531     
00545     void insert(const RemoteIndex& index) throw(InvalidPosition);
00546     
00547     
00562     void insert(const RemoteIndex& index, const GlobalIndex& global) throw(InvalidPosition);
00563 
00571     bool remove(const GlobalIndex& global) throw(InvalidPosition);
00572 
00585     void repairLocalIndexPointers() throw(InvalidIndexSetState);
00586 
00587 
00588     RemoteIndexListModifier(const RemoteIndexListModifier<T,mode>&);
00589 
00590   private:
00591     
00597     RemoteIndexListModifier(const ParallelIndexSet& indexSet,
00598                             RemoteIndexList& rList);
00599 
00600     typedef SLList<GlobalIndex,Allocator> GlobalList;
00601     typedef typename GlobalList::ModifyIterator GlobalModifyIterator;
00602     RemoteIndices<ParallelIndexSet>* remoteIndices_;
00603     RemoteIndexList* rList_;
00604     const ParallelIndexSet* indexSet_;
00605     GlobalList* glist_;
00606     ModifyIterator iter_;
00607     GlobalModifyIterator giter_;
00608     ConstIterator end_;
00609     bool first_;
00610     GlobalIndex last_;
00611   };
00612   
00617   template<class T>
00618   class CollectiveIterator
00619   {
00620     
00624     typedef T ParallelIndexSet;
00625 
00629     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00630 
00634     typedef typename ParallelIndexSet::LocalIndex LocalIndex;
00635     
00639     typedef typename LocalIndex::Attribute Attribute;
00640 
00642     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00643     
00645     typedef Dune::PoolAllocator<RemoteIndex,1> Allocator;
00646 
00648     typedef Dune::SLList<RemoteIndex,Allocator> RemoteIndexList;    
00649 
00651     typedef std::map<int,std::pair<typename RemoteIndexList::const_iterator,
00652                                    const typename RemoteIndexList::const_iterator> >
00653     Map;
00654 
00655   public:
00656     
00658     typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
00659     RemoteIndexMap;
00660     
00666     inline CollectiveIterator(const RemoteIndexMap& map_, bool send);
00667     
00676     inline void advance(const GlobalIndex& global);
00677     
00681     inline bool empty();
00682     
00689     class iterator
00690     {
00691     public:
00692       typedef typename Map::iterator RealIterator;
00693       typedef typename Map::iterator ConstRealIterator;
00694       
00695       
00697       iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex& index)
00698         : iter_(iter), end_(end), index_(index)
00699       {
00700         // Move to the first valid entry
00701         while(iter_!=end_ && iter_->second.first->localIndexPair().global()!=index_)
00702           ++iter_;
00703       }
00704       
00706       iterator(const iterator& other)
00707         : iter_(other.iter_), end_(other.end_), index_(other.index_)
00708       { }
00709          
00711       iterator& operator++()
00712       { 
00713         ++iter_;
00714         // If entry is not valid move on
00715         while(iter_!=end_ && iter_->second.first->localIndexPair().global()!=index_)
00716           ++iter_;
00717                   
00718         return *this;
00719       }
00720       
00722       const RemoteIndex& operator*()const
00723       {
00724         return *(iter_->second.first);
00725       }
00726       
00728       int process() const
00729       {
00730         return iter_->first;
00731       }
00732       
00734       const RemoteIndex* operator->()const
00735       {
00736         return iter_->second.first.operator->();
00737       }
00738       
00740       bool operator==(const iterator& other)
00741       {
00742         return other.iter_==iter_;
00743       }
00744 
00746       bool operator!=(const iterator& other)
00747       {
00748         return other.iter_!=iter_;
00749       }
00750 
00751     private:
00752       iterator();
00753       
00754       RealIterator iter_;
00755       RealIterator end_;
00756       GlobalIndex index_;
00757     };
00758     
00759     iterator begin();
00760     
00761     iterator end();
00762     
00763   private:
00764     
00765     Map map_;
00766     GlobalIndex index_;
00767   };
00768     
00769   template<typename TG, typename TA>
00770   MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::getType()
00771   {
00772     if(type==MPI_DATATYPE_NULL){
00773       int length[4];
00774       MPI_Aint disp[4];
00775       MPI_Datatype types[4] = {MPI_LB, MPITraits<TG>::getType(), 
00776                                MPITraits<ParallelLocalIndex<TA> >::getType(), MPI_UB};
00777       IndexPair<TG,ParallelLocalIndex<TA> > rep[2];
00778       length[0]=length[1]=length[2]=length[3]=1;
00779       MPI_Address(rep, disp); // lower bound of the datatype
00780       MPI_Address(&(rep[0].global_), disp+1);
00781       MPI_Address(&(rep[0].local_), disp+2);
00782       MPI_Address(rep+1, disp+3); // upper bound of the datatype
00783       for(int i=3; i >= 0; --i)
00784         disp[i] -= disp[0];
00785       MPI_Type_struct(4, length, disp, types, &type);
00786       MPI_Type_commit(&type);
00787     }
00788     return type;
00789   }
00790 
00791   template<typename TG, typename TA>
00792   MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::type=MPI_DATATYPE_NULL;
00793 
00794   template<typename T1, typename T2>
00795   RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute, const PairType* local)
00796     : localIndex_(local), attribute_(attribute)
00797   {}
00798 
00799   template<typename T1, typename T2>
00800   RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute)
00801     : localIndex_(0), attribute_(attribute)
00802   {}
00803 
00804   template<typename T1, typename T2>
00805   RemoteIndex<T1,T2>::RemoteIndex()
00806     : localIndex_(0), attribute_()
00807   {}
00808   template<typename T1, typename T2>
00809   inline bool RemoteIndex<T1,T2>::operator==(const RemoteIndex& ri) const
00810   {
00811     return localIndex_==ri.localIndex_ && attribute_==ri.attribute;
00812   }
00813   
00814   template<typename T1, typename T2>
00815   inline bool RemoteIndex<T1,T2>::operator!=(const RemoteIndex& ri) const
00816   {
00817     return localIndex_!=ri.localIndex_ || attribute_!=ri.attribute_;
00818   }
00819   
00820   template<typename T1, typename T2>
00821   inline const T2 RemoteIndex<T1,T2>::attribute() const
00822   {
00823     return T2(attribute_);
00824   }
00825 
00826   template<typename T1, typename T2>
00827   inline const IndexPair<T1,ParallelLocalIndex<T2> >& RemoteIndex<T1,T2>::localIndexPair() const
00828   {
00829     return *localIndex_;
00830   }
00831   
00832   template<typename T>
00833   inline RemoteIndices<T>::RemoteIndices(const ParallelIndexSet& source,
00834                                          const ParallelIndexSet& destination,
00835                                          const MPI_Comm& comm,
00836                                          const std::vector<int>& neighbours)
00837     : source_(&source), target_(&destination), comm_(comm), neighbourIds(neighbours),
00838       sourceSeqNo_(-1), destSeqNo_(-1), publicIgnored(false), firstBuild(true)
00839   {}
00840 
00841   template<typename T>
00842   RemoteIndices<T>::RemoteIndices()
00843     :source_(0), target_(0), sourceSeqNo_(-1), 
00844      destSeqNo_(-1), publicIgnored(false), firstBuild(true)
00845   {}
00846   
00847   template<class T>
00848   void RemoteIndices<T>::setIndexSets(const ParallelIndexSet& source,
00849                                       const ParallelIndexSet& destination,
00850                                       const MPI_Comm& comm,
00851                                       const std::vector<int>& neighbours)
00852   {
00853     free();
00854     source_ = &source;
00855     target_ = &destination;
00856     comm_ = comm;
00857     firstBuild = true;
00858     neighbourIds=neighbours;
00859   }
00860   
00861   template<typename T>
00862   RemoteIndices<T>::~RemoteIndices()
00863   {
00864     free();
00865   }
00866 
00867   template<typename T>
00868   template<bool ignorePublic>
00869   inline void RemoteIndices<T>::packEntries(IndexPair<GlobalIndex,LocalIndex>** pairs,
00870                                                 const ParallelIndexSet& indexSet,
00871                                                 char* p_out, MPI_Datatype type, 
00872                                                 int bufferSize,
00873                                                 int *position, int n)
00874   {
00875     // fill with own indices
00876     typedef typename ParallelIndexSet::const_iterator const_iterator;
00877     typedef IndexPair<GlobalIndex,LocalIndex> PairType;
00878     const const_iterator end = indexSet.end();
00879 
00880     //Now pack the source indices
00881     int i=0;
00882     for(const_iterator index = indexSet.begin(); index != end; ++index)
00883       if(ignorePublic || index->local().isPublic()){
00884         
00885         MPI_Pack(const_cast<PairType*>(&(*index)), 1, 
00886                  type,
00887                  p_out, bufferSize, position, comm_);
00888         pairs[i++] = const_cast<PairType*>(&(*index));
00889       
00890       }
00891     assert(i==n);
00892   }
00893   
00894   template<typename T>
00895   inline int RemoteIndices<T>::noPublic(const ParallelIndexSet& indexSet)
00896   {
00897     typedef typename ParallelIndexSet::const_iterator const_iterator;
00898     
00899     int noPublic=0;
00900     
00901     const const_iterator end=indexSet.end();
00902     for(const_iterator index=indexSet.begin(); index!=end; ++index)
00903       if(index->local().isPublic())
00904         noPublic++;
00905     
00906     return noPublic;
00907     
00908   }
00909     
00910 
00911   template<typename T>
00912   inline void RemoteIndices<T>::unpackCreateRemote(char* p_in, PairType** sourcePairs,
00913                                                    PairType** destPairs, int remoteProc, 
00914                                                    int sourcePublish, int destPublish,
00915                                                    int bufferSize, bool sendTwo)
00916   {
00917     
00918       // unpack the number of indices we received
00919       int noRemoteSource=-1, noRemoteDest=-1;
00920       char twoIndexSets=0;
00921       int position=0;
00922       // Did we receive two index sets?
00923       MPI_Unpack(p_in, bufferSize, &position, &twoIndexSets, 1, MPI_CHAR, comm_);
00924       // The number of source indices received
00925       MPI_Unpack(p_in, bufferSize, &position, &noRemoteSource, 1, MPI_INT, comm_);
00926       // The number of destination indices received
00927       MPI_Unpack(p_in, bufferSize, &position, &noRemoteDest, 1, MPI_INT, comm_);      
00928 
00929 
00930       // Indices for which we receive
00931       RemoteIndexList* receive= new RemoteIndexList();
00932       // Indices for which we send
00933       RemoteIndexList* send=0;
00934 
00935       MPI_Datatype type= MPITraits<PairType>::getType();
00936       
00937       if(!twoIndexSets){
00938         if(sendTwo){
00939           send = new RemoteIndexList();
00940           // Create both remote index sets simultaneously
00941           unpackIndices(*send, *receive, noRemoteSource, sourcePairs, sourcePublish,
00942                         destPairs, destPublish, p_in, type, &position, bufferSize);
00943         }else{
00944           // we only need one list
00945           unpackIndices(*receive, noRemoteSource, sourcePairs, sourcePublish,
00946                         p_in, type, &position, bufferSize);
00947           send=receive;
00948         }
00949       }else{
00950         
00951         int oldPos=position;
00952         // Two index sets received
00953         unpackIndices(*receive, noRemoteSource, destPairs, destPublish,
00954                       p_in, type, &position, bufferSize);
00955         if(!sendTwo)
00956           //unpack source entries again as destination entries
00957           position=oldPos;
00958         
00959         send = new RemoteIndexList();
00960         unpackIndices(*send, noRemoteDest, sourcePairs, sourcePublish, 
00961                       p_in, type, &position, bufferSize);
00962       }
00963       
00964       if(receive->empty() && send->empty()){
00965         if(send==receive){
00966           delete send;
00967         }else{
00968           delete send;
00969           delete receive;
00970         }
00971       }else{
00972         remoteIndices_.insert(std::make_pair(remoteProc, 
00973                                              std::make_pair(send,receive)));
00974       }
00975   }
00976   
00977 
00978   template<typename T>
00979   template<bool ignorePublic>
00980   inline void RemoteIndices<T>::buildRemote()
00981   {
00982     // Processor configuration
00983     int rank, procs;
00984     MPI_Comm_rank(comm_, &rank);
00985     MPI_Comm_size(comm_, &procs);
00986     
00987     // number of local indices to publish
00988     // The indices of the destination will be send.
00989     int sourcePublish, destPublish;
00990     
00991     // Do we need to send two index sets?
00992     char sendTwo = (source_ != target_);
00993     
00994     if(procs==0 && !sendTwo)
00995       // Nothing to communicate
00996       return;   
00997 
00998     sourcePublish = (ignorePublic)? source_->size() : noPublic(*source_);
00999     
01000     if(sendTwo)
01001       destPublish = (ignorePublic)? target_->size() : noPublic(*target_);
01002     else  
01003       // we only need to send one set of indices
01004       destPublish = 0;
01005         
01006     int maxPublish, publish=sourcePublish+destPublish;
01007 
01008     // Calucate maximum number of indices send
01009     MPI_Allreduce(&publish, &maxPublish, 1, MPI_INT, MPI_MAX, comm_);
01010 
01011     // allocate buffers
01012     typedef IndexPair<GlobalIndex,LocalIndex> PairType;
01013     
01014     PairType** destPairs;
01015     PairType** sourcePairs = new PairType*[sourcePublish];
01016     
01017     if(sendTwo)
01018       destPairs = new PairType*[destPublish];
01019     else
01020       destPairs=sourcePairs;
01021     
01022     char** buffer = new char*[2];
01023     int bufferSize;    
01024     int position=0;
01025     int intSize;
01026     int charSize;
01027 
01028     // calculate buffer size
01029     MPI_Datatype type = MPITraits<PairType>::getType();
01030     
01031     MPI_Pack_size(maxPublish, type, comm_,
01032                   &bufferSize);
01033     MPI_Pack_size(1, MPI_INT, comm_,
01034                   &intSize);
01035     MPI_Pack_size(1, MPI_CHAR, comm_,
01036                   &charSize);
01037     // Our message will contain the following:
01038     // a bool wether two index sets where sent
01039     // the size of the source and the dest indexset, 
01040     // then the source and destination indices
01041     bufferSize += 2 * intSize + charSize;
01042 
01043     buffer[0] = new char[bufferSize];
01044     buffer[1] = new char[bufferSize];
01045     
01046     
01047     // pack entries into buffer[0], p_out below!
01048     MPI_Pack(&sendTwo, 1, MPI_CHAR, buffer[0], bufferSize, &position,
01049                  comm_);
01050         
01051     // The number of indices we send for each index set
01052     MPI_Pack(&sourcePublish, 1, MPI_INT, buffer[0], bufferSize, &position, 
01053              comm_);
01054     MPI_Pack(&destPublish, 1, MPI_INT, buffer[0], bufferSize, &position,
01055              comm_);
01056         
01057     // Now pack the source indices and setup the destination pairs
01058     packEntries<ignorePublic>(sourcePairs, *source_, buffer[0], type, 
01059                               bufferSize, &position, sourcePublish);
01060     // If necessary send the dest indices and setup the source pairs
01061     if(sendTwo)
01062       packEntries<ignorePublic>(destPairs, *target_, buffer[0], type,
01063                                 bufferSize, &position, destPublish);
01064 
01065 
01066     // Update remote indices for ourself
01067     if(sendTwo)
01068       unpackCreateRemote(buffer[0], sourcePairs, destPairs, rank, sourcePublish, 
01069                          destPublish, bufferSize, sendTwo);
01070 
01071     if(neighbourIds.size()==0)
01072       {
01073         // send mesages in ring
01074         for(int proc=1; proc<procs; proc++){
01075           // pointers to the current input and output buffers
01076           char* p_out = buffer[1-(proc%2)];
01077           char* p_in = buffer[proc%2];
01078           
01079           MPI_Status status;
01080           if(rank%2==0){
01081             MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
01082                       commTag_, comm_);
01083             MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
01084                      commTag_, comm_, &status);
01085           }else{
01086             MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
01087                      commTag_, comm_, &status);
01088             MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
01089                       commTag_, comm_);
01090           }
01091 
01092 
01093           // The process these indices are from
01094           int remoteProc = (rank+procs-proc)%procs;
01095           
01096           unpackCreateRemote(p_in, sourcePairs, destPairs, remoteProc, sourcePublish, 
01097                              destPublish, bufferSize, sendTwo);
01098           
01099         }
01100         
01101       }
01102     else
01103       {
01104         MPI_Request* requests=new MPI_Request[neighbourIds.size()];
01105         MPI_Request* req=requests;
01106         
01107         // setup sends
01108         for(std::vector<int>::iterator neighbour=neighbourIds.begin();
01109             neighbour!= neighbourIds.end(); ++neighbour){
01110           // Only send the information to the neighbouring processors
01111           MPI_Issend(buffer[0], position , MPI_PACKED, *neighbour, commTag_, comm_, req++);
01112         }
01113         
01114         //Test for received messages
01115         
01116         typedef typename std::vector<int>::size_type size_type;
01117         for(size_type received=0; received <neighbourIds.size(); ++received)
01118           {
01119             MPI_Status status;
01120             // probe for next message
01121             MPI_Probe(MPI_ANY_SOURCE, commTag_, comm_, &status);
01122             int remoteProc=status.MPI_SOURCE;
01123             int size;
01124             MPI_Get_count(&status, MPI_PACKED, &size);
01125             // receive message
01126             MPI_Recv(buffer[1], size, MPI_PACKED, remoteProc,
01127                      commTag_, comm_, &status);
01128               
01129             unpackCreateRemote(buffer[1], sourcePairs, destPairs, remoteProc, sourcePublish, 
01130                                destPublish, bufferSize, sendTwo);
01131           }
01132         // wait for completion of pending requests
01133         MPI_Status* statuses = new MPI_Status[neighbourIds.size()];
01134         
01135         if(MPI_ERR_IN_STATUS==MPI_Waitall(neighbourIds.size(), requests, statuses)){
01136           for(size_type i=0; i < neighbourIds.size(); ++i)
01137             if(statuses[i].MPI_ERROR!=MPI_SUCCESS){
01138               std::cerr<<rank<<": MPI_Error occurred while receiving message."<<std::endl;
01139               MPI_Abort(comm_, 999);
01140           }
01141         }
01142         delete[] requests;
01143         delete[] statuses;
01144       }
01145     
01146     MPI_Barrier(comm_);
01147     
01148     // delete allocated memory
01149     if(destPairs!=sourcePairs)
01150       delete[] destPairs;
01151     
01152     delete[] sourcePairs;
01153     delete[] buffer[0];
01154     delete[] buffer[1];
01155     delete[] buffer;
01156   }
01157   
01158   template<typename T>
01159   inline void RemoteIndices<T>::unpackIndices(RemoteIndexList& remote,
01160                                                   int remoteEntries,
01161                                                   PairType** local,
01162                                                   int localEntries,
01163                                                   char* p_in,
01164                                                   MPI_Datatype type,
01165                                                   int* position,
01166                                                   int bufferSize)
01167   {
01168     if(remoteEntries==0)
01169       return;
01170   
01171     PairType index;
01172     MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01173                type, comm_);
01174         
01175     int n_in=0, localIndex=0;
01176         
01177     //Check if we know the global index
01178     while(localIndex<localEntries){
01179       if(local[localIndex]->global()==index.global()){
01180         remote.push_back(RemoteIndex(index.local().attribute(), 
01181                                             local[localIndex]));
01182         localIndex++;
01183         // unpack next remote index
01184         if((++n_in) < remoteEntries){
01185           MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01186                      type, comm_);
01187         }else{
01188           // No more received indices
01189           break;
01190         }
01191         continue;
01192       }
01193 
01194       if (local[localIndex]->global()<index.global()){
01195         // compare with next entry in our list
01196         ++localIndex;
01197       }else{
01198         // We do not know the index, unpack next
01199         if((++n_in) < remoteEntries){
01200           MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01201                      type, comm_);
01202         }else
01203           // No more received indices
01204           break;
01205       }
01206     }
01207     
01208     // Unpack the other received indices without doing anything
01209     while(++n_in < remoteEntries)
01210       MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01211                type, comm_);
01212   }
01213   
01214     
01215   template<typename T>
01216   inline void RemoteIndices<T>::unpackIndices(RemoteIndexList& send,
01217                                                   RemoteIndexList& receive,
01218                                                   int remoteEntries,
01219                                                   PairType** localSource,
01220                                                   int localSourceEntries,
01221                                                   PairType** localDest,
01222                                                   int localDestEntries,
01223                                                   char* p_in,
01224                                                   MPI_Datatype type,
01225                                                   int* position,
01226                                                   int bufferSize)
01227   {             
01228     int n_in=0, sourceIndex=0, destIndex=0;
01229         
01230     //Check if we know the global index
01231     while(n_in<remoteEntries && (sourceIndex<localSourceEntries || destIndex<localDestEntries)){
01232       // Unpack next index
01233       PairType index;
01234       MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01235                  type, comm_);
01236       n_in++;
01237       
01238       // Advance until global index in localSource and localDest are >= than the one in the unpacked index
01239       while(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()<index.global())
01240         sourceIndex++;
01241       
01242       while(destIndex<localDestEntries && localDest[destIndex]->global()<index.global())
01243         destIndex++;
01244 
01245       // Add a remote index if we found the global index.
01246       if(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()==index.global())
01247           send.push_back(RemoteIndex(index.local().attribute(), 
01248                                               localSource[sourceIndex]));
01249 
01250       if(destIndex < localDestEntries && localDest[destIndex]->global() == index.global())
01251           receive.push_back(RemoteIndex(index.local().attribute(), 
01252                                         localDest[sourceIndex]));
01253     }
01254     
01255   }
01256     
01257   template<typename T>
01258   inline void RemoteIndices<T>::free()
01259   {
01260     typedef typename RemoteIndexMap::iterator Iterator;
01261     Iterator lend = remoteIndices_.end();
01262     for(Iterator lists=remoteIndices_.begin(); lists != lend; ++lists){
01263       if(lists->second.first==lists->second.second){
01264         // there is only one remote index list.
01265         delete lists->second.first;
01266       }else{
01267         delete lists->second.first;
01268         delete lists->second.second;
01269       }
01270     }
01271     remoteIndices_.clear();
01272     firstBuild=true;
01273   }
01274 
01275   template<typename T>
01276   inline int RemoteIndices<T>::neighbours() const
01277   {
01278     return remoteIndices_.size();
01279   }
01280   
01281   template<typename T>
01282   template<bool ignorePublic>
01283   inline void RemoteIndices<T>::rebuild()
01284   {
01285     // Test wether a rebuild is Needed.
01286     if(firstBuild || 
01287        ignorePublic!=publicIgnored || !
01288        isSynced()){
01289       free();
01290       
01291       buildRemote<ignorePublic>();
01292 
01293       sourceSeqNo_ = source_->seqNo();
01294       destSeqNo_ = target_->seqNo();
01295       firstBuild=false;
01296       publicIgnored=ignorePublic;
01297     }
01298     
01299     
01300   }
01301   
01302   template<typename T>
01303   inline bool RemoteIndices<T>::isSynced() const
01304   {
01305     return sourceSeqNo_==source_->seqNo() && destSeqNo_ ==target_->seqNo();
01306   }
01307 
01308   template<typename T>
01309   template<bool mode, bool send>
01310   RemoteIndexListModifier<T,mode> RemoteIndices<T>::getModifier(int process)
01311   {
01312 
01313     // The user are on their own now!
01314     // We assume they know what they are doing and just set the
01315     // remote indices to synced status.
01316     sourceSeqNo_ = source_->seqNo();
01317     destSeqNo_ = target_->seqNo();
01318 
01319     typename RemoteIndexMap::iterator found = remoteIndices_.find(process);
01320     
01321     if(found == remoteIndices_.end())
01322     {
01323       if(source_ != target_)
01324         remoteIndices_.insert(std::make_pair(process,
01325                                              std::make_pair(new RemoteIndexList(), 
01326                                                             new RemoteIndexList())));
01327       else{
01328         RemoteIndexList* rlist = new RemoteIndexList();
01329         remoteIndices_.insert(std::make_pair(process,
01330                                              std::make_pair(rlist, rlist)));
01331         
01332         found = remoteIndices_.find(process);
01333       }
01334     }
01335     
01336     firstBuild = false;
01337     
01338     if(send)
01339       return RemoteIndexListModifier<T,mode>(*source_, *(found->second.first));
01340     else
01341       return RemoteIndexListModifier<T,mode>(*target_, *(found->second.second));
01342   }
01343   
01344   template<typename T>
01345   inline typename std::map<int, std::pair<SLList<typename RemoteIndices<T>::RemoteIndex,
01346                                                  PoolAllocator<typename RemoteIndices<T>::RemoteIndex,1> >*,
01347                                           SLList<typename RemoteIndices<T>::RemoteIndex,
01348                                                  PoolAllocator<typename RemoteIndices<T>::RemoteIndex,1> >*> >::const_iterator
01349   RemoteIndices<T>::begin() const
01350   {
01351     return remoteIndices_.begin();
01352   }
01353   
01354   template<typename T>
01355   inline typename std::map<int, std::pair<SLList<typename RemoteIndices<T>::RemoteIndex,
01356                                                  PoolAllocator<typename RemoteIndices<T>::RemoteIndex,1> >*,
01357                                           SLList<typename RemoteIndices<T>::RemoteIndex,
01358                                                  PoolAllocator<typename RemoteIndices<T>::RemoteIndex,1> >*> >::const_iterator
01359   RemoteIndices<T>::end() const
01360   {
01361     return remoteIndices_.end();
01362   }
01363 
01364 
01365   template<typename T>
01366   bool RemoteIndices<T>::operator==(const RemoteIndices& ri)
01367   {
01368     if(neighbours()!=ri.neighbours())
01369       return false;
01370     
01371     typedef RemoteIndexList RList;
01372     typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
01373     
01374     const const_iterator rend = remoteIndices_.end();
01375 
01376     for(const_iterator rindex = remoteIndices_.begin(), rindex1=ri.remoteIndices_.begin(); rindex!=rend; ++rindex, ++rindex1){
01377       if(rindex->first != rindex1->first)
01378         return false;
01379       if(*(rindex->second.first) != *(rindex1->second.first))
01380         return false;
01381       if(*(rindex->second.second) != *(rindex1->second.second))
01382         return false;
01383     }
01384     return true;
01385   }
01386   
01387   template<class T, bool mode>
01388   RemoteIndexListModifier<T,mode>::RemoteIndexListModifier(const ParallelIndexSet& indexSet,
01389                                                            RemoteIndexList& rList)
01390     : rList_(&rList), indexSet_(&indexSet), glist_(new GlobalList()), iter_(rList.beginModify()), end_(rList.end()), first_(true)
01391   {
01392     if(MODIFYINDEXSET){
01393       assert(indexSet_);
01394       for(ConstIterator iter=iter_; iter != end_; ++iter)
01395         glist_->push_back(iter->localIndexPair().global());
01396       giter_ = glist_->beginModify();
01397     }
01398   }
01399 
01400   template<typename T, bool mode>
01401   RemoteIndexListModifier<T,mode>::RemoteIndexListModifier(const RemoteIndexListModifier<T,mode>& other)
01402     : rList_(other.rList_), indexSet_(other.indexSet_), 
01403       glist_(other.glist_), iter_(other.iter_), giter_(other.giter_), end_(other.end_), 
01404       first_(other.first_), last_(other.last_)
01405   {}
01406     
01407   template<typename T, bool mode>
01408   inline void RemoteIndexListModifier<T,mode>::repairLocalIndexPointers() throw(InvalidIndexSetState)
01409   { 
01410     if(MODIFYINDEXSET){
01411       // repair pointers to local index set.
01412 #ifdef DUNE_ISTL_WITH_CHECKING
01413       if(indexSet_->state()!=GROUND)
01414         DUNE_THROW(InvalidIndexSetState, "Index has to be in ground mode for repairing pointers to indices");
01415 #endif
01416       typedef typename ParallelIndexSet::const_iterator IndexIterator;
01417       typedef typename GlobalList::const_iterator GlobalIterator;
01418       typedef typename RemoteIndexList::iterator Iterator;
01419       GlobalIterator giter = glist_->begin();
01420       IndexIterator index = indexSet_->begin();
01421       
01422       for(Iterator iter=rList_->begin(); iter != end_; ++iter){
01423         while(index->global()<*giter){
01424           ++index;
01425 #ifdef DUNE_ISTL_WITH_CHECKING
01426           if(index == indexSet_.end())
01427             DUNE_THROW(ISTLError, "No such global index in set!");
01428 #endif
01429         }
01430 
01431 #ifdef DUNE_ISTL_WITH_CHECKING
01432           if(index->global() != *giter)
01433             DUNE_THROW(ISTLError, "No such global index in set!");
01434 #endif
01435           iter->localIndex_ = &(*index);
01436       }
01437     }
01438   }
01439   
01440   template<typename T, bool mode>
01441   inline void RemoteIndexListModifier<T,mode>::insert(const RemoteIndex& index) throw(InvalidPosition)
01442   {
01443     dune_static_assert(!mode,"Not allowed if the mode indicates that new indices"
01444                        "might be added to the underlying index set. Use "
01445                        "insert(const RemoteIndex&, const GlobalIndex&) instead");
01446     
01447 #ifdef DUNE_ISTL_WITH_CHECKING
01448     if(!first_ && index.localIndexPair().global()<last_)
01449       DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
01450 #endif
01451     // Move to the correct position
01452     while(iter_ != end_ && iter_->localIndexPair().global() < index.localIndexPair().global()){
01453       ++iter_;
01454     }
01455     
01456     // No duplicate entries allowed
01457     assert(iter_==end_ || iter_->localIndexPair().global() != index.localIndexPair().global());
01458     iter_.insert(index);
01459     last_ = index.localIndexPair().global();
01460     first_ = false;
01461   }
01462   
01463   template<typename T, bool mode>
01464   inline void RemoteIndexListModifier<T,mode>::insert(const RemoteIndex& index, const GlobalIndex& global) throw(InvalidPosition)
01465   {
01466      dune_static_assert(mode,"Not allowed if the mode indicates that no new indices"
01467                        "might be added to the underlying index set. Use "
01468                        "insert(const RemoteIndex&) instead");
01469 #ifdef DUNE_ISTL_WITH_CHECKING
01470     if(!first_ && global<last_)
01471       DUNE_THROW(InvalidPosition, "Modification of remote indices have to occur with ascending global index.");
01472 #endif
01473     // Move to the correct position
01474     while(iter_ != end_ && *giter_ < global){
01475       ++giter_;
01476       ++iter_;
01477     }
01478     
01479     // No duplicate entries allowed
01480     assert(iter_->localIndexPair().global() != global);
01481     iter_.insert(index);
01482     giter_.insert(global);
01483     
01484     last_ = global;
01485     first_ = false;
01486   }
01487   
01488   template<typename T, bool mode>
01489   bool RemoteIndexListModifier<T,mode>::remove(const GlobalIndex& global) throw(InvalidPosition)
01490   {
01491 #ifdef DUNE_ISTL_WITH_CHECKING
01492     if(!first_ && global<last_)
01493       DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
01494 #endif
01495 
01496     bool found= false;
01497     
01498     if(MODIFYINDEXSET){
01499     // Move to the correct position
01500       while(iter_!=end_ && *giter_< global){
01501         ++giter_;
01502         ++iter_;
01503       }
01504       if(*giter_ == global){
01505         giter_.remove();
01506         iter_.remove();
01507         found=true;
01508       }
01509     }else{
01510       while(iter_!=end_ && iter_->localIndexPair().global() < global)
01511         ++iter_;
01512         
01513       if(iter_->localIndexPair().global()==global){
01514         iter_.remove();
01515         found = true;
01516       }
01517     }
01518     
01519     last_ = global;
01520     first_ = false;
01521     return found;
01522   }
01523 
01524   template<typename T>
01525   template<bool send>
01526   inline CollectiveIterator<T> RemoteIndices<T>::iterator() const
01527   {
01528     return CollectiveIterator<T>(remoteIndices_, send);
01529   }
01530 
01531   template<typename T>
01532   inline MPI_Comm RemoteIndices<T>::communicator() const
01533   {
01534     return comm_;
01535     
01536   }
01537   
01538   template<typename T>
01539   CollectiveIterator<T>::CollectiveIterator(const RemoteIndexMap& pmap, bool send)
01540   {
01541     typedef typename RemoteIndexMap::const_iterator const_iterator;
01542     typedef typename RemoteIndexMap::iterator iterator;
01543     
01544     const const_iterator end=pmap.end();
01545     for(const_iterator process=pmap.begin(); process != end; ++process){
01546       const RemoteIndexList* list = send? process->second.first : process->second.second;
01547       typedef typename RemoteIndexList::const_iterator iterator;
01548       map_.insert(std::make_pair(process->first, 
01549                                  std::pair<iterator, const iterator>(list->begin(), list->end())));
01550     }
01551   }
01552 
01553   template<typename T>
01554   inline void CollectiveIterator<T>::advance(const GlobalIndex& index)
01555   {
01556     typedef typename Map::iterator iterator;
01557     typedef typename Map::const_iterator const_iterator;
01558     const const_iterator end = map_.end();
01559     
01560     for(iterator iter = map_.begin(); iter != end;){
01561       // Step the iterator until we are >= index
01562       typename RemoteIndexList::const_iterator current = iter->second.first;
01563       typename RemoteIndexList::const_iterator rend = iter->second.second;
01564       RemoteIndex remoteIndex;
01565       if(current != rend)
01566         remoteIndex = *current;
01567       
01568       while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
01569         ++(iter->second.first);
01570       
01571       // erase from the map if there are no more entries.
01572       if(iter->second.first == iter->second.second)
01573         map_.erase(iter++);
01574       else{
01575         ++iter;
01576       }
01577     }
01578     index_=index;
01579   }
01580   
01581   template<typename T>
01582   inline bool CollectiveIterator<T>::empty()
01583   {
01584     return map_.empty();
01585   }
01586     
01587   template<typename T>
01588   inline typename CollectiveIterator<T>::iterator
01589   CollectiveIterator<T>::begin()
01590   {
01591     return iterator(map_.begin(), map_.end(), index_);
01592   }
01593    
01594   template<typename T>
01595   inline typename CollectiveIterator<T>::iterator
01596   CollectiveIterator<T>::end()
01597   {
01598     return iterator(map_.end(), map_.end(), index_);
01599   }
01600   
01601   template<typename TG, typename TA>
01602   inline std::ostream& operator<<(std::ostream& os, const RemoteIndex<TG,TA>& index)
01603   {
01604     os<<"[global="<<index.localIndexPair().global()<<",attribute="<<index.attribute()<<"]";
01605     return os;
01606   }
01607   
01608   template<typename T>
01609   inline std::ostream& operator<<(std::ostream& os, const RemoteIndices<T>& indices)
01610   {
01611     int rank;
01612     MPI_Comm_rank(indices.comm_, &rank);
01613 
01614     typedef typename RemoteIndices<T>::RemoteIndexList RList;
01615     typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
01616 
01617     const const_iterator rend = indices.remoteIndices_.end();
01618 
01619     for(const_iterator rindex = indices.remoteIndices_.begin(); rindex!=rend; ++rindex){
01620       os<<rank<<": Prozess "<<rindex->first<<":";
01621       
01622       if(!rindex->second.first->empty()){
01623         os<<" send:";
01624 
01625         const typename RList::const_iterator send= rindex->second.first->end();
01626       
01627         for(typename RList::const_iterator index = rindex->second.first->begin(); 
01628             index != send; ++index)
01629           os<<*index<<" ";
01630         os<<std::endl;
01631       }
01632       if(!rindex->second.second->empty()){
01633         os<<rank<<": Prozess "<<rindex->first<<": "<<"receive: ";
01634         
01635         const typename RList::const_iterator rend= rindex->second.second->end();
01636         
01637         for(typename RList::const_iterator index = rindex->second.second->begin(); 
01638             index != rend; ++index)
01639           os<<*index<<" ";
01640       }
01641       os<<std::endl<<std::flush;
01642     }
01643     return os;
01644   }
01646 }
01647 
01648 #endif
01649 #endif

Generated on Thu Apr 2 10:40:15 2009 for dune-istl by  doxygen 1.5.6