remoteindices.hh

Go to the documentation of this file.
00001 // $Id: remoteindices.hh 1093 2009-10-02 15:11:15Z 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<set>
00013 #include<utility>
00014 #include<iostream>
00015 #include<algorithm>
00016 #include<iterator>
00017 #if HAVE_MPI
00018 #include"mpitraits.hh"
00019 #include"mpi.h"
00020 
00021 namespace Dune{
00032 
00033   template<typename TG, typename TA>
00034   class MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >
00035   {
00036   public:
00037     inline static MPI_Datatype getType();
00038   private:
00039     static MPI_Datatype type;
00040   };
00041   
00042   
00043   template<typename T, typename A>
00044   class RemoteIndices;
00045   
00046   template<typename T1, typename T2>
00047   class RemoteIndex;
00048   
00049   template<typename T>
00050   class IndicesSyncer;
00051   
00052   template<typename T1, typename T2>
00053   std::ostream& operator<<(std::ostream& os, const RemoteIndex<T1,T2>& index);
00054   
00055 
00056   template<typename T, typename A, bool mode>
00057   class RemoteIndexListModifier;
00058 
00062   template<typename T1, typename T2>
00063   class RemoteIndex
00064   {
00065     template<typename T>
00066     friend class IndicesSyncer;
00067     
00068     template<typename T, typename A, typename A1>
00069     friend void repairLocalIndexPointers(std::map<int,SLList<typename T::GlobalIndex,A> >&, RemoteIndices<T,A1>&,
00070                                          const T&); 
00071 
00072     template<typename T, typename A, bool mode>
00073     friend class RemoteIndexListModifier;
00074 
00075   public:    
00080     typedef T1 GlobalIndex;
00089     typedef T2 Attribute;
00090 
00094     typedef IndexPair<GlobalIndex,ParallelLocalIndex<Attribute> >
00095     PairType;
00096     
00101     const Attribute attribute() const;
00102 
00108     const PairType& localIndexPair() const;
00109 
00113     RemoteIndex();
00114     
00115     
00121     RemoteIndex(const T2& attribute, 
00122                 const PairType* local);
00123 
00124     
00130     RemoteIndex(const T2& attribute);
00131 
00132     bool operator==(const RemoteIndex& ri) const;
00133     
00134     bool operator!=(const RemoteIndex& ri) const;
00135   private:
00137     const PairType* localIndex_;
00138 
00140     char attribute_;
00141   };
00142   
00143   template<class T, class A>
00144   std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices);
00145   
00146   class InterfaceBuilder;
00147 
00148   template<class T, class A>
00149   class CollectiveIterator;
00150 
00151   template<class T>
00152     class IndicesSyncer;
00153   
00154   // forward declaration needed for friend declaration.
00155   template<typename T1, typename T2>
00156   class OwnerOverlapCopyCommunication;
00157   
00158 
00175   template<class T, class A=std::allocator<RemoteIndex<typename T::GlobalIndex,
00176                                                        typename T::LocalIndex::Attribute> > >
00177   class RemoteIndices
00178   {
00179     friend class InterfaceBuilder;
00180     friend class IndicesSyncer<T>;
00181     template<typename T1,typename A1,typename A2>
00182     friend void repairLocalIndexPointers(std::map<int,SLList<typename T1::GlobalIndex,A1> >&, 
00183                                          RemoteIndices<T1,A2>&,
00184                                          const T1&);
00185     
00186     template<class G, class T1, class T2>
00187     friend void fillIndexSetHoles(const G& graph, Dune::OwnerOverlapCopyCommunication<T1,T2>& oocomm);
00188     friend std::ostream& operator<<<>(std::ostream&, const RemoteIndices<T>&);
00189     
00190   public:
00191    
00195     typedef T ParallelIndexSet;
00196  
00199     typedef CollectiveIterator<T,A> CollectiveIteratorT;
00200     
00204     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00205     
00206     
00210     typedef typename ParallelIndexSet::LocalIndex LocalIndex;
00211 
00215     typedef typename LocalIndex::Attribute Attribute;
00216  
00220     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00221 
00222     
00226     typedef typename A::template rebind<RemoteIndex>::other Allocator;
00227 
00229     typedef Dune::SLList<RemoteIndex,Allocator>
00230     RemoteIndexList;
00231     
00233     typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
00234     RemoteIndexMap;    
00235     
00236     typedef typename RemoteIndexMap::const_iterator const_iterator;
00237     
00251     inline RemoteIndices(const ParallelIndexSet& source, const ParallelIndexSet& destination, 
00252                          const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>());
00253 
00254     RemoteIndices();
00255     
00272     void setIndexSets(const ParallelIndexSet& source, const ParallelIndexSet& destination, 
00273                       const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>());
00274 
00275     template<typename C>
00276     void setNeighbours(const C& neighbours)
00277     {
00278       typedef typename C::const_iterator Iter;
00279       neighbourIds.clear();
00280       neighbourIds.insert(neighbours.begin(), neighbours.end());
00281         
00282     }
00283     
00287     ~RemoteIndices();
00288     
00298     template<bool ignorePublic>
00299     void rebuild();
00300 
00301     bool operator==(const RemoteIndices& ri);
00302     
00310     inline bool isSynced() const;
00311 
00315     inline MPI_Comm communicator() const;
00316 
00331     template<bool mode, bool send>
00332     inline RemoteIndexListModifier<T,A,mode> getModifier(int process);
00333     
00338     inline const_iterator begin() const;
00339     
00344     inline const_iterator end() const;
00345     
00349     template<bool send>
00350     inline CollectiveIteratorT iterator() const;
00351 
00355     inline void free();
00356     
00361     inline int neighbours() const;
00362 
00364     inline const ParallelIndexSet& sourceIndexSet() const;
00365     
00367     inline const ParallelIndexSet& destinationIndexSet() const;
00368 
00369   private:    
00371     const ParallelIndexSet* source_;
00372 
00374     const ParallelIndexSet* target_;
00375 
00377     MPI_Comm comm_;
00378 
00381     std::set<int> neighbourIds;
00382     
00384     const static int commTag_=333;
00385     
00390     int sourceSeqNo_;
00391     
00396     int destSeqNo_;
00397 
00401     bool publicIgnored;
00402 
00406     bool firstBuild;
00407         
00409     typedef IndexPair<GlobalIndex, LocalIndex> 
00410     PairType;
00411     
00418     RemoteIndexMap remoteIndices_;
00419     
00426     template<bool ignorePublic>
00427     inline void buildRemote();
00428 
00434     inline int noPublic(const ParallelIndexSet& indexSet);
00435     
00447     template<bool ignorePublic>
00448     inline void packEntries(PairType** myPairs, const ParallelIndexSet& indexSet,
00449                             char* p_out, MPI_Datatype type, int bufferSize, 
00450                             int* position, int n);
00451     
00465     inline void unpackIndices(RemoteIndexList& remote, int remoteEntries,
00466                               PairType** local, int localEntries, char* p_in, 
00467                               MPI_Datatype type, int* positon, int bufferSize);
00468 
00469     inline void unpackIndices(RemoteIndexList& send, RemoteIndexList& receive,
00470                               int remoteEntries, PairType** localSource,
00471                               int localSourceEntries, PairType** localDest,
00472                               int localDestEntries, char* p_in,
00473                               MPI_Datatype type, int* position, int bufferSize);
00474 
00475     void unpackCreateRemote(char* p_in, PairType** sourcePairs, PairType** DestPairs, 
00476                             int remoteProc,  int sourcePublish, int destPublish, 
00477                             int bufferSize, bool sendTwo);
00478   };
00479   
00497   template<class T, class A, bool mode>
00498   class RemoteIndexListModifier
00499   {
00500 
00501     template<typename T1, typename A1>
00502     friend class RemoteIndices;
00503   
00504   public:
00505     class InvalidPosition : public ISTLError
00506     {};
00507 
00508     enum {
00517       MODIFYINDEXSET=mode
00518     };
00519   
00523     typedef T ParallelIndexSet;
00524 
00528     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00529 
00533     typedef typename ParallelIndexSet::LocalIndex LocalIndex;
00534     
00538     typedef typename LocalIndex::Attribute Attribute;
00539     
00543     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00544     
00548     typedef A Allocator;
00549 
00551     typedef Dune::SLList<RemoteIndex,Allocator>
00552     RemoteIndexList;
00553 
00557     typedef SLListModifyIterator<RemoteIndex,Allocator> ModifyIterator;
00558     
00562     typedef typename RemoteIndexList::const_iterator ConstIterator;
00563     
00577     void insert(const RemoteIndex& index) throw(InvalidPosition);
00578     
00579     
00594     void insert(const RemoteIndex& index, const GlobalIndex& global) throw(InvalidPosition);
00595 
00603     bool remove(const GlobalIndex& global) throw(InvalidPosition);
00604 
00617     void repairLocalIndexPointers() throw(InvalidIndexSetState);
00618 
00619 
00620     RemoteIndexListModifier(const RemoteIndexListModifier&);
00621     
00626     RemoteIndexListModifier()
00627       : glist_()
00628     {}
00629     
00631     ~RemoteIndexListModifier()
00632     {
00633       if(glist_)
00634         delete glist_;
00635     }
00636     
00637   private:
00638     
00644     RemoteIndexListModifier(const ParallelIndexSet& indexSet,
00645                             RemoteIndexList& rList);
00646 
00647     typedef SLList<GlobalIndex,Allocator> GlobalList;
00648     typedef typename GlobalList::ModifyIterator GlobalModifyIterator;
00649     RemoteIndices<ParallelIndexSet>* remoteIndices_;
00650     RemoteIndexList* rList_;
00651     const ParallelIndexSet* indexSet_;
00652     GlobalList* glist_;
00653     ModifyIterator iter_;
00654     GlobalModifyIterator giter_;
00655     ConstIterator end_;
00656     bool first_;
00657     GlobalIndex last_;
00658   };
00659   
00664   template<class T, class A>
00665   class CollectiveIterator
00666   {
00667     
00671     typedef T ParallelIndexSet;
00672 
00676     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00677 
00681     typedef typename ParallelIndexSet::LocalIndex LocalIndex;
00682     
00686     typedef typename LocalIndex::Attribute Attribute;
00687 
00689     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00690     
00692     typedef typename A::template rebind<RemoteIndex>::other Allocator;
00693 
00695     typedef Dune::SLList<RemoteIndex,Allocator> RemoteIndexList;    
00696 
00698     typedef std::map<int,std::pair<typename RemoteIndexList::const_iterator,
00699                                    const typename RemoteIndexList::const_iterator> >
00700     Map;
00701 
00702   public:
00703     
00705     typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
00706     RemoteIndexMap;
00707     
00713     inline CollectiveIterator(const RemoteIndexMap& map_, bool send);
00714     
00723     inline void advance(const GlobalIndex& global);
00724     
00728     inline bool empty();
00729     
00736     class iterator
00737     {
00738     public:
00739       typedef typename Map::iterator RealIterator;
00740       typedef typename Map::iterator ConstRealIterator;
00741       
00742       
00744       iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex& index)
00745         : iter_(iter), end_(end), index_(index)
00746       {
00747         // Move to the first valid entry
00748         while(iter_!=end_ && iter_->second.first->localIndexPair().global()!=index_)
00749           ++iter_;
00750       }
00751       
00753       iterator(const iterator& other)
00754         : iter_(other.iter_), end_(other.end_), index_(other.index_)
00755       { }
00756          
00758       iterator& operator++()
00759       { 
00760         ++iter_;
00761         // If entry is not valid move on
00762         while(iter_!=end_ && iter_->second.first->localIndexPair().global()!=index_)
00763           ++iter_;
00764                   
00765         return *this;
00766       }
00767       
00769       const RemoteIndex& operator*()const
00770       {
00771         return *(iter_->second.first);
00772       }
00773       
00775       int process() const
00776       {
00777         return iter_->first;
00778       }
00779       
00781       const RemoteIndex* operator->()const
00782       {
00783         return iter_->second.first.operator->();
00784       }
00785       
00787       bool operator==(const iterator& other)
00788       {
00789         return other.iter_==iter_;
00790       }
00791 
00793       bool operator!=(const iterator& other)
00794       {
00795         return other.iter_!=iter_;
00796       }
00797 
00798     private:
00799       iterator();
00800       
00801       RealIterator iter_;
00802       RealIterator end_;
00803       GlobalIndex index_;
00804     };
00805     
00806     iterator begin();
00807     
00808     iterator end();
00809     
00810   private:
00811     
00812     Map map_;
00813     GlobalIndex index_;
00814   };
00815     
00816   template<typename TG, typename TA>
00817   MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::getType()
00818   {
00819     if(type==MPI_DATATYPE_NULL){
00820       int length[4];
00821       MPI_Aint disp[4];
00822       MPI_Datatype types[4] = {MPI_LB, MPITraits<TG>::getType(), 
00823                                MPITraits<ParallelLocalIndex<TA> >::getType(), MPI_UB};
00824       IndexPair<TG,ParallelLocalIndex<TA> > rep[2];
00825       length[0]=length[1]=length[2]=length[3]=1;
00826       MPI_Address(rep, disp); // lower bound of the datatype
00827       MPI_Address(&(rep[0].global_), disp+1);
00828       MPI_Address(&(rep[0].local_), disp+2);
00829       MPI_Address(rep+1, disp+3); // upper bound of the datatype
00830       for(int i=3; i >= 0; --i)
00831         disp[i] -= disp[0];
00832       MPI_Type_struct(4, length, disp, types, &type);
00833       MPI_Type_commit(&type);
00834     }
00835     return type;
00836   }
00837 
00838   template<typename TG, typename TA>
00839   MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::type=MPI_DATATYPE_NULL;
00840 
00841   template<typename T1, typename T2>
00842   RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute, const PairType* local)
00843     : localIndex_(local), attribute_(attribute)
00844   {}
00845 
00846   template<typename T1, typename T2>
00847   RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute)
00848     : localIndex_(0), attribute_(attribute)
00849   {}
00850 
00851   template<typename T1, typename T2>
00852   RemoteIndex<T1,T2>::RemoteIndex()
00853     : localIndex_(0), attribute_()
00854   {}
00855   template<typename T1, typename T2>
00856   inline bool RemoteIndex<T1,T2>::operator==(const RemoteIndex& ri) const
00857   {
00858     return localIndex_==ri.localIndex_ && attribute_==ri.attribute;
00859   }
00860   
00861   template<typename T1, typename T2>
00862   inline bool RemoteIndex<T1,T2>::operator!=(const RemoteIndex& ri) const
00863   {
00864     return localIndex_!=ri.localIndex_ || attribute_!=ri.attribute_;
00865   }
00866   
00867   template<typename T1, typename T2>
00868   inline const T2 RemoteIndex<T1,T2>::attribute() const
00869   {
00870     return T2(attribute_);
00871   }
00872 
00873   template<typename T1, typename T2>
00874   inline const IndexPair<T1,ParallelLocalIndex<T2> >& RemoteIndex<T1,T2>::localIndexPair() const
00875   {
00876     return *localIndex_;
00877   }
00878   
00879   template<typename T, typename A>
00880   inline RemoteIndices<T,A>::RemoteIndices(const ParallelIndexSet& source,
00881                                          const ParallelIndexSet& destination,
00882                                          const MPI_Comm& comm,
00883                                          const std::vector<int>& neighbours)
00884     : source_(&source), target_(&destination), comm_(comm),
00885       sourceSeqNo_(-1), destSeqNo_(-1), publicIgnored(false), firstBuild(true)
00886   {
00887     setNeighbours(neighbours);
00888   }
00889 
00890   template<typename T, typename A>
00891   RemoteIndices<T,A>::RemoteIndices()
00892     :source_(0), target_(0), sourceSeqNo_(-1), 
00893      destSeqNo_(-1), publicIgnored(false), firstBuild(true)
00894   {}
00895   
00896   template<class T, typename A>
00897   void RemoteIndices<T,A>::setIndexSets(const ParallelIndexSet& source,
00898                                       const ParallelIndexSet& destination,
00899                                       const MPI_Comm& comm,
00900                                       const std::vector<int>& neighbours)
00901   {
00902     free();
00903     source_ = &source;
00904     target_ = &destination;
00905     comm_ = comm;
00906     firstBuild = true;
00907     setNeighbours(neighbours);
00908   }
00909 
00910   template<typename T, typename A>
00911   const typename RemoteIndices<T,A>::ParallelIndexSet& 
00912   RemoteIndices<T,A>::sourceIndexSet() const
00913   {
00914     return *source_;
00915   }
00916   
00917     
00918   template<typename T, typename A>
00919   const typename RemoteIndices<T,A>::ParallelIndexSet& 
00920   RemoteIndices<T,A>::destinationIndexSet() const
00921   {
00922     return *target_;
00923   }
00924   
00925   
00926   template<typename T, typename A>
00927   RemoteIndices<T,A>::~RemoteIndices()
00928   {
00929     free();
00930   }
00931 
00932   template<typename T, typename A>
00933   template<bool ignorePublic>
00934   inline void RemoteIndices<T,A>::packEntries(IndexPair<GlobalIndex,LocalIndex>** pairs,
00935                                                 const ParallelIndexSet& indexSet,
00936                                                 char* p_out, MPI_Datatype type, 
00937                                                 int bufferSize,
00938                                                 int *position, int n)
00939   {
00940     // fill with own indices
00941     typedef typename ParallelIndexSet::const_iterator const_iterator;
00942     typedef IndexPair<GlobalIndex,LocalIndex> PairType;
00943     const const_iterator end = indexSet.end();
00944 
00945     //Now pack the source indices
00946     int i=0;
00947     for(const_iterator index = indexSet.begin(); index != end; ++index)
00948       if(ignorePublic || index->local().isPublic()){
00949         
00950         MPI_Pack(const_cast<PairType*>(&(*index)), 1, 
00951                  type,
00952                  p_out, bufferSize, position, comm_);
00953         pairs[i++] = const_cast<PairType*>(&(*index));
00954       
00955       }
00956     assert(i==n);
00957   }
00958   
00959   template<typename T, typename A>
00960   inline int RemoteIndices<T,A>::noPublic(const ParallelIndexSet& indexSet)
00961   {
00962     typedef typename ParallelIndexSet::const_iterator const_iterator;
00963     
00964     int noPublic=0;
00965     
00966     const const_iterator end=indexSet.end();
00967     for(const_iterator index=indexSet.begin(); index!=end; ++index)
00968       if(index->local().isPublic())
00969         noPublic++;
00970     
00971     return noPublic;
00972     
00973   }
00974     
00975 
00976   template<typename T, typename A>
00977   inline void RemoteIndices<T,A>::unpackCreateRemote(char* p_in, PairType** sourcePairs,
00978                                                    PairType** destPairs, int remoteProc, 
00979                                                    int sourcePublish, int destPublish,
00980                                                    int bufferSize, bool sendTwo)
00981   {
00982     
00983       // unpack the number of indices we received
00984       int noRemoteSource=-1, noRemoteDest=-1;
00985       char twoIndexSets=0;
00986       int position=0;
00987       // Did we receive two index sets?
00988       MPI_Unpack(p_in, bufferSize, &position, &twoIndexSets, 1, MPI_CHAR, comm_);
00989       // The number of source indices received
00990       MPI_Unpack(p_in, bufferSize, &position, &noRemoteSource, 1, MPI_INT, comm_);
00991       // The number of destination indices received
00992       MPI_Unpack(p_in, bufferSize, &position, &noRemoteDest, 1, MPI_INT, comm_);      
00993 
00994 
00995       // Indices for which we receive
00996       RemoteIndexList* receive= new RemoteIndexList();
00997       // Indices for which we send
00998       RemoteIndexList* send=0;
00999 
01000       MPI_Datatype type= MPITraits<PairType>::getType();
01001       
01002       if(!twoIndexSets){
01003         if(sendTwo){
01004           send = new RemoteIndexList();
01005           // Create both remote index sets simultaneously
01006           unpackIndices(*send, *receive, noRemoteSource, sourcePairs, sourcePublish,
01007                         destPairs, destPublish, p_in, type, &position, bufferSize);
01008         }else{
01009           // we only need one list
01010           unpackIndices(*receive, noRemoteSource, sourcePairs, sourcePublish,
01011                         p_in, type, &position, bufferSize);
01012           send=receive;
01013         }
01014       }else{
01015         
01016         int oldPos=position;
01017         // Two index sets received
01018         unpackIndices(*receive, noRemoteSource, destPairs, destPublish,
01019                       p_in, type, &position, bufferSize);
01020         if(!sendTwo)
01021           //unpack source entries again as destination entries
01022           position=oldPos;
01023         
01024         send = new RemoteIndexList();
01025         unpackIndices(*send, noRemoteDest, sourcePairs, sourcePublish, 
01026                       p_in, type, &position, bufferSize);
01027       }
01028       
01029       if(receive->empty() && send->empty()){
01030         if(send==receive){
01031           delete send;
01032         }else{
01033           delete send;
01034           delete receive;
01035         }
01036       }else{
01037         remoteIndices_.insert(std::make_pair(remoteProc, 
01038                                              std::make_pair(send,receive)));
01039       }
01040   }
01041   
01042 
01043   template<typename T, typename A>
01044   template<bool ignorePublic>
01045   inline void RemoteIndices<T,A>::buildRemote()
01046   {
01047     // Processor configuration
01048     int rank, procs;
01049     MPI_Comm_rank(comm_, &rank);
01050     MPI_Comm_size(comm_, &procs);
01051     
01052     // number of local indices to publish
01053     // The indices of the destination will be send.
01054     int sourcePublish, destPublish;
01055     
01056     // Do we need to send two index sets?
01057     char sendTwo = (source_ != target_);
01058     
01059     if(procs==0 && !sendTwo)
01060       // Nothing to communicate
01061       return;   
01062 
01063     sourcePublish = (ignorePublic)? source_->size() : noPublic(*source_);
01064     
01065     if(sendTwo)
01066       destPublish = (ignorePublic)? target_->size() : noPublic(*target_);
01067     else  
01068       // we only need to send one set of indices
01069       destPublish = 0;
01070         
01071     int maxPublish, publish=sourcePublish+destPublish;
01072 
01073     // Calucate maximum number of indices send
01074     MPI_Allreduce(&publish, &maxPublish, 1, MPI_INT, MPI_MAX, comm_);
01075 
01076     // allocate buffers
01077     typedef IndexPair<GlobalIndex,LocalIndex> PairType;
01078     
01079     PairType** destPairs;
01080     PairType** sourcePairs = new PairType*[sourcePublish>0?sourcePublish:1];
01081     
01082     if(sendTwo)
01083       destPairs = new PairType*[destPublish>0?destPublish:1];
01084     else
01085       destPairs=sourcePairs;
01086     
01087     char** buffer = new char*[2];
01088     int bufferSize;    
01089     int position=0;
01090     int intSize;
01091     int charSize;
01092 
01093     // calculate buffer size
01094     MPI_Datatype type = MPITraits<PairType>::getType();
01095     
01096     MPI_Pack_size(maxPublish, type, comm_,
01097                   &bufferSize);
01098     MPI_Pack_size(1, MPI_INT, comm_,
01099                   &intSize);
01100     MPI_Pack_size(1, MPI_CHAR, comm_,
01101                   &charSize);
01102     // Our message will contain the following:
01103     // a bool wether two index sets where sent
01104     // the size of the source and the dest indexset, 
01105     // then the source and destination indices
01106     bufferSize += 2 * intSize + charSize;
01107     
01108     if(bufferSize<=0) bufferSize=1;
01109     
01110     buffer[0] = new char[bufferSize];
01111     buffer[1] = new char[bufferSize];
01112     
01113     
01114     // pack entries into buffer[0], p_out below!
01115     MPI_Pack(&sendTwo, 1, MPI_CHAR, buffer[0], bufferSize, &position,
01116                  comm_);
01117         
01118     // The number of indices we send for each index set
01119     MPI_Pack(&sourcePublish, 1, MPI_INT, buffer[0], bufferSize, &position, 
01120              comm_);
01121     MPI_Pack(&destPublish, 1, MPI_INT, buffer[0], bufferSize, &position,
01122              comm_);
01123         
01124     // Now pack the source indices and setup the destination pairs
01125     packEntries<ignorePublic>(sourcePairs, *source_, buffer[0], type, 
01126                               bufferSize, &position, sourcePublish);
01127     // If necessary send the dest indices and setup the source pairs
01128     if(sendTwo)
01129       packEntries<ignorePublic>(destPairs, *target_, buffer[0], type,
01130                                 bufferSize, &position, destPublish);
01131 
01132 
01133     // Update remote indices for ourself
01134     if(sendTwo)
01135       unpackCreateRemote(buffer[0], sourcePairs, destPairs, rank, sourcePublish, 
01136                          destPublish, bufferSize, sendTwo);
01137 
01138     neighbourIds.erase(rank);
01139 
01140     if(neighbourIds.size()==0)
01141       {
01142         // send mesages in ring
01143         for(int proc=1; proc<procs; proc++){
01144           // pointers to the current input and output buffers
01145           char* p_out = buffer[1-(proc%2)];
01146           char* p_in = buffer[proc%2];
01147           
01148           MPI_Status status;
01149           if(rank%2==0){
01150             MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
01151                       commTag_, comm_);
01152             MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
01153                      commTag_, comm_, &status);
01154           }else{
01155             MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
01156                      commTag_, comm_, &status);
01157             MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
01158                       commTag_, comm_);
01159           }
01160 
01161 
01162           // The process these indices are from
01163           int remoteProc = (rank+procs-proc)%procs;
01164           
01165           unpackCreateRemote(p_in, sourcePairs, destPairs, remoteProc, sourcePublish, 
01166                              destPublish, bufferSize, sendTwo);
01167           
01168         }
01169         
01170       }
01171     else
01172       {
01173         MPI_Request* requests=new MPI_Request[neighbourIds.size()];
01174         MPI_Request* req=requests;
01175         
01176         typedef typename std::set<int>::size_type size_type;
01177         size_type noNeighbours=neighbourIds.size();
01178 
01179         // setup sends  
01180         for(std::set<int>::iterator neighbour=neighbourIds.begin();
01181             neighbour!= neighbourIds.end(); ++neighbour){
01182             // Only send the information to the neighbouring processors
01183             MPI_Issend(buffer[0], position , MPI_PACKED, *neighbour, commTag_, comm_, req++);
01184         }
01185         
01186         //Test for received messages
01187         
01188         for(size_type received=0; received <noNeighbours; ++received)
01189           {
01190             MPI_Status status;
01191             // probe for next message
01192             MPI_Probe(MPI_ANY_SOURCE, commTag_, comm_, &status);
01193             int remoteProc=status.MPI_SOURCE;
01194             int size;
01195             MPI_Get_count(&status, MPI_PACKED, &size);
01196             // receive message
01197             MPI_Recv(buffer[1], size, MPI_PACKED, remoteProc,
01198                      commTag_, comm_, &status);
01199               
01200             unpackCreateRemote(buffer[1], sourcePairs, destPairs, remoteProc, sourcePublish, 
01201                                destPublish, bufferSize, sendTwo);
01202           }
01203         // wait for completion of pending requests
01204         MPI_Status* statuses = new MPI_Status[neighbourIds.size()];
01205         
01206         if(MPI_ERR_IN_STATUS==MPI_Waitall(neighbourIds.size(), requests, statuses)){
01207           for(size_type i=0; i < neighbourIds.size(); ++i)
01208             if(statuses[i].MPI_ERROR!=MPI_SUCCESS){
01209               std::cerr<<rank<<": MPI_Error occurred while receiving message."<<std::endl;
01210               MPI_Abort(comm_, 999);
01211           }
01212         }
01213         delete[] requests;
01214         delete[] statuses;
01215       }
01216     
01217     
01218     // delete allocated memory
01219     if(destPairs!=sourcePairs)
01220       delete[] destPairs;
01221     
01222     delete[] sourcePairs;
01223     delete[] buffer[0];
01224     delete[] buffer[1];
01225     delete[] buffer;
01226   }
01227   
01228   template<typename T, typename A>
01229   inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& remote,
01230                                                   int remoteEntries,
01231                                                   PairType** local,
01232                                                   int localEntries,
01233                                                   char* p_in,
01234                                                   MPI_Datatype type,
01235                                                   int* position,
01236                                                   int bufferSize)
01237   {
01238     if(remoteEntries==0)
01239       return;
01240   
01241     PairType index(-1);
01242     MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01243                type, comm_);
01244         
01245     int n_in=0, localIndex=0;
01246         
01247     //Check if we know the global index
01248     while(localIndex<localEntries){
01249       if(local[localIndex]->global()==index.global()){
01250         remote.push_back(RemoteIndex(index.local().attribute(), 
01251                                             local[localIndex]));
01252         localIndex++;
01253         // unpack next remote index
01254         if((++n_in) < remoteEntries){
01255           MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01256                      type, comm_);
01257         }else{
01258           // No more received indices
01259           break;
01260         }
01261         continue;
01262       }
01263       
01264       if (local[localIndex]->global()<index.global()){
01265         // compare with next entry in our list
01266         ++localIndex;
01267       }else{
01268         // We do not know the index, unpack next
01269         if((++n_in) < remoteEntries){
01270           MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01271                      type, comm_);
01272         }else
01273           // No more received indices
01274           break;
01275       }
01276     }
01277     
01278     // Unpack the other received indices without doing anything
01279     while(++n_in < remoteEntries)
01280       MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01281                type, comm_);
01282   }
01283   
01284     
01285   template<typename T, typename A>
01286   inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& send,
01287                                                   RemoteIndexList& receive,
01288                                                   int remoteEntries,
01289                                                   PairType** localSource,
01290                                                   int localSourceEntries,
01291                                                   PairType** localDest,
01292                                                   int localDestEntries,
01293                                                   char* p_in,
01294                                                   MPI_Datatype type,
01295                                                   int* position,
01296                                                   int bufferSize)
01297   {             
01298     int n_in=0, sourceIndex=0, destIndex=0;
01299         
01300     //Check if we know the global index
01301     while(n_in<remoteEntries && (sourceIndex<localSourceEntries || destIndex<localDestEntries)){
01302       // Unpack next index
01303       PairType index;
01304       MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01305                  type, comm_);
01306       n_in++;
01307       
01308       // Advance until global index in localSource and localDest are >= than the one in the unpacked index
01309       while(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()<index.global())
01310         sourceIndex++;
01311       
01312       while(destIndex<localDestEntries && localDest[destIndex]->global()<index.global())
01313         destIndex++;
01314 
01315       // Add a remote index if we found the global index.
01316       if(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()==index.global())
01317           send.push_back(RemoteIndex(index.local().attribute(), 
01318                                               localSource[sourceIndex]));
01319 
01320       if(destIndex < localDestEntries && localDest[destIndex]->global() == index.global())
01321           receive.push_back(RemoteIndex(index.local().attribute(), 
01322                                         localDest[sourceIndex]));
01323     }
01324     
01325   }
01326     
01327   template<typename T, typename A>
01328   inline void RemoteIndices<T,A>::free()
01329   {
01330     typedef typename RemoteIndexMap::iterator Iterator;
01331     Iterator lend = remoteIndices_.end();
01332     for(Iterator lists=remoteIndices_.begin(); lists != lend; ++lists){
01333       if(lists->second.first==lists->second.second){
01334         // there is only one remote index list.
01335         delete lists->second.first;
01336       }else{
01337         delete lists->second.first;
01338         delete lists->second.second;
01339       }
01340     }
01341     remoteIndices_.clear();
01342     firstBuild=true;
01343   }
01344 
01345   template<typename T, typename A>
01346   inline int RemoteIndices<T,A>::neighbours() const
01347   {
01348     return remoteIndices_.size();
01349   }
01350   
01351   template<typename T, typename A>
01352   template<bool ignorePublic>
01353   inline void RemoteIndices<T,A>::rebuild()
01354   {
01355     // Test wether a rebuild is Needed.
01356     if(firstBuild || 
01357        ignorePublic!=publicIgnored || !
01358        isSynced()){
01359       free();
01360       
01361       buildRemote<ignorePublic>();
01362 
01363       sourceSeqNo_ = source_->seqNo();
01364       destSeqNo_ = target_->seqNo();
01365       firstBuild=false;
01366       publicIgnored=ignorePublic;
01367     }
01368     
01369     
01370   }
01371   
01372   template<typename T, typename A>
01373   inline bool RemoteIndices<T,A>::isSynced() const
01374   {
01375     return sourceSeqNo_==source_->seqNo() && destSeqNo_ ==target_->seqNo();
01376   }
01377 
01378   template<typename T, typename A>
01379   template<bool mode, bool send>
01380   RemoteIndexListModifier<T,A,mode> RemoteIndices<T,A>::getModifier(int process)
01381   {
01382 
01383     // The user are on their own now!
01384     // We assume they know what they are doing and just set the
01385     // remote indices to synced status.
01386     sourceSeqNo_ = source_->seqNo();
01387     destSeqNo_ = target_->seqNo();
01388 
01389     typename RemoteIndexMap::iterator found = remoteIndices_.find(process);
01390     
01391     if(found == remoteIndices_.end())
01392     {
01393       if(source_ != target_)
01394         remoteIndices_.insert(std::make_pair(process,
01395                                              std::make_pair(new RemoteIndexList(), 
01396                                                             new RemoteIndexList())));
01397       else{
01398         RemoteIndexList* rlist = new RemoteIndexList();
01399         remoteIndices_.insert(std::make_pair(process,
01400                                              std::make_pair(rlist, rlist)));
01401         
01402         found = remoteIndices_.find(process);
01403       }
01404     }
01405     
01406     firstBuild = false;
01407     
01408     if(send)
01409       return RemoteIndexListModifier<T,A,mode>(*source_, *(found->second.first));
01410     else
01411       return RemoteIndexListModifier<T,A,mode>(*target_, *(found->second.second));
01412   }
01413   
01414   template<typename T, typename A>
01415   inline typename std::map<int, std::pair<SLList<typename RemoteIndices<T,A>::RemoteIndex,
01416                                                  typename RemoteIndices<T,A>::Allocator >*,
01417                                           SLList<typename RemoteIndices<T,A>::RemoteIndex,
01418                                                  typename RemoteIndices<T,A>::Allocator >*> >::const_iterator
01419   RemoteIndices<T,A>::begin() const
01420   {
01421     return remoteIndices_.begin();
01422   }
01423   
01424   template<typename T, typename A>
01425   inline typename std::map<int, std::pair<SLList<typename RemoteIndices<T,A>::RemoteIndex,
01426                                                  typename RemoteIndices<T,A>::Allocator >*,
01427                                           SLList<typename RemoteIndices<T,A>::RemoteIndex,
01428                                                  typename RemoteIndices<T,A>::Allocator >*> >::const_iterator
01429   RemoteIndices<T,A>::end() const
01430   {
01431     return remoteIndices_.end();
01432   }
01433 
01434 
01435   template<typename T, typename A>
01436   bool RemoteIndices<T,A>::operator==(const RemoteIndices& ri)
01437   {
01438     if(neighbours()!=ri.neighbours())
01439       return false;
01440     
01441     typedef RemoteIndexList RList;
01442     typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
01443     
01444     const const_iterator rend = remoteIndices_.end();
01445 
01446     for(const_iterator rindex = remoteIndices_.begin(), rindex1=ri.remoteIndices_.begin(); rindex!=rend; ++rindex, ++rindex1){
01447       if(rindex->first != rindex1->first)
01448         return false;
01449       if(*(rindex->second.first) != *(rindex1->second.first))
01450         return false;
01451       if(*(rindex->second.second) != *(rindex1->second.second))
01452         return false;
01453     }
01454     return true;
01455   }
01456   
01457   template<class T, class A, bool mode>
01458   RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(const ParallelIndexSet& indexSet,
01459                                                            RemoteIndexList& rList)
01460     : rList_(&rList), indexSet_(&indexSet), glist_(new GlobalList()), iter_(rList.beginModify()), end_(rList.end()), first_(true)
01461   {
01462     if(MODIFYINDEXSET){
01463       assert(indexSet_);
01464       for(ConstIterator iter=iter_; iter != end_; ++iter)
01465         glist_->push_back(iter->localIndexPair().global());
01466       giter_ = glist_->beginModify();
01467     }
01468   }
01469 
01470   template<typename T, typename A, bool mode>
01471   RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(const RemoteIndexListModifier<T,A,mode>& other)
01472     : rList_(other.rList_), indexSet_(other.indexSet_), 
01473       glist_(other.glist_), iter_(other.iter_), giter_(other.giter_), end_(other.end_), 
01474       first_(other.first_), last_(other.last_)
01475   {}
01476     
01477   template<typename T, typename A, bool mode>
01478   inline void RemoteIndexListModifier<T,A,mode>::repairLocalIndexPointers() throw(InvalidIndexSetState)
01479   { 
01480     if(MODIFYINDEXSET){
01481       // repair pointers to local index set.
01482 #ifdef DUNE_ISTL_WITH_CHECKING
01483       if(indexSet_->state()!=GROUND)
01484         DUNE_THROW(InvalidIndexSetState, "Index has to be in ground mode for repairing pointers to indices");
01485 #endif
01486       typedef typename ParallelIndexSet::const_iterator IndexIterator;
01487       typedef typename GlobalList::const_iterator GlobalIterator;
01488       typedef typename RemoteIndexList::iterator Iterator;
01489       GlobalIterator giter = glist_->begin();
01490       IndexIterator index = indexSet_->begin();
01491       
01492       for(Iterator iter=rList_->begin(); iter != end_; ++iter){
01493         while(index->global()<*giter){
01494           ++index;
01495 #ifdef DUNE_ISTL_WITH_CHECKING
01496           if(index == indexSet_.end())
01497             DUNE_THROW(ISTLError, "No such global index in set!");
01498 #endif
01499         }
01500 
01501 #ifdef DUNE_ISTL_WITH_CHECKING
01502           if(index->global() != *giter)
01503             DUNE_THROW(ISTLError, "No such global index in set!");
01504 #endif
01505           iter->localIndex_ = &(*index);
01506       }
01507     }
01508   }
01509   
01510   template<typename T, typename A, bool mode>
01511   inline void RemoteIndexListModifier<T,A,mode>::insert(const RemoteIndex& index) throw(InvalidPosition)
01512   {
01513     dune_static_assert(!mode,"Not allowed if the mode indicates that new indices"
01514                        "might be added to the underlying index set. Use "
01515                        "insert(const RemoteIndex&, const GlobalIndex&) instead");
01516     
01517 #ifdef DUNE_ISTL_WITH_CHECKING
01518     if(!first_ && index.localIndexPair().global()<last_)
01519       DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
01520 #endif
01521     // Move to the correct position
01522     while(iter_ != end_ && iter_->localIndexPair().global() < index.localIndexPair().global()){
01523       ++iter_;
01524     }
01525     
01526     // No duplicate entries allowed
01527     assert(iter_==end_ || iter_->localIndexPair().global() != index.localIndexPair().global());
01528     iter_.insert(index);
01529     last_ = index.localIndexPair().global();
01530     first_ = false;
01531   }
01532   
01533   template<typename T, typename A, bool mode>
01534   inline void RemoteIndexListModifier<T,A,mode>::insert(const RemoteIndex& index, const GlobalIndex& global) throw(InvalidPosition)
01535   {
01536      dune_static_assert(mode,"Not allowed if the mode indicates that no new indices"
01537                        "might be added to the underlying index set. Use "
01538                        "insert(const RemoteIndex&) instead");
01539 #ifdef DUNE_ISTL_WITH_CHECKING
01540     if(!first_ && global<last_)
01541       DUNE_THROW(InvalidPosition, "Modification of remote indices have to occur with ascending global index.");
01542 #endif
01543     // Move to the correct position
01544     while(iter_ != end_ && *giter_ < global){
01545       ++giter_;
01546       ++iter_;
01547     }
01548     
01549     // No duplicate entries allowed
01550     assert(iter_->localIndexPair().global() != global);
01551     iter_.insert(index);
01552     giter_.insert(global);
01553     
01554     last_ = global;
01555     first_ = false;
01556   }
01557   
01558   template<typename T, typename A, bool mode>
01559   bool RemoteIndexListModifier<T,A,mode>::remove(const GlobalIndex& global) throw(InvalidPosition)
01560   {
01561 #ifdef DUNE_ISTL_WITH_CHECKING
01562     if(!first_ && global<last_)
01563       DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
01564 #endif
01565 
01566     bool found= false;
01567     
01568     if(MODIFYINDEXSET){
01569     // Move to the correct position
01570       while(iter_!=end_ && *giter_< global){
01571         ++giter_;
01572         ++iter_;
01573       }
01574       if(*giter_ == global){
01575         giter_.remove();
01576         iter_.remove();
01577         found=true;
01578       }
01579     }else{
01580       while(iter_!=end_ && iter_->localIndexPair().global() < global)
01581         ++iter_;
01582         
01583       if(iter_->localIndexPair().global()==global){
01584         iter_.remove();
01585         found = true;
01586       }
01587     }
01588     
01589     last_ = global;
01590     first_ = false;
01591     return found;
01592   }
01593 
01594   template<typename T, typename A>
01595   template<bool send>
01596   inline typename RemoteIndices<T,A>::CollectiveIteratorT RemoteIndices<T,A>::iterator() const
01597   {
01598     return CollectiveIterator<T,A>(remoteIndices_, send);
01599   }
01600 
01601   template<typename T, typename A>
01602   inline MPI_Comm RemoteIndices<T,A>::communicator() const
01603   {
01604     return comm_;
01605     
01606   }
01607   
01608   template<typename T, typename A>
01609   CollectiveIterator<T,A>::CollectiveIterator(const RemoteIndexMap& pmap, bool send)
01610   {
01611     typedef typename RemoteIndexMap::const_iterator const_iterator;
01612     typedef typename RemoteIndexMap::iterator iterator;
01613     
01614     const const_iterator end=pmap.end();
01615     for(const_iterator process=pmap.begin(); process != end; ++process){
01616       const RemoteIndexList* list = send? process->second.first : process->second.second;
01617       typedef typename RemoteIndexList::const_iterator iterator;
01618       map_.insert(std::make_pair(process->first, 
01619                                  std::pair<iterator, const iterator>(list->begin(), list->end())));
01620     }
01621   }
01622 
01623   template<typename T, typename A>
01624   inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index)
01625   {
01626     typedef typename Map::iterator iterator;
01627     typedef typename Map::const_iterator const_iterator;
01628     const const_iterator end = map_.end();
01629     
01630     for(iterator iter = map_.begin(); iter != end;){
01631       // Step the iterator until we are >= index
01632       typename RemoteIndexList::const_iterator current = iter->second.first;
01633       typename RemoteIndexList::const_iterator rend = iter->second.second;
01634       RemoteIndex remoteIndex;
01635       if(current != rend)
01636         remoteIndex = *current;
01637       
01638       while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
01639         ++(iter->second.first);
01640       
01641       // erase from the map if there are no more entries.
01642       if(iter->second.first == iter->second.second)
01643         map_.erase(iter++);
01644       else{
01645         ++iter;
01646       }
01647     }
01648     index_=index;
01649   }
01650   
01651   template<typename T, typename A>
01652   inline bool CollectiveIterator<T,A>::empty()
01653   {
01654     return map_.empty();
01655   }
01656     
01657   template<typename T, typename A>
01658   inline typename CollectiveIterator<T,A>::iterator
01659   CollectiveIterator<T,A>::begin()
01660   {
01661     return iterator(map_.begin(), map_.end(), index_);
01662   }
01663    
01664   template<typename T, typename A>
01665   inline typename CollectiveIterator<T,A>::iterator
01666   CollectiveIterator<T,A>::end()
01667   {
01668     return iterator(map_.end(), map_.end(), index_);
01669   }
01670   
01671   template<typename TG, typename TA>
01672   inline std::ostream& operator<<(std::ostream& os, const RemoteIndex<TG,TA>& index)
01673   {
01674     os<<"[global="<<index.localIndexPair().global()<<",attribute="<<index.attribute()<<"]";
01675     return os;
01676   }
01677   
01678   template<typename T, typename A>
01679   inline std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices)
01680   {
01681     int rank;
01682     MPI_Comm_rank(indices.comm_, &rank);
01683 
01684     typedef typename RemoteIndices<T,A>::RemoteIndexList RList;
01685     typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
01686 
01687     const const_iterator rend = indices.remoteIndices_.end();
01688 
01689     for(const_iterator rindex = indices.remoteIndices_.begin(); rindex!=rend; ++rindex){
01690       os<<rank<<": Prozess "<<rindex->first<<":";
01691       
01692       if(!rindex->second.first->empty()){
01693         os<<" send:";
01694 
01695         const typename RList::const_iterator send= rindex->second.first->end();
01696       
01697         for(typename RList::const_iterator index = rindex->second.first->begin(); 
01698             index != send; ++index)
01699           os<<*index<<" ";
01700         os<<std::endl;
01701       }
01702       if(!rindex->second.second->empty()){
01703         os<<rank<<": Prozess "<<rindex->first<<": "<<"receive: ";
01704         
01705         const typename RList::const_iterator rend= rindex->second.second->end();
01706         
01707         for(typename RList::const_iterator index = rindex->second.second->begin(); 
01708             index != rend; ++index)
01709           os<<*index<<" ";
01710       }
01711       os<<std::endl<<std::flush;
01712     }
01713     return os;
01714   }
01716 }
01717 
01718 #endif
01719 #endif
Generated on Sat Apr 24 11:13:47 2010 for dune-istl by  doxygen 1.6.3