dune-common 2.1.1
remoteindices.hh
Go to the documentation of this file.
00001 // $Id$
00002 #ifndef DUNE_REMOTEINDICES_HH
00003 #define DUNE_REMOTEINDICES_HH
00004 
00005 #include"indexset.hh"
00006 #include <dune/common/exceptions.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 <dune/common/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 
00059  
00063   template<typename T1, typename T2>
00064   class RemoteIndex
00065   {
00066     template<typename T>
00067     friend class IndicesSyncer;
00068 
00069     template<typename T, typename A, typename A1>
00070     friend void repairLocalIndexPointers(std::map<int,SLList<std::pair<typename T::GlobalIndex, typename T::LocalIndex::Attribute>,A> >&, 
00071                                          RemoteIndices<T,A1>&,
00072                                          const T&);
00073 
00074     template<typename T, typename A, bool mode>
00075     friend class RemoteIndexListModifier;
00076 
00077   public:    
00082     typedef T1 GlobalIndex;
00091     typedef T2 Attribute;
00092 
00096     typedef IndexPair<GlobalIndex,ParallelLocalIndex<Attribute> >
00097     PairType;
00098     
00103     const Attribute attribute() const;
00104 
00110     const PairType& localIndexPair() const;
00111 
00115     RemoteIndex();
00116     
00117     
00123     RemoteIndex(const T2& attribute, 
00124                 const PairType* local);
00125 
00126     
00132     RemoteIndex(const T2& attribute);
00133 
00134     bool operator==(const RemoteIndex& ri) const;
00135     
00136     bool operator!=(const RemoteIndex& ri) const;
00137   private:
00139     const PairType* localIndex_;
00140 
00142     char attribute_;
00143   };
00144   
00145   template<class T, class A>
00146   std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices);
00147   
00148   class InterfaceBuilder;
00149 
00150   template<class T, class A>
00151   class CollectiveIterator;
00152 
00153   template<class T>
00154     class IndicesSyncer;
00155   
00156   // forward declaration needed for friend declaration.
00157   template<typename T1, typename T2>
00158   class OwnerOverlapCopyCommunication;
00159   
00160 
00177   template<class T, class A=std::allocator<RemoteIndex<typename T::GlobalIndex,
00178                                                        typename T::LocalIndex::Attribute> > >
00179   class RemoteIndices
00180   {
00181     friend class InterfaceBuilder;
00182     friend class IndicesSyncer<T>;
00183     template<typename T1, typename A2, typename A1>
00184     friend void repairLocalIndexPointers(std::map<int,SLList<std::pair<typename T1::GlobalIndex, typename T1::LocalIndex::Attribute>,A2> >&, 
00185                                          RemoteIndices<T1,A1>&,
00186                                          const T1&);
00187     
00188     template<class G, class T1, class T2>
00189     friend void fillIndexSetHoles(const G& graph, Dune::OwnerOverlapCopyCommunication<T1,T2>& oocomm);
00190     friend std::ostream& operator<<<>(std::ostream&, const RemoteIndices<T>&);
00191     
00192   public:
00193    
00197     typedef T ParallelIndexSet;
00198  
00201     typedef CollectiveIterator<T,A> CollectiveIteratorT;
00202     
00206     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00207     
00208     
00212     typedef typename ParallelIndexSet::LocalIndex LocalIndex;
00213 
00217     typedef typename LocalIndex::Attribute Attribute;
00218  
00222     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00223 
00224     
00228     typedef typename A::template rebind<RemoteIndex>::other Allocator;
00229 
00231     typedef Dune::SLList<RemoteIndex,Allocator>
00232     RemoteIndexList;
00233     
00235     typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
00236     RemoteIndexMap;    
00237     
00238     typedef typename RemoteIndexMap::const_iterator const_iterator;
00239     
00257     inline RemoteIndices(const ParallelIndexSet& source, const ParallelIndexSet& destination, 
00258                          const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>(), bool includeSelf=false);
00259 
00260     RemoteIndices();
00261     
00269     void setIncludeSelf(bool includeSelf);
00270     
00287     void setIndexSets(const ParallelIndexSet& source, const ParallelIndexSet& destination, 
00288                       const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>());
00289 
00290     template<typename C>
00291     void setNeighbours(const C& neighbours)
00292     {
00293       typedef typename C::const_iterator Iter;
00294       neighbourIds.clear();
00295       neighbourIds.insert(neighbours.begin(), neighbours.end());
00296         
00297     }
00298     
00302     ~RemoteIndices();
00303     
00313     template<bool ignorePublic>
00314     void rebuild();
00315 
00316     bool operator==(const RemoteIndices& ri);
00317     
00325     inline bool isSynced() const;
00326 
00330     inline MPI_Comm communicator() const;
00331 
00346     template<bool mode, bool send>
00347     inline RemoteIndexListModifier<T,A,mode> getModifier(int process);
00348 
00355     inline const_iterator find(int proc) const;
00356 
00361     inline const_iterator begin() const;
00362     
00367     inline const_iterator end() const;
00368     
00372     template<bool send>
00373     inline CollectiveIteratorT iterator() const;
00374 
00378     inline void free();
00379     
00384     inline int neighbours() const;
00385 
00387     inline const ParallelIndexSet& sourceIndexSet() const;
00388     
00390     inline const ParallelIndexSet& destinationIndexSet() const;
00391 
00392   private:
00394     RemoteIndices(const RemoteIndices&)
00395     {} 
00396 
00398     const ParallelIndexSet* source_;
00399 
00401     const ParallelIndexSet* target_;
00402 
00404     MPI_Comm comm_;
00405 
00408     std::set<int> neighbourIds;
00409     
00411     const static int commTag_=333;
00412     
00417     int sourceSeqNo_;
00418     
00423     int destSeqNo_;
00424 
00428     bool publicIgnored;
00429 
00433     bool firstBuild;
00434     
00435     /*
00436      * @brief If true, sending from indices of the processor to other 
00437      * indices on the same processor is enabled even if the same indexset is used 
00438      * on both the
00439      * sending and receiving side.
00440      */
00441     bool includeSelf;
00442     
00444     typedef IndexPair<GlobalIndex, LocalIndex> 
00445     PairType;
00446     
00453     RemoteIndexMap remoteIndices_;
00454     
00465     template<bool ignorePublic>
00466     inline void buildRemote(bool includeSelf);
00467 
00473     inline int noPublic(const ParallelIndexSet& indexSet);
00474     
00486     template<bool ignorePublic>
00487     inline void packEntries(PairType** myPairs, const ParallelIndexSet& indexSet,
00488                             char* p_out, MPI_Datatype type, int bufferSize, 
00489                             int* position, int n);
00490     
00504     inline void unpackIndices(RemoteIndexList& remote, int remoteEntries,
00505                               PairType** local, int localEntries, char* p_in, 
00506                               MPI_Datatype type, int* positon, int bufferSize,
00507                               bool fromOurself);
00508 
00509     inline void unpackIndices(RemoteIndexList& send, RemoteIndexList& receive,
00510                               int remoteEntries, PairType** localSource,
00511                               int localSourceEntries, PairType** localDest,
00512                               int localDestEntries, char* p_in,
00513                               MPI_Datatype type, int* position, int bufferSize);
00514 
00515     void unpackCreateRemote(char* p_in, PairType** sourcePairs, PairType** DestPairs, 
00516                             int remoteProc,  int sourcePublish, int destPublish, 
00517                             int bufferSize, bool sendTwo, bool fromOurSelf=false);
00518   };
00519   
00537   template<class T, class A, bool mode>
00538   class RemoteIndexListModifier
00539   {
00540 
00541     template<typename T1, typename A1>
00542     friend class RemoteIndices;
00543   
00544   public:
00545     class InvalidPosition : public RangeError
00546     {};
00547 
00548     enum {
00557       MODIFYINDEXSET=mode
00558     };
00559   
00563     typedef T ParallelIndexSet;
00564 
00568     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00569 
00573     typedef typename ParallelIndexSet::LocalIndex LocalIndex;
00574     
00578     typedef typename LocalIndex::Attribute Attribute;
00579     
00583     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00584     
00588     typedef A Allocator;
00589 
00591     typedef Dune::SLList<RemoteIndex,Allocator>
00592     RemoteIndexList;
00593 
00597     typedef SLListModifyIterator<RemoteIndex,Allocator> ModifyIterator;
00598     
00602     typedef typename RemoteIndexList::const_iterator ConstIterator;
00603     
00617     void insert(const RemoteIndex& index) throw(InvalidPosition);
00618     
00619     
00634     void insert(const RemoteIndex& index, const GlobalIndex& global) throw(InvalidPosition);
00635 
00643     bool remove(const GlobalIndex& global) throw(InvalidPosition);
00644 
00657     void repairLocalIndexPointers() throw(InvalidIndexSetState);
00658 
00659 
00660     RemoteIndexListModifier(const RemoteIndexListModifier&);
00661     
00666     RemoteIndexListModifier()
00667       : glist_()
00668     {}
00669     
00671     ~RemoteIndexListModifier()
00672     {
00673       if(glist_)
00674         delete glist_;
00675     }
00676     
00677   private:
00678     
00684     RemoteIndexListModifier(const ParallelIndexSet& indexSet,
00685                             RemoteIndexList& rList);
00686 
00687     typedef SLList<GlobalIndex,Allocator> GlobalList;
00688     typedef typename GlobalList::ModifyIterator GlobalModifyIterator;
00689     RemoteIndices<ParallelIndexSet>* remoteIndices_;
00690     RemoteIndexList* rList_;
00691     const ParallelIndexSet* indexSet_;
00692     GlobalList* glist_;
00693     ModifyIterator iter_;
00694     GlobalModifyIterator giter_;
00695     ConstIterator end_;
00696     bool first_;
00697     GlobalIndex last_;
00698   };
00699   
00704   template<class T, class A>
00705   class CollectiveIterator
00706   {
00707     
00711     typedef T ParallelIndexSet;
00712 
00716     typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
00717 
00721     typedef typename ParallelIndexSet::LocalIndex LocalIndex;
00722     
00726     typedef typename LocalIndex::Attribute Attribute;
00727 
00729     typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
00730     
00732     typedef typename A::template rebind<RemoteIndex>::other Allocator;
00733 
00735     typedef Dune::SLList<RemoteIndex,Allocator> RemoteIndexList;    
00736 
00738     typedef std::map<int,std::pair<typename RemoteIndexList::const_iterator,
00739                                    const typename RemoteIndexList::const_iterator> >
00740     Map;
00741 
00742   public:
00743     
00745     typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
00746     RemoteIndexMap;
00747     
00753     inline CollectiveIterator(const RemoteIndexMap& map_, bool send);
00754     
00763     inline void advance(const GlobalIndex& global);
00764 
00774     inline void advance(const GlobalIndex& global, const Attribute& attribute);
00775 
00776     CollectiveIterator& operator++();
00777     
00781     inline bool empty();
00782     
00789     class iterator
00790     {
00791     public:
00792       typedef typename Map::iterator RealIterator;
00793       typedef typename Map::iterator ConstRealIterator;
00794       
00795       
00797       iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex& index)
00798         : iter_(iter), end_(end), index_(index), hasAttribute(false)
00799       {
00800         // Move to the first valid entry
00801         while(iter_!=end_ && iter_->second.first->localIndexPair().global()!=index_)
00802           ++iter_;
00803       }
00804       
00805       iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex index,
00806                Attribute attribute)
00807         : iter_(iter), end_(end), index_(index), attribute_(attribute), hasAttribute(true)
00808       {
00809         // Move to the first valid entry or the end
00810         while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_
00811                               || iter_->second.first->localIndexPair().local().attribute()!=attribute))
00812           ++iter_;
00813       }
00815       iterator(const iterator& other)
00816         : iter_(other.iter_), end_(other.end_), index_(other.index_)
00817       { }
00818          
00820       iterator& operator++()
00821       { 
00822         ++iter_;
00823         // If entry is not valid move on
00824         while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_ ||
00825               (hasAttribute &&
00826                iter_->second.first->localIndexPair().local().attribute()!=attribute_)))
00827           ++iter_;
00828         assert(iter_==end_ || 
00829                (iter_->second.first->localIndexPair().global()==index_));
00830         assert(iter_==end_ || !hasAttribute ||
00831                (iter_->second.first->localIndexPair().local().attribute()==attribute_));
00832         return *this;
00833       }
00834       
00836       const RemoteIndex& operator*()const
00837       {
00838         return *(iter_->second.first);
00839       }
00840       
00842       int process() const
00843       {
00844         return iter_->first;
00845       }
00846       
00848       const RemoteIndex* operator->()const
00849       {
00850         return iter_->second.first.operator->();
00851       }
00852       
00854       bool operator==(const iterator& other)
00855       {
00856         return other.iter_==iter_;
00857       }
00858 
00860       bool operator!=(const iterator& other)
00861       {
00862         return other.iter_!=iter_;
00863       }
00864 
00865     private:
00866       iterator();
00867       
00868       RealIterator iter_;
00869       RealIterator end_;
00870       GlobalIndex index_;
00871       Attribute attribute_;
00872       bool hasAttribute;
00873     };
00874     
00875     iterator begin();
00876     
00877     iterator end();
00878     
00879   private:
00880     
00881     Map map_;
00882     GlobalIndex index_;
00883     Attribute attribute_;
00884     bool noattribute;
00885   };
00886     
00887   template<typename TG, typename TA>
00888   MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::getType()
00889   {
00890     if(type==MPI_DATATYPE_NULL){
00891       int length[4];
00892       MPI_Aint disp[4];
00893       MPI_Datatype types[4] = {MPI_LB, MPITraits<TG>::getType(), 
00894                                MPITraits<ParallelLocalIndex<TA> >::getType(), MPI_UB};
00895       IndexPair<TG,ParallelLocalIndex<TA> > rep[2];
00896       length[0]=length[1]=length[2]=length[3]=1;
00897       MPI_Address(rep, disp); // lower bound of the datatype
00898       MPI_Address(&(rep[0].global_), disp+1);
00899       MPI_Address(&(rep[0].local_), disp+2);
00900       MPI_Address(rep+1, disp+3); // upper bound of the datatype
00901       for(int i=3; i >= 0; --i)
00902         disp[i] -= disp[0];
00903       MPI_Type_struct(4, length, disp, types, &type);
00904       MPI_Type_commit(&type);
00905     }
00906     return type;
00907   }
00908 
00909   template<typename TG, typename TA>
00910   MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::type=MPI_DATATYPE_NULL;
00911 
00912   template<typename T1, typename T2>
00913   RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute, const PairType* local)
00914     : localIndex_(local), attribute_(attribute)
00915   {}
00916 
00917   template<typename T1, typename T2>
00918   RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute)
00919     : localIndex_(0), attribute_(attribute)
00920   {}
00921 
00922   template<typename T1, typename T2>
00923   RemoteIndex<T1,T2>::RemoteIndex()
00924     : localIndex_(0), attribute_()
00925   {}
00926   template<typename T1, typename T2>
00927   inline bool RemoteIndex<T1,T2>::operator==(const RemoteIndex& ri) const
00928   {
00929     return localIndex_==ri.localIndex_ && attribute_==ri.attribute;
00930   }
00931   
00932   template<typename T1, typename T2>
00933   inline bool RemoteIndex<T1,T2>::operator!=(const RemoteIndex& ri) const
00934   {
00935     return localIndex_!=ri.localIndex_ || attribute_!=ri.attribute_;
00936   }
00937   
00938   template<typename T1, typename T2>
00939   inline const T2 RemoteIndex<T1,T2>::attribute() const
00940   {
00941     return T2(attribute_);
00942   }
00943 
00944   template<typename T1, typename T2>
00945   inline const IndexPair<T1,ParallelLocalIndex<T2> >& RemoteIndex<T1,T2>::localIndexPair() const
00946   {
00947     return *localIndex_;
00948   }
00949   
00950   template<typename T, typename A>
00951   inline RemoteIndices<T,A>::RemoteIndices(const ParallelIndexSet& source,
00952                                          const ParallelIndexSet& destination,
00953                                          const MPI_Comm& comm,
00954                                            const std::vector<int>& neighbours,
00955                                            bool includeSelf_)
00956     : source_(&source), target_(&destination), comm_(comm),
00957       sourceSeqNo_(-1), destSeqNo_(-1), publicIgnored(false), firstBuild(true),
00958       includeSelf(includeSelf_)
00959   {
00960     setNeighbours(neighbours);
00961   }
00962 
00963   template<typename T, typename A>
00964   void RemoteIndices<T,A>::setIncludeSelf(bool b)
00965   {
00966     includeSelf=b;
00967   }
00968   
00969   template<typename T, typename A>
00970   RemoteIndices<T,A>::RemoteIndices()
00971     :source_(0), target_(0), sourceSeqNo_(-1), 
00972      destSeqNo_(-1), publicIgnored(false), firstBuild(true)
00973   {}
00974   
00975   template<class T, typename A>
00976   void RemoteIndices<T,A>::setIndexSets(const ParallelIndexSet& source,
00977                                       const ParallelIndexSet& destination,
00978                                       const MPI_Comm& comm,
00979                                       const std::vector<int>& neighbours)
00980   {
00981     free();
00982     source_ = &source;
00983     target_ = &destination;
00984     comm_ = comm;
00985     firstBuild = true;
00986     setNeighbours(neighbours);
00987   }
00988 
00989   template<typename T, typename A>
00990   const typename RemoteIndices<T,A>::ParallelIndexSet& 
00991   RemoteIndices<T,A>::sourceIndexSet() const
00992   {
00993     return *source_;
00994   }
00995   
00996     
00997   template<typename T, typename A>
00998   const typename RemoteIndices<T,A>::ParallelIndexSet& 
00999   RemoteIndices<T,A>::destinationIndexSet() const
01000   {
01001     return *target_;
01002   }
01003   
01004   
01005   template<typename T, typename A>
01006   RemoteIndices<T,A>::~RemoteIndices()
01007   {
01008     free();
01009   }
01010 
01011   template<typename T, typename A>
01012   template<bool ignorePublic>
01013   inline void RemoteIndices<T,A>::packEntries(IndexPair<GlobalIndex,LocalIndex>** pairs,
01014                                                 const ParallelIndexSet& indexSet,
01015                                                 char* p_out, MPI_Datatype type, 
01016                                                 int bufferSize,
01017                                                 int *position, int n)
01018   {
01019     // fill with own indices
01020     typedef typename ParallelIndexSet::const_iterator const_iterator;
01021     typedef IndexPair<GlobalIndex,LocalIndex> PairType;
01022     const const_iterator end = indexSet.end();
01023 
01024     //Now pack the source indices
01025     int i=0;
01026     for(const_iterator index = indexSet.begin(); index != end; ++index)
01027       if(ignorePublic || index->local().isPublic()){
01028         
01029         MPI_Pack(const_cast<PairType*>(&(*index)), 1, 
01030                  type,
01031                  p_out, bufferSize, position, comm_);
01032         pairs[i++] = const_cast<PairType*>(&(*index));
01033       
01034       }
01035     assert(i==n);
01036   }
01037   
01038   template<typename T, typename A>
01039   inline int RemoteIndices<T,A>::noPublic(const ParallelIndexSet& indexSet)
01040   {
01041     typedef typename ParallelIndexSet::const_iterator const_iterator;
01042     
01043     int noPublic=0;
01044     
01045     const const_iterator end=indexSet.end();
01046     for(const_iterator index=indexSet.begin(); index!=end; ++index)
01047       if(index->local().isPublic())
01048         noPublic++;
01049     
01050     return noPublic;
01051     
01052   }
01053     
01054 
01055   template<typename T, typename A>
01056   inline void RemoteIndices<T,A>::unpackCreateRemote(char* p_in, PairType** sourcePairs,
01057                                                      PairType** destPairs, int remoteProc, 
01058                                                      int sourcePublish, int destPublish,
01059                                                      int bufferSize, bool sendTwo,
01060                                                      bool fromOurSelf)
01061   {
01062     
01063       // unpack the number of indices we received
01064       int noRemoteSource=-1, noRemoteDest=-1;
01065       char twoIndexSets=0;
01066       int position=0;
01067       // Did we receive two index sets?
01068       MPI_Unpack(p_in, bufferSize, &position, &twoIndexSets, 1, MPI_CHAR, comm_);
01069       // The number of source indices received
01070       MPI_Unpack(p_in, bufferSize, &position, &noRemoteSource, 1, MPI_INT, comm_);
01071       // The number of destination indices received
01072       MPI_Unpack(p_in, bufferSize, &position, &noRemoteDest, 1, MPI_INT, comm_);      
01073 
01074 
01075       // Indices for which we receive
01076       RemoteIndexList* receive= new RemoteIndexList();
01077       // Indices for which we send
01078       RemoteIndexList* send=0;
01079 
01080       MPI_Datatype type= MPITraits<PairType>::getType();
01081       
01082       if(!twoIndexSets){
01083         if(sendTwo){
01084           send = new RemoteIndexList();
01085           // Create both remote index sets simultaneously
01086           unpackIndices(*send, *receive, noRemoteSource, sourcePairs, sourcePublish,
01087                         destPairs, destPublish, p_in, type, &position, bufferSize);
01088         }else{
01089           // we only need one list
01090           unpackIndices(*receive, noRemoteSource, sourcePairs, sourcePublish,
01091                         p_in, type, &position, bufferSize, fromOurSelf);
01092           send=receive;
01093         }
01094       }else{
01095         
01096         int oldPos=position;
01097         // Two index sets received
01098         unpackIndices(*receive, noRemoteSource, destPairs, destPublish,
01099                       p_in, type, &position, bufferSize, fromOurSelf);
01100         if(!sendTwo)
01101           //unpack source entries again as destination entries
01102           position=oldPos;
01103         
01104         send = new RemoteIndexList();
01105         unpackIndices(*send, noRemoteDest, sourcePairs, sourcePublish, 
01106                       p_in, type, &position, bufferSize, fromOurSelf);
01107       }
01108       
01109       if(receive->empty() && send->empty()){
01110         if(send==receive){
01111           delete send;
01112         }else{
01113           delete send;
01114           delete receive;
01115         }
01116       }else{
01117         remoteIndices_.insert(std::make_pair(remoteProc, 
01118                                              std::make_pair(send,receive)));
01119       }
01120   }
01121   
01122 
01123   template<typename T, typename A>
01124   template<bool ignorePublic>
01125   inline void RemoteIndices<T,A>::buildRemote(bool includeSelf)
01126   {
01127     // Processor configuration
01128     int rank, procs;
01129     MPI_Comm_rank(comm_, &rank);
01130     MPI_Comm_size(comm_, &procs);
01131     
01132     // number of local indices to publish
01133     // The indices of the destination will be send.
01134     int sourcePublish, destPublish;
01135     
01136     // Do we need to send two index sets?
01137     char sendTwo = (source_ != target_);
01138     
01139     if(procs==1 && !(sendTwo || includeSelf))
01140       // Nothing to communicate
01141       return;   
01142 
01143     sourcePublish = (ignorePublic)? source_->size() : noPublic(*source_);
01144     
01145     if(sendTwo)
01146       destPublish = (ignorePublic)? target_->size() : noPublic(*target_);
01147     else  
01148       // we only need to send one set of indices
01149       destPublish = 0;
01150         
01151     int maxPublish, publish=sourcePublish+destPublish;
01152 
01153     // Calucate maximum number of indices send
01154     MPI_Allreduce(&publish, &maxPublish, 1, MPI_INT, MPI_MAX, comm_);
01155 
01156     // allocate buffers
01157     typedef IndexPair<GlobalIndex,LocalIndex> PairType;
01158     
01159     PairType** destPairs;
01160     PairType** sourcePairs = new PairType*[sourcePublish>0?sourcePublish:1];
01161     
01162     if(sendTwo)
01163       destPairs = new PairType*[destPublish>0?destPublish:1];
01164     else
01165       destPairs=sourcePairs;
01166     
01167     char** buffer = new char*[2];
01168     int bufferSize;    
01169     int position=0;
01170     int intSize;
01171     int charSize;
01172 
01173     // calculate buffer size
01174     MPI_Datatype type = MPITraits<PairType>::getType();
01175     
01176     MPI_Pack_size(maxPublish, type, comm_,
01177                   &bufferSize);
01178     MPI_Pack_size(1, MPI_INT, comm_,
01179                   &intSize);
01180     MPI_Pack_size(1, MPI_CHAR, comm_,
01181                   &charSize);
01182     // Our message will contain the following:
01183     // a bool wether two index sets where sent
01184     // the size of the source and the dest indexset, 
01185     // then the source and destination indices
01186     bufferSize += 2 * intSize + charSize;
01187     
01188     if(bufferSize<=0) bufferSize=1;
01189     
01190     buffer[0] = new char[bufferSize];
01191     buffer[1] = new char[bufferSize];
01192     
01193     
01194     // pack entries into buffer[0], p_out below!
01195     MPI_Pack(&sendTwo, 1, MPI_CHAR, buffer[0], bufferSize, &position,
01196                  comm_);
01197         
01198     // The number of indices we send for each index set
01199     MPI_Pack(&sourcePublish, 1, MPI_INT, buffer[0], bufferSize, &position, 
01200              comm_);
01201     MPI_Pack(&destPublish, 1, MPI_INT, buffer[0], bufferSize, &position,
01202              comm_);
01203         
01204     // Now pack the source indices and setup the destination pairs
01205     packEntries<ignorePublic>(sourcePairs, *source_, buffer[0], type, 
01206                               bufferSize, &position, sourcePublish);
01207     // If necessary send the dest indices and setup the source pairs
01208     if(sendTwo)
01209       packEntries<ignorePublic>(destPairs, *target_, buffer[0], type,
01210                                 bufferSize, &position, destPublish);
01211 
01212 
01213     // Update remote indices for ourself
01214     if(sendTwo|| includeSelf)
01215       unpackCreateRemote(buffer[0], sourcePairs, destPairs, rank, sourcePublish, 
01216                          destPublish, bufferSize, sendTwo, includeSelf);
01217 
01218     neighbourIds.erase(rank);
01219 
01220     if(neighbourIds.size()==0)
01221       {
01222         std::cout<<rank<<": Sending messages in a ring"<<std::endl;
01223         // send mesages in ring
01224         for(int proc=1; proc<procs; proc++){
01225           // pointers to the current input and output buffers
01226           char* p_out = buffer[1-(proc%2)];
01227           char* p_in = buffer[proc%2];
01228           
01229           MPI_Status status;
01230           if(rank%2==0){
01231             MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
01232                       commTag_, comm_);
01233             MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
01234                      commTag_, comm_, &status);
01235           }else{
01236             MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
01237                      commTag_, comm_, &status);
01238             MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
01239                       commTag_, comm_);
01240           }
01241 
01242 
01243           // The process these indices are from
01244           int remoteProc = (rank+procs-proc)%procs;
01245           
01246           unpackCreateRemote(p_in, sourcePairs, destPairs, remoteProc, sourcePublish, 
01247                              destPublish, bufferSize, sendTwo);
01248           
01249         }
01250         
01251       }
01252     else
01253       {
01254         MPI_Request* requests=new MPI_Request[neighbourIds.size()];
01255         MPI_Request* req=requests;
01256         
01257         typedef typename std::set<int>::size_type size_type;
01258         size_type noNeighbours=neighbourIds.size();
01259 
01260         // setup sends  
01261         for(std::set<int>::iterator neighbour=neighbourIds.begin();
01262             neighbour!= neighbourIds.end(); ++neighbour){
01263             // Only send the information to the neighbouring processors
01264             MPI_Issend(buffer[0], position , MPI_PACKED, *neighbour, commTag_, comm_, req++);
01265         }
01266         
01267         //Test for received messages
01268         
01269         for(size_type received=0; received <noNeighbours; ++received)
01270           {
01271             MPI_Status status;
01272             // probe for next message
01273             MPI_Probe(MPI_ANY_SOURCE, commTag_, comm_, &status);
01274             int remoteProc=status.MPI_SOURCE;
01275             int size;
01276             MPI_Get_count(&status, MPI_PACKED, &size);
01277             // receive message
01278             MPI_Recv(buffer[1], size, MPI_PACKED, remoteProc,
01279                      commTag_, comm_, &status);
01280               
01281             unpackCreateRemote(buffer[1], sourcePairs, destPairs, remoteProc, sourcePublish, 
01282                                destPublish, bufferSize, sendTwo);
01283           }
01284         // wait for completion of pending requests
01285         MPI_Status* statuses = new MPI_Status[neighbourIds.size()];
01286         
01287         if(MPI_ERR_IN_STATUS==MPI_Waitall(neighbourIds.size(), requests, statuses)){
01288           for(size_type i=0; i < neighbourIds.size(); ++i)
01289             if(statuses[i].MPI_ERROR!=MPI_SUCCESS){
01290               std::cerr<<rank<<": MPI_Error occurred while receiving message."<<std::endl;
01291               MPI_Abort(comm_, 999);
01292           }
01293         }
01294         delete[] requests;
01295         delete[] statuses;
01296       }
01297     
01298     
01299     // delete allocated memory
01300     if(destPairs!=sourcePairs)
01301       delete[] destPairs;
01302     
01303     delete[] sourcePairs;
01304     delete[] buffer[0];
01305     delete[] buffer[1];
01306     delete[] buffer;
01307   }
01308   
01309   template<typename T, typename A>
01310   inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& remote,
01311                                                 int remoteEntries,
01312                                                 PairType** local,
01313                                                 int localEntries,
01314                                                 char* p_in,
01315                                                 MPI_Datatype type,
01316                                                 int* position,
01317                                                 int bufferSize,
01318                                                 bool fromOurSelf)
01319   {
01320     if(remoteEntries==0)
01321       return;
01322   
01323     PairType index(-1);
01324     MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01325                type, comm_);
01326     GlobalIndex oldGlobal=index.global();
01327     int n_in=0, localIndex=0;
01328         
01329     //Check if we know the global index
01330     while(localIndex<localEntries){
01331       if(local[localIndex]->global()==index.global()){
01332           int oldLocalIndex=localIndex;
01333           
01334           while(localIndex<localEntries && 
01335                 local[localIndex]->global()==index.global()){
01336             if(!fromOurSelf || index.local().attribute() != 
01337                local[localIndex]->local().attribute())
01338               // if index is from us it has to have a different attribute
01339               remote.push_back(RemoteIndex(index.local().attribute(), 
01340                                            local[localIndex]));
01341             localIndex++;
01342           }
01343           
01344         // unpack next remote index
01345         if((++n_in) < remoteEntries){
01346           MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01347                      type, comm_);
01348           if(index.global()==oldGlobal)
01349             // Restart comparison for the same global indices
01350             localIndex=oldLocalIndex;
01351           else
01352             oldGlobal=index.global();
01353         }else{
01354           // No more received indices
01355           break;
01356         }
01357         continue;
01358       }
01359       
01360       if (local[localIndex]->global()<index.global()){
01361         // compare with next entry in our list
01362         ++localIndex;
01363       }else{
01364         // We do not know the index, unpack next
01365         if((++n_in) < remoteEntries){
01366           MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01367                      type, comm_);
01368           oldGlobal=index.global();
01369         }else
01370           // No more received indices
01371           break;
01372       }
01373     }
01374     
01375     // Unpack the other received indices without doing anything
01376     while(++n_in < remoteEntries)
01377       MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01378                type, comm_);
01379   }
01380   
01381     
01382   template<typename T, typename A>
01383   inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& send,
01384                                                   RemoteIndexList& receive,
01385                                                   int remoteEntries,
01386                                                   PairType** localSource,
01387                                                   int localSourceEntries,
01388                                                   PairType** localDest,
01389                                                   int localDestEntries,
01390                                                   char* p_in,
01391                                                   MPI_Datatype type,
01392                                                   int* position,
01393                                                   int bufferSize)
01394   {             
01395     int n_in=0, sourceIndex=0, destIndex=0;
01396         
01397     //Check if we know the global index
01398     while(n_in<remoteEntries && (sourceIndex<localSourceEntries || destIndex<localDestEntries)){
01399       // Unpack next index
01400       PairType index;
01401       MPI_Unpack(p_in, bufferSize, position, &index, 1, 
01402                  type, comm_);
01403       n_in++;
01404       
01405       // Advance until global index in localSource and localDest are >= than the one in the unpacked index
01406       while(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()<index.global())
01407         sourceIndex++;
01408       
01409       while(destIndex<localDestEntries && localDest[destIndex]->global()<index.global())
01410         destIndex++;
01411 
01412       // Add a remote index if we found the global index.
01413       if(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()==index.global())
01414           send.push_back(RemoteIndex(index.local().attribute(), 
01415                                               localSource[sourceIndex]));
01416 
01417       if(destIndex < localDestEntries && localDest[destIndex]->global() == index.global())
01418           receive.push_back(RemoteIndex(index.local().attribute(), 
01419                                         localDest[sourceIndex]));
01420     }
01421     
01422   }
01423     
01424   template<typename T, typename A>
01425   inline void RemoteIndices<T,A>::free()
01426   {
01427     typedef typename RemoteIndexMap::iterator Iterator;
01428     Iterator lend = remoteIndices_.end();
01429     for(Iterator lists=remoteIndices_.begin(); lists != lend; ++lists){
01430       if(lists->second.first==lists->second.second){
01431         // there is only one remote index list.
01432         delete lists->second.first;
01433       }else{
01434         delete lists->second.first;
01435         delete lists->second.second;
01436       }
01437     }
01438     remoteIndices_.clear();
01439     firstBuild=true;
01440   }
01441 
01442   template<typename T, typename A>
01443   inline int RemoteIndices<T,A>::neighbours() const
01444   {
01445     return remoteIndices_.size();
01446   }
01447   
01448   template<typename T, typename A>
01449   template<bool ignorePublic>
01450   inline void RemoteIndices<T,A>::rebuild()
01451   {
01452     // Test wether a rebuild is Needed.
01453     if(firstBuild || 
01454        ignorePublic!=publicIgnored || !
01455        isSynced()){
01456       free();
01457       
01458       buildRemote<ignorePublic>(includeSelf);
01459 
01460       sourceSeqNo_ = source_->seqNo();
01461       destSeqNo_ = target_->seqNo();
01462       firstBuild=false;
01463       publicIgnored=ignorePublic;
01464     }
01465     
01466     
01467   }
01468   
01469   template<typename T, typename A>
01470   inline bool RemoteIndices<T,A>::isSynced() const
01471   {
01472     return sourceSeqNo_==source_->seqNo() && destSeqNo_ ==target_->seqNo();
01473   }
01474 
01475   template<typename T, typename A>
01476   template<bool mode, bool send>
01477   RemoteIndexListModifier<T,A,mode> RemoteIndices<T,A>::getModifier(int process)
01478   {
01479 
01480     // The user are on their own now!
01481     // We assume they know what they are doing and just set the
01482     // remote indices to synced status.
01483     sourceSeqNo_ = source_->seqNo();
01484     destSeqNo_ = target_->seqNo();
01485 
01486     typename RemoteIndexMap::iterator found = remoteIndices_.find(process);
01487     
01488     if(found == remoteIndices_.end())
01489     {
01490       if(source_ != target_)
01491         remoteIndices_.insert(std::make_pair(process,
01492                                              std::make_pair(new RemoteIndexList(), 
01493                                                             new RemoteIndexList())));
01494       else{
01495         RemoteIndexList* rlist = new RemoteIndexList();
01496         remoteIndices_.insert(std::make_pair(process,
01497                                              std::make_pair(rlist, rlist)));
01498         
01499         found = remoteIndices_.find(process);
01500       }
01501     }
01502     
01503     firstBuild = false;
01504     
01505     if(send)
01506       return RemoteIndexListModifier<T,A,mode>(*source_, *(found->second.first));
01507     else
01508       return RemoteIndexListModifier<T,A,mode>(*target_, *(found->second.second));
01509   }
01510 
01511   template<typename T, typename A>
01512   inline typename RemoteIndices<T,A>::const_iterator
01513   RemoteIndices<T,A>::find(int proc) const
01514   {
01515     return remoteIndices_.find(proc);
01516   }
01517 
01518   template<typename T, typename A>
01519   inline typename RemoteIndices<T,A>::const_iterator
01520   RemoteIndices<T,A>::begin() const
01521   {
01522     return remoteIndices_.begin();
01523   }
01524   
01525   template<typename T, typename A>
01526   inline typename RemoteIndices<T,A>::const_iterator
01527   RemoteIndices<T,A>::end() const
01528   {
01529     return remoteIndices_.end();
01530   }
01531 
01532 
01533   template<typename T, typename A>
01534   bool RemoteIndices<T,A>::operator==(const RemoteIndices& ri)
01535   {
01536     if(neighbours()!=ri.neighbours())
01537       return false;
01538     
01539     typedef RemoteIndexList RList;
01540     typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
01541     
01542     const const_iterator rend = remoteIndices_.end();
01543 
01544     for(const_iterator rindex = remoteIndices_.begin(), rindex1=ri.remoteIndices_.begin(); rindex!=rend; ++rindex, ++rindex1){
01545       if(rindex->first != rindex1->first)
01546         return false;
01547       if(*(rindex->second.first) != *(rindex1->second.first))
01548         return false;
01549       if(*(rindex->second.second) != *(rindex1->second.second))
01550         return false;
01551     }
01552     return true;
01553   }
01554   
01555   template<class T, class A, bool mode>
01556   RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(const ParallelIndexSet& indexSet,
01557                                                            RemoteIndexList& rList)
01558     : rList_(&rList), indexSet_(&indexSet), glist_(new GlobalList()), iter_(rList.beginModify()), end_(rList.end()), first_(true)
01559   {
01560     if(MODIFYINDEXSET){
01561       assert(indexSet_);
01562       for(ConstIterator iter=iter_; iter != end_; ++iter)
01563         glist_->push_back(iter->localIndexPair().global());
01564       giter_ = glist_->beginModify();
01565     }
01566   }
01567 
01568   template<typename T, typename A, bool mode>
01569   RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(const RemoteIndexListModifier<T,A,mode>& other)
01570     : rList_(other.rList_), indexSet_(other.indexSet_), 
01571       glist_(other.glist_), iter_(other.iter_), giter_(other.giter_), end_(other.end_), 
01572       first_(other.first_), last_(other.last_)
01573   {}
01574     
01575   template<typename T, typename A, bool mode>
01576   inline void RemoteIndexListModifier<T,A,mode>::repairLocalIndexPointers() throw(InvalidIndexSetState)
01577   { 
01578     if(MODIFYINDEXSET){
01579       // repair pointers to local index set.
01580 #ifdef DUNE_ISTL_WITH_CHECKING
01581       if(indexSet_->state()!=GROUND)
01582         DUNE_THROW(InvalidIndexSetState, "Index has to be in ground mode for repairing pointers to indices");
01583 #endif
01584       typedef typename ParallelIndexSet::const_iterator IndexIterator;
01585       typedef typename GlobalList::const_iterator GlobalIterator;
01586       typedef typename RemoteIndexList::iterator Iterator;
01587       GlobalIterator giter = glist_->begin();
01588       IndexIterator index = indexSet_->begin();
01589       
01590       for(Iterator iter=rList_->begin(); iter != end_; ++iter){
01591         while(index->global()<*giter){
01592           ++index;
01593 #ifdef DUNE_ISTL_WITH_CHECKING
01594           if(index == indexSet_.end())
01595             DUNE_THROW(InvalidPosition, "No such global index in set!");
01596 #endif
01597         }
01598 
01599 #ifdef DUNE_ISTL_WITH_CHECKING
01600           if(index->global() != *giter)
01601             DUNE_THROW(InvalidPosition, "No such global index in set!");
01602 #endif
01603           iter->localIndex_ = &(*index);
01604       }
01605     }
01606   }
01607   
01608   template<typename T, typename A, bool mode>
01609   inline void RemoteIndexListModifier<T,A,mode>::insert(const RemoteIndex& index) throw(InvalidPosition)
01610   {
01611     dune_static_assert(!mode,"Not allowed if the mode indicates that new indices"
01612                        "might be added to the underlying index set. Use "
01613                        "insert(const RemoteIndex&, const GlobalIndex&) instead");
01614     
01615 #ifdef DUNE_ISTL_WITH_CHECKING
01616     if(!first_ && index.localIndexPair().global()<last_)
01617       DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
01618 #endif
01619     // Move to the correct position
01620     while(iter_ != end_ && iter_->localIndexPair().global() < index.localIndexPair().global()){
01621       ++iter_;
01622     }
01623     
01624     // No duplicate entries allowed
01625     assert(iter_==end_ || iter_->localIndexPair().global() != index.localIndexPair().global());
01626     iter_.insert(index);
01627     last_ = index.localIndexPair().global();
01628     first_ = false;
01629   }
01630   
01631   template<typename T, typename A, bool mode>
01632   inline void RemoteIndexListModifier<T,A,mode>::insert(const RemoteIndex& index, const GlobalIndex& global) throw(InvalidPosition)
01633   {
01634      dune_static_assert(mode,"Not allowed if the mode indicates that no new indices"
01635                        "might be added to the underlying index set. Use "
01636                        "insert(const RemoteIndex&) instead");
01637 #ifdef DUNE_ISTL_WITH_CHECKING
01638     if(!first_ && global<last_)
01639       DUNE_THROW(InvalidPosition, "Modification of remote indices have to occur with ascending global index.");
01640 #endif
01641     // Move to the correct position
01642     while(iter_ != end_ && *giter_ < global){
01643       ++giter_;
01644       ++iter_;
01645     }
01646     
01647     // No duplicate entries allowed
01648     assert(iter_->localIndexPair().global() != global);
01649     iter_.insert(index);
01650     giter_.insert(global);
01651     
01652     last_ = global;
01653     first_ = false;
01654   }
01655   
01656   template<typename T, typename A, bool mode>
01657   bool RemoteIndexListModifier<T,A,mode>::remove(const GlobalIndex& global) throw(InvalidPosition)
01658   {
01659 #ifdef DUNE_ISTL_WITH_CHECKING
01660     if(!first_ && global<last_)
01661       DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
01662 #endif
01663 
01664     bool found= false;
01665     
01666     if(MODIFYINDEXSET){
01667     // Move to the correct position
01668       while(iter_!=end_ && *giter_< global){
01669         ++giter_;
01670         ++iter_;
01671       }
01672       if(*giter_ == global){
01673         giter_.remove();
01674         iter_.remove();
01675         found=true;
01676       }
01677     }else{
01678       while(iter_!=end_ && iter_->localIndexPair().global() < global)
01679         ++iter_;
01680         
01681       if(iter_->localIndexPair().global()==global){
01682         iter_.remove();
01683         found = true;
01684       }
01685     }
01686     
01687     last_ = global;
01688     first_ = false;
01689     return found;
01690   }
01691 
01692   template<typename T, typename A>
01693   template<bool send>
01694   inline typename RemoteIndices<T,A>::CollectiveIteratorT RemoteIndices<T,A>::iterator() const
01695   {
01696     return CollectiveIterator<T,A>(remoteIndices_, send);
01697   }
01698 
01699   template<typename T, typename A>
01700   inline MPI_Comm RemoteIndices<T,A>::communicator() const
01701   {
01702     return comm_;
01703     
01704   }
01705   
01706   template<typename T, typename A>
01707   CollectiveIterator<T,A>::CollectiveIterator(const RemoteIndexMap& pmap, bool send)
01708   {
01709     typedef typename RemoteIndexMap::const_iterator const_iterator;
01710     typedef typename RemoteIndexMap::iterator iterator;
01711     
01712     const const_iterator end=pmap.end();
01713     for(const_iterator process=pmap.begin(); process != end; ++process){
01714       const RemoteIndexList* list = send? process->second.first : process->second.second;
01715       typedef typename RemoteIndexList::const_iterator iterator;
01716       map_.insert(std::make_pair(process->first, 
01717                                  std::pair<iterator, const iterator>(list->begin(), list->end())));
01718     }
01719   }
01720 
01721   template<typename T, typename A>
01722   inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index)
01723   {
01724     typedef typename Map::iterator iterator;
01725     typedef typename Map::const_iterator const_iterator;
01726     const const_iterator end = map_.end();
01727     
01728     for(iterator iter = map_.begin(); iter != end;){
01729       // Step the iterator until we are >= index
01730       typename RemoteIndexList::const_iterator current = iter->second.first;
01731       typename RemoteIndexList::const_iterator rend = iter->second.second;
01732       RemoteIndex remoteIndex;
01733       if(current != rend)
01734         remoteIndex = *current;
01735       
01736       while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
01737         ++(iter->second.first);
01738       
01739       // erase from the map if there are no more entries.
01740       if(iter->second.first == iter->second.second)
01741         map_.erase(iter++);
01742       else{
01743         ++iter;
01744       }
01745     }
01746     index_=index;
01747     noattribute=true;
01748   }
01749 
01750   template<typename T, typename A>
01751   inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index,
01752                                                const Attribute& attribute)
01753   {
01754     typedef typename Map::iterator iterator;
01755     typedef typename Map::const_iterator const_iterator;
01756     const const_iterator end = map_.end();
01757     
01758     for(iterator iter = map_.begin(); iter != end;){
01759       // Step the iterator until we are >= index
01760       typename RemoteIndexList::const_iterator current = iter->second.first;
01761       typename RemoteIndexList::const_iterator rend = iter->second.second;
01762       RemoteIndex remoteIndex;
01763       if(current != rend)
01764         remoteIndex = *current;
01765 
01766       // Move to global index or bigger
01767       while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
01768         ++(iter->second.first);
01769       
01770       // move to attribute or bigger
01771       while(iter->second.first!=iter->second.second 
01772             && iter->second.first->localIndexPair().global()==index
01773             && iter->second.first->localIndexPair().local().attribute()<attribute)
01774         ++(iter->second.first);
01775 
01776       // erase from the map if there are no more entries.
01777       if(iter->second.first == iter->second.second)
01778         map_.erase(iter++);
01779       else{
01780         ++iter;
01781       }
01782     }
01783     index_=index;
01784     attribute_=attribute;
01785     noattribute=false;
01786   }
01787 
01788   template<typename T, typename A>
01789   inline  CollectiveIterator<T,A>& CollectiveIterator<T,A>::operator++()
01790   {
01791     typedef typename Map::iterator iterator;
01792     typedef typename Map::const_iterator const_iterator;
01793     const const_iterator end = map_.end();
01794     
01795     for(iterator iter = map_.begin(); iter != end;){
01796       // Step the iterator until we are >= index
01797       typename RemoteIndexList::const_iterator current = iter->second.first;
01798       typename RemoteIndexList::const_iterator rend = iter->second.second;
01799 
01800       // move all iterators pointing to the current global index to next value
01801       if(iter->second.first->localIndexPair().global()==index_ &&
01802          (noattribute || iter->second.first->localIndexPair().local().attribute() == attribute_))
01803         ++(iter->second.first);
01804       
01805       // erase from the map if there are no more entries.
01806       if(iter->second.first == iter->second.second)
01807         map_.erase(iter++);
01808       else{
01809         ++iter;
01810       }
01811     }
01812     return *this;
01813   }
01814 
01815   template<typename T, typename A>
01816   inline bool CollectiveIterator<T,A>::empty()
01817   {
01818     return map_.empty();
01819   }
01820     
01821   template<typename T, typename A>
01822   inline typename CollectiveIterator<T,A>::iterator
01823   CollectiveIterator<T,A>::begin()
01824   {
01825     if(noattribute)
01826       return iterator(map_.begin(), map_.end(), index_);
01827     else
01828       return iterator(map_.begin(), map_.end(), index_,
01829                       attribute_);
01830   }
01831    
01832   template<typename T, typename A>
01833   inline typename CollectiveIterator<T,A>::iterator
01834   CollectiveIterator<T,A>::end()
01835   {
01836     return iterator(map_.end(), map_.end(), index_);
01837   }
01838   
01839   template<typename TG, typename TA>
01840   inline std::ostream& operator<<(std::ostream& os, const RemoteIndex<TG,TA>& index)
01841   {
01842     os<<"[global="<<index.localIndexPair().global()<<", remote attribute="<<index.attribute()<<" local attribute="<<index.localIndexPair().local().attribute()<<"]";
01843     return os;
01844   }
01845   
01846   template<typename T, typename A>
01847   inline std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices)
01848   {
01849     int rank;
01850     MPI_Comm_rank(indices.comm_, &rank);
01851 
01852     typedef typename RemoteIndices<T,A>::RemoteIndexList RList;
01853     typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
01854 
01855     const const_iterator rend = indices.remoteIndices_.end();
01856 
01857     for(const_iterator rindex = indices.remoteIndices_.begin(); rindex!=rend; ++rindex){
01858       os<<rank<<": Prozess "<<rindex->first<<":";
01859       
01860       if(!rindex->second.first->empty()){
01861         os<<" send:";
01862 
01863         const typename RList::const_iterator send= rindex->second.first->end();
01864       
01865         for(typename RList::const_iterator index = rindex->second.first->begin(); 
01866             index != send; ++index)
01867           os<<*index<<" ";
01868         os<<std::endl;
01869       }
01870       if(!rindex->second.second->empty()){
01871         os<<rank<<": Prozess "<<rindex->first<<": "<<"receive: ";
01872         
01873         const typename RList::const_iterator rend= rindex->second.second->end();
01874         
01875         for(typename RList::const_iterator index = rindex->second.second->begin(); 
01876             index != rend; ++index)
01877           os<<*index<<" ";
01878       }
01879       os<<std::endl<<std::flush;
01880     }
01881     return os;
01882   }
01884 }
01885 
01886 #endif
01887 #endif