dune-common  2.3.0
remoteindices.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 // $Id$
4 #ifndef DUNE_REMOTEINDICES_HH
5 #define DUNE_REMOTEINDICES_HH
6 
7 #include "indexset.hh"
8 #include "plocalindex.hh"
11 #include <dune/common/sllist.hh>
14 #include <map>
15 #include <set>
16 #include <utility>
17 #include <iostream>
18 #include <algorithm>
19 #include <iterator>
20 #if HAVE_MPI
21 #include "mpitraits.hh"
22 #include <mpi.h>
23 
24 namespace Dune {
35 
36  template<typename TG, typename TA>
38  {
39  public:
40  inline static MPI_Datatype getType();
41  private:
42  static MPI_Datatype type;
43  };
44 
45 
46  template<typename T, typename A>
47  class RemoteIndices;
48 
49  template<typename T1, typename T2>
50  class RemoteIndex;
51 
52  template<typename T>
53  class IndicesSyncer;
54 
55  template<typename T1, typename T2>
56  std::ostream& operator<<(std::ostream& os, const RemoteIndex<T1,T2>& index);
57 
58 
59  template<typename T, typename A, bool mode>
61 
62 
66  template<typename T1, typename T2>
68  {
69  template<typename T>
70  friend class IndicesSyncer;
71 
72  template<typename T, typename A, typename A1>
73  friend void repairLocalIndexPointers(std::map<int,SLList<std::pair<typename T::GlobalIndex, typename T::LocalIndex::Attribute>,A> >&,
75  const T&);
76 
77  template<typename T, typename A, bool mode>
79 
80  public:
85  typedef T1 GlobalIndex;
94  typedef T2 Attribute;
95 
101 
106  const Attribute attribute() const;
107 
113  const PairType& localIndexPair() const;
114 
118  RemoteIndex();
119 
120 
126  RemoteIndex(const T2& attribute,
127  const PairType* local);
128 
129 
135  RemoteIndex(const T2& attribute);
136 
137  bool operator==(const RemoteIndex& ri) const;
138 
139  bool operator!=(const RemoteIndex& ri) const;
140  private:
142  const PairType* localIndex_;
143 
145  char attribute_;
146  };
147 
148  template<class T, class A>
149  std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices);
150 
151  class InterfaceBuilder;
152 
153  template<class T, class A>
154  class CollectiveIterator;
155 
156  template<class T>
157  class IndicesSyncer;
158 
159  // forward declaration needed for friend declaration.
160  template<typename T1, typename T2>
161  class OwnerOverlapCopyCommunication;
162 
163 
180  template<class T, class A=std::allocator<RemoteIndex<typename T::GlobalIndex,
181  typename T::LocalIndex::Attribute> > >
183  {
184  friend class InterfaceBuilder;
185  friend class IndicesSyncer<T>;
186  template<typename T1, typename A2, typename A1>
187  friend void repairLocalIndexPointers(std::map<int,SLList<std::pair<typename T1::GlobalIndex, typename T1::LocalIndex::Attribute>,A2> >&,
189  const T1&);
190 
191  template<class G, class T1, class T2>
192  friend void fillIndexSetHoles(const G& graph, Dune::OwnerOverlapCopyCommunication<T1,T2>& oocomm);
193  friend std::ostream& operator<<<>(std::ostream&, const RemoteIndices<T>&);
194 
195  public:
196 
200  typedef T ParallelIndexSet;
201 
205 
210 
211 
216 
220  typedef typename LocalIndex::Attribute Attribute;
221 
226 
227 
231  typedef typename A::template rebind<RemoteIndex>::other Allocator;
232 
236 
238  typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
240 
241  typedef typename RemoteIndexMap::const_iterator const_iterator;
242 
260  inline RemoteIndices(const ParallelIndexSet& source, const ParallelIndexSet& destination,
261  const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>(), bool includeSelf=false);
262 
263  RemoteIndices();
264 
272  void setIncludeSelf(bool includeSelf);
273 
290  void setIndexSets(const ParallelIndexSet& source, const ParallelIndexSet& destination,
291  const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>());
292 
293  template<typename C>
294  void setNeighbours(const C& neighbours)
295  {
296  neighbourIds.clear();
297  neighbourIds.insert(neighbours.begin(), neighbours.end());
298 
299  }
300 
301  const std::set<int>& getNeighbours() const
302  {
303  return neighbourIds;
304  }
305 
309  ~RemoteIndices();
310 
320  template<bool ignorePublic>
321  void rebuild();
322 
323  bool operator==(const RemoteIndices& ri);
324 
332  inline bool isSynced() const;
333 
337  inline MPI_Comm communicator() const;
338 
353  template<bool mode, bool send>
355 
362  inline const_iterator find(int proc) const;
363 
368  inline const_iterator begin() const;
369 
374  inline const_iterator end() const;
375 
379  template<bool send>
380  inline CollectiveIteratorT iterator() const;
381 
385  inline void free();
386 
391  inline int neighbours() const;
392 
394  inline const ParallelIndexSet& sourceIndexSet() const;
395 
397  inline const ParallelIndexSet& destinationIndexSet() const;
398 
399  private:
402  {}
403 
405  const ParallelIndexSet* source_;
406 
408  const ParallelIndexSet* target_;
409 
411  MPI_Comm comm_;
412 
415  std::set<int> neighbourIds;
416 
418  const static int commTag_=333;
419 
424  int sourceSeqNo_;
425 
430  int destSeqNo_;
431 
435  bool publicIgnored;
436 
440  bool firstBuild;
441 
442  /*
443  * @brief If true, sending from indices of the processor to other
444  * indices on the same processor is enabled even if the same indexset is used
445  * on both the
446  * sending and receiving side.
447  */
448  bool includeSelf;
449 
451  typedef IndexPair<GlobalIndex, LocalIndex>
452  PairType;
453 
460  RemoteIndexMap remoteIndices_;
461 
472  template<bool ignorePublic>
473  inline void buildRemote(bool includeSelf);
474 
480  inline int noPublic(const ParallelIndexSet& indexSet);
481 
493  template<bool ignorePublic>
494  inline void packEntries(PairType** myPairs, const ParallelIndexSet& indexSet,
495  char* p_out, MPI_Datatype type, int bufferSize,
496  int* position, int n);
497 
511  inline void unpackIndices(RemoteIndexList& remote, int remoteEntries,
512  PairType** local, int localEntries, char* p_in,
513  MPI_Datatype type, int* positon, int bufferSize,
514  bool fromOurself);
515 
516  inline void unpackIndices(RemoteIndexList& send, RemoteIndexList& receive,
517  int remoteEntries, PairType** localSource,
518  int localSourceEntries, PairType** localDest,
519  int localDestEntries, char* p_in,
520  MPI_Datatype type, int* position, int bufferSize);
521 
522  void unpackCreateRemote(char* p_in, PairType** sourcePairs, PairType** DestPairs,
523  int remoteProc, int sourcePublish, int destPublish,
524  int bufferSize, bool sendTwo, bool fromOurSelf=false);
525  };
526 
544  template<class T, class A, bool mode>
546  {
547 
548  template<typename T1, typename A1>
549  friend class RemoteIndices;
550 
551  public:
553  {};
554 
555  enum {
565  };
566 
570  typedef T ParallelIndexSet;
571 
576 
581 
585  typedef typename LocalIndex::Attribute Attribute;
586 
591 
595  typedef A Allocator;
596 
600 
605 
610 
624  void insert(const RemoteIndex& index) throw(InvalidPosition);
625 
626 
641  void insert(const RemoteIndex& index, const GlobalIndex& global) throw(InvalidPosition);
642 
650  bool remove(const GlobalIndex& global) throw(InvalidPosition);
651 
665 
666 
668 
673  RemoteIndexListModifier()
674  : glist_()
675  {}
676 
677  private:
678 
685  RemoteIndexList& rList);
686 
687  typedef SLList<GlobalIndex,Allocator> GlobalList;
688  typedef typename GlobalList::ModifyIterator GlobalModifyIterator;
689  RemoteIndexList* rList_;
690  const ParallelIndexSet* indexSet_;
691  GlobalList glist_;
692  ModifyIterator iter_;
693  GlobalModifyIterator giter_;
694  ConstIterator end_;
695  bool first_;
696  GlobalIndex last_;
697  };
698 
703  template<class T, class A>
705  {
706 
710  typedef T ParallelIndexSet;
711 
715  typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
716 
720  typedef typename ParallelIndexSet::LocalIndex LocalIndex;
721 
725  typedef typename LocalIndex::Attribute Attribute;
726 
729 
731  typedef typename A::template rebind<RemoteIndex>::other Allocator;
732 
735 
737  typedef std::map<int,std::pair<typename RemoteIndexList::const_iterator,
738  const typename RemoteIndexList::const_iterator> >
739  Map;
740 
741  public:
742 
744  typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
746 
752  inline CollectiveIterator(const RemoteIndexMap& map_, bool send);
753 
762  inline void advance(const GlobalIndex& global);
763 
773  inline void advance(const GlobalIndex& global, const Attribute& attribute);
774 
776 
780  inline bool empty();
781 
788  class iterator
789  {
790  public:
791  typedef typename Map::iterator RealIterator;
792  typedef typename Map::iterator ConstRealIterator;
793 
794 
796  iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex& index)
797  : iter_(iter), end_(end), index_(index), hasAttribute(false)
798  {
799  // Move to the first valid entry
800  while(iter_!=end_ && iter_->second.first->localIndexPair().global()!=index_)
801  ++iter_;
802  }
803 
804  iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex index,
805  Attribute attribute)
806  : iter_(iter), end_(end), index_(index), attribute_(attribute), hasAttribute(true)
807  {
808  // Move to the first valid entry or the end
809  while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_
810  || iter_->second.first->localIndexPair().local().attribute()!=attribute))
811  ++iter_;
812  }
814  iterator(const iterator& other)
815  : iter_(other.iter_), end_(other.end_), index_(other.index_)
816  { }
817 
820  {
821  ++iter_;
822  // If entry is not valid move on
823  while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_ ||
824  (hasAttribute &&
825  iter_->second.first->localIndexPair().local().attribute()!=attribute_)))
826  ++iter_;
827  assert(iter_==end_ ||
828  (iter_->second.first->localIndexPair().global()==index_));
829  assert(iter_==end_ || !hasAttribute ||
830  (iter_->second.first->localIndexPair().local().attribute()==attribute_));
831  return *this;
832  }
833 
835  const RemoteIndex& operator*() const
836  {
837  return *(iter_->second.first);
838  }
839 
841  int process() const
842  {
843  return iter_->first;
844  }
845 
847  const RemoteIndex* operator->() const
848  {
849  return iter_->second.first.operator->();
850  }
851 
853  bool operator==(const iterator& other)
854  {
855  return other.iter_==iter_;
856  }
857 
859  bool operator!=(const iterator& other)
860  {
861  return other.iter_!=iter_;
862  }
863 
864  private:
865  iterator();
866 
867  RealIterator iter_;
868  RealIterator end_;
869  GlobalIndex index_;
870  Attribute attribute_;
871  bool hasAttribute;
872  };
873 
874  iterator begin();
875 
876  iterator end();
877 
878  private:
879 
880  Map map_;
881  GlobalIndex index_;
882  Attribute attribute_;
883  bool noattribute;
884  };
885 
886  template<typename TG, typename TA>
888  {
889  if(type==MPI_DATATYPE_NULL) {
890  int length[4];
891  MPI_Aint disp[4];
892  MPI_Datatype types[4] = {MPI_LB, MPITraits<TG>::getType(),
893  MPITraits<ParallelLocalIndex<TA> >::getType(), MPI_UB};
895  length[0]=length[1]=length[2]=length[3]=1;
896  MPI_Address(rep, disp); // lower bound of the datatype
897  MPI_Address(&(rep[0].global_), disp+1);
898  MPI_Address(&(rep[0].local_), disp+2);
899  MPI_Address(rep+1, disp+3); // upper bound of the datatype
900  for(int i=3; i >= 0; --i)
901  disp[i] -= disp[0];
902  MPI_Type_struct(4, length, disp, types, &type);
903  MPI_Type_commit(&type);
904  }
905  return type;
906  }
907 
908  template<typename TG, typename TA>
909  MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::type=MPI_DATATYPE_NULL;
910 
911  template<typename T1, typename T2>
912  RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute, const PairType* local)
913  : localIndex_(local), attribute_(attribute)
914  {}
915 
916  template<typename T1, typename T2>
917  RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute)
918  : localIndex_(0), attribute_(attribute)
919  {}
920 
921  template<typename T1, typename T2>
923  : localIndex_(0), attribute_()
924  {}
925  template<typename T1, typename T2>
926  inline bool RemoteIndex<T1,T2>::operator==(const RemoteIndex& ri) const
927  {
928  return localIndex_==ri.localIndex_ && attribute_==ri.attribute;
929  }
930 
931  template<typename T1, typename T2>
932  inline bool RemoteIndex<T1,T2>::operator!=(const RemoteIndex& ri) const
933  {
934  return localIndex_!=ri.localIndex_ || attribute_!=ri.attribute_;
935  }
936 
937  template<typename T1, typename T2>
938  inline const T2 RemoteIndex<T1,T2>::attribute() const
939  {
940  return T2(attribute_);
941  }
942 
943  template<typename T1, typename T2>
945  {
946  return *localIndex_;
947  }
948 
949  template<typename T, typename A>
951  const ParallelIndexSet& destination,
952  const MPI_Comm& comm,
953  const std::vector<int>& neighbours,
954  bool includeSelf_)
955  : source_(&source), target_(&destination), comm_(comm),
956  sourceSeqNo_(-1), destSeqNo_(-1), publicIgnored(false), firstBuild(true),
957  includeSelf(includeSelf_)
958  {
959  setNeighbours(neighbours);
960  }
961 
962  template<typename T, typename A>
964  {
965  includeSelf=b;
966  }
967 
968  template<typename T, typename A>
970  : source_(0), target_(0), sourceSeqNo_(-1),
971  destSeqNo_(-1), publicIgnored(false), firstBuild(true)
972  {}
973 
974  template<class T, typename A>
976  const ParallelIndexSet& destination,
977  const MPI_Comm& comm,
978  const std::vector<int>& neighbours)
979  {
980  free();
981  source_ = &source;
982  target_ = &destination;
983  comm_ = comm;
984  firstBuild = true;
985  setNeighbours(neighbours);
986  }
987 
988  template<typename T, typename A>
991  {
992  return *source_;
993  }
994 
995 
996  template<typename T, typename A>
999  {
1000  return *target_;
1001  }
1002 
1003 
1004  template<typename T, typename A>
1006  {
1007  free();
1008  }
1009 
1010  template<typename T, typename A>
1011  template<bool ignorePublic>
1013  const ParallelIndexSet& indexSet,
1014  char* p_out, MPI_Datatype type,
1015  int bufferSize,
1016  int *position, int n)
1017  {
1018  // fill with own indices
1021  const const_iterator end = indexSet.end();
1022 
1023  //Now pack the source indices
1024  int i=0;
1025  for(const_iterator index = indexSet.begin(); index != end; ++index)
1026  if(ignorePublic || index->local().isPublic()) {
1027 
1028  MPI_Pack(const_cast<PairType*>(&(*index)), 1,
1029  type,
1030  p_out, bufferSize, position, comm_);
1031  pairs[i++] = const_cast<PairType*>(&(*index));
1032 
1033  }
1034  assert(i==n);
1035  }
1036 
1037  template<typename T, typename A>
1038  inline int RemoteIndices<T,A>::noPublic(const ParallelIndexSet& indexSet)
1039  {
1040  typedef typename ParallelIndexSet::const_iterator const_iterator;
1041 
1042  int noPublic=0;
1043 
1044  const const_iterator end=indexSet.end();
1045  for(const_iterator index=indexSet.begin(); index!=end; ++index)
1046  if(index->local().isPublic())
1047  noPublic++;
1048 
1049  return noPublic;
1050 
1051  }
1052 
1053 
1054  template<typename T, typename A>
1055  inline void RemoteIndices<T,A>::unpackCreateRemote(char* p_in, PairType** sourcePairs,
1056  PairType** destPairs, int remoteProc,
1057  int sourcePublish, int destPublish,
1058  int bufferSize, bool sendTwo,
1059  bool fromOurSelf)
1060  {
1061 
1062  // unpack the number of indices we received
1063  int noRemoteSource=-1, noRemoteDest=-1;
1064  char twoIndexSets=0;
1065  int position=0;
1066  // Did we receive two index sets?
1067  MPI_Unpack(p_in, bufferSize, &position, &twoIndexSets, 1, MPI_CHAR, comm_);
1068  // The number of source indices received
1069  MPI_Unpack(p_in, bufferSize, &position, &noRemoteSource, 1, MPI_INT, comm_);
1070  // The number of destination indices received
1071  MPI_Unpack(p_in, bufferSize, &position, &noRemoteDest, 1, MPI_INT, comm_);
1072 
1073 
1074  // Indices for which we receive
1075  RemoteIndexList* receive= new RemoteIndexList();
1076  // Indices for which we send
1077  RemoteIndexList* send=0;
1078 
1079  MPI_Datatype type= MPITraits<PairType>::getType();
1080 
1081  if(!twoIndexSets) {
1082  if(sendTwo) {
1083  send = new RemoteIndexList();
1084  // Create both remote index sets simultaneously
1085  unpackIndices(*send, *receive, noRemoteSource, sourcePairs, sourcePublish,
1086  destPairs, destPublish, p_in, type, &position, bufferSize);
1087  }else{
1088  // we only need one list
1089  unpackIndices(*receive, noRemoteSource, sourcePairs, sourcePublish,
1090  p_in, type, &position, bufferSize, fromOurSelf);
1091  send=receive;
1092  }
1093  }else{
1094 
1095  int oldPos=position;
1096  // Two index sets received
1097  unpackIndices(*receive, noRemoteSource, destPairs, destPublish,
1098  p_in, type, &position, bufferSize, fromOurSelf);
1099  if(!sendTwo)
1100  //unpack source entries again as destination entries
1101  position=oldPos;
1102 
1103  send = new RemoteIndexList();
1104  unpackIndices(*send, noRemoteDest, sourcePairs, sourcePublish,
1105  p_in, type, &position, bufferSize, fromOurSelf);
1106  }
1107 
1108  if(receive->empty() && send->empty()) {
1109  if(send==receive) {
1110  delete send;
1111  }else{
1112  delete send;
1113  delete receive;
1114  }
1115  }else{
1116  remoteIndices_.insert(std::make_pair(remoteProc,
1117  std::make_pair(send,receive)));
1118  }
1119  }
1120 
1121 
1122  template<typename T, typename A>
1123  template<bool ignorePublic>
1124  inline void RemoteIndices<T,A>::buildRemote(bool includeSelf)
1125  {
1126  // Processor configuration
1127  int rank, procs;
1128  MPI_Comm_rank(comm_, &rank);
1129  MPI_Comm_size(comm_, &procs);
1130 
1131  // number of local indices to publish
1132  // The indices of the destination will be send.
1133  int sourcePublish, destPublish;
1134 
1135  // Do we need to send two index sets?
1136  char sendTwo = (source_ != target_);
1137 
1138  if(procs==1 && !(sendTwo || includeSelf))
1139  // Nothing to communicate
1140  return;
1141 
1142  sourcePublish = (ignorePublic) ? source_->size() : noPublic(*source_);
1143 
1144  if(sendTwo)
1145  destPublish = (ignorePublic) ? target_->size() : noPublic(*target_);
1146  else
1147  // we only need to send one set of indices
1148  destPublish = 0;
1149 
1150  int maxPublish, publish=sourcePublish+destPublish;
1151 
1152  // Calucate maximum number of indices send
1153  MPI_Allreduce(&publish, &maxPublish, 1, MPI_INT, MPI_MAX, comm_);
1154 
1155  // allocate buffers
1156  typedef IndexPair<GlobalIndex,LocalIndex> PairType;
1157 
1158  PairType** destPairs;
1159  PairType** sourcePairs = new PairType*[sourcePublish>0 ? sourcePublish : 1];
1160 
1161  if(sendTwo)
1162  destPairs = new PairType*[destPublish>0 ? destPublish : 1];
1163  else
1164  destPairs=sourcePairs;
1165 
1166  char** buffer = new char*[2];
1167  int bufferSize;
1168  int position=0;
1169  int intSize;
1170  int charSize;
1171 
1172  // calculate buffer size
1173  MPI_Datatype type = MPITraits<PairType>::getType();
1174 
1175  MPI_Pack_size(maxPublish, type, comm_,
1176  &bufferSize);
1177  MPI_Pack_size(1, MPI_INT, comm_,
1178  &intSize);
1179  MPI_Pack_size(1, MPI_CHAR, comm_,
1180  &charSize);
1181  // Our message will contain the following:
1182  // a bool wether two index sets where sent
1183  // the size of the source and the dest indexset,
1184  // then the source and destination indices
1185  bufferSize += 2 * intSize + charSize;
1186 
1187  if(bufferSize<=0) bufferSize=1;
1188 
1189  buffer[0] = new char[bufferSize];
1190  buffer[1] = new char[bufferSize];
1191 
1192 
1193  // pack entries into buffer[0], p_out below!
1194  MPI_Pack(&sendTwo, 1, MPI_CHAR, buffer[0], bufferSize, &position,
1195  comm_);
1196 
1197  // The number of indices we send for each index set
1198  MPI_Pack(&sourcePublish, 1, MPI_INT, buffer[0], bufferSize, &position,
1199  comm_);
1200  MPI_Pack(&destPublish, 1, MPI_INT, buffer[0], bufferSize, &position,
1201  comm_);
1202 
1203  // Now pack the source indices and setup the destination pairs
1204  packEntries<ignorePublic>(sourcePairs, *source_, buffer[0], type,
1205  bufferSize, &position, sourcePublish);
1206  // If necessary send the dest indices and setup the source pairs
1207  if(sendTwo)
1208  packEntries<ignorePublic>(destPairs, *target_, buffer[0], type,
1209  bufferSize, &position, destPublish);
1210 
1211 
1212  // Update remote indices for ourself
1213  if(sendTwo|| includeSelf)
1214  unpackCreateRemote(buffer[0], sourcePairs, destPairs, rank, sourcePublish,
1215  destPublish, bufferSize, sendTwo, includeSelf);
1216 
1217  neighbourIds.erase(rank);
1218 
1219  if(neighbourIds.size()==0)
1220  {
1221  Dune::dvverb<<rank<<": Sending messages in a ring"<<std::endl;
1222  // send messages in ring
1223  for(int proc=1; proc<procs; proc++) {
1224  // pointers to the current input and output buffers
1225  char* p_out = buffer[1-(proc%2)];
1226  char* p_in = buffer[proc%2];
1227 
1228  MPI_Status status;
1229  if(rank%2==0) {
1230  MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
1231  commTag_, comm_);
1232  MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
1233  commTag_, comm_, &status);
1234  }else{
1235  MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
1236  commTag_, comm_, &status);
1237  MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
1238  commTag_, comm_);
1239  }
1240 
1241 
1242  // The process these indices are from
1243  int remoteProc = (rank+procs-proc)%procs;
1244 
1245  unpackCreateRemote(p_in, sourcePairs, destPairs, remoteProc, sourcePublish,
1246  destPublish, bufferSize, sendTwo);
1247 
1248  }
1249 
1250  }
1251  else
1252  {
1253  MPI_Request* requests=new MPI_Request[neighbourIds.size()];
1254  MPI_Request* req=requests;
1255 
1256  typedef typename std::set<int>::size_type size_type;
1257  size_type noNeighbours=neighbourIds.size();
1258 
1259  // setup sends
1260  for(std::set<int>::iterator neighbour=neighbourIds.begin();
1261  neighbour!= neighbourIds.end(); ++neighbour) {
1262  // Only send the information to the neighbouring processors
1263  MPI_Issend(buffer[0], position , MPI_PACKED, *neighbour, commTag_, comm_, req++);
1264  }
1265 
1266  //Test for received messages
1267 
1268  for(size_type received=0; received <noNeighbours; ++received)
1269  {
1270  MPI_Status status;
1271  // probe for next message
1272  MPI_Probe(MPI_ANY_SOURCE, commTag_, comm_, &status);
1273  int remoteProc=status.MPI_SOURCE;
1274  int size;
1275  MPI_Get_count(&status, MPI_PACKED, &size);
1276  // receive message
1277  MPI_Recv(buffer[1], size, MPI_PACKED, remoteProc,
1278  commTag_, comm_, &status);
1279 
1280  unpackCreateRemote(buffer[1], sourcePairs, destPairs, remoteProc, sourcePublish,
1281  destPublish, bufferSize, sendTwo);
1282  }
1283  // wait for completion of pending requests
1284  MPI_Status* statuses = new MPI_Status[neighbourIds.size()];
1285 
1286  if(MPI_ERR_IN_STATUS==MPI_Waitall(neighbourIds.size(), requests, statuses)) {
1287  for(size_type i=0; i < neighbourIds.size(); ++i)
1288  if(statuses[i].MPI_ERROR!=MPI_SUCCESS) {
1289  std::cerr<<rank<<": MPI_Error occurred while receiving message."<<std::endl;
1290  MPI_Abort(comm_, 999);
1291  }
1292  }
1293  delete[] requests;
1294  delete[] statuses;
1295  }
1296 
1297 
1298  // delete allocated memory
1299  if(destPairs!=sourcePairs)
1300  delete[] destPairs;
1301 
1302  delete[] sourcePairs;
1303  delete[] buffer[0];
1304  delete[] buffer[1];
1305  delete[] buffer;
1306  }
1307 
1308  template<typename T, typename A>
1309  inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& remote,
1310  int remoteEntries,
1311  PairType** local,
1312  int localEntries,
1313  char* p_in,
1314  MPI_Datatype type,
1315  int* position,
1316  int bufferSize,
1317  bool fromOurSelf)
1318  {
1319  if(remoteEntries==0)
1320  return;
1321 
1322  PairType index(-1);
1323  MPI_Unpack(p_in, bufferSize, position, &index, 1,
1324  type, comm_);
1325  GlobalIndex oldGlobal=index.global();
1326  int n_in=0, localIndex=0;
1327 
1328  //Check if we know the global index
1329  while(localIndex<localEntries) {
1330  if(local[localIndex]->global()==index.global()) {
1331  int oldLocalIndex=localIndex;
1332 
1333  while(localIndex<localEntries &&
1334  local[localIndex]->global()==index.global()) {
1335  if(!fromOurSelf || index.local().attribute() !=
1336  local[localIndex]->local().attribute())
1337  // if index is from us it has to have a different attribute
1338  remote.push_back(RemoteIndex(index.local().attribute(),
1339  local[localIndex]));
1340  localIndex++;
1341  }
1342 
1343  // unpack next remote index
1344  if((++n_in) < remoteEntries) {
1345  MPI_Unpack(p_in, bufferSize, position, &index, 1,
1346  type, comm_);
1347  if(index.global()==oldGlobal)
1348  // Restart comparison for the same global indices
1349  localIndex=oldLocalIndex;
1350  else
1351  oldGlobal=index.global();
1352  }else{
1353  // No more received indices
1354  break;
1355  }
1356  continue;
1357  }
1358 
1359  if (local[localIndex]->global()<index.global()) {
1360  // compare with next entry in our list
1361  ++localIndex;
1362  }else{
1363  // We do not know the index, unpack next
1364  if((++n_in) < remoteEntries) {
1365  MPI_Unpack(p_in, bufferSize, position, &index, 1,
1366  type, comm_);
1367  oldGlobal=index.global();
1368  }else
1369  // No more received indices
1370  break;
1371  }
1372  }
1373 
1374  // Unpack the other received indices without doing anything
1375  while(++n_in < remoteEntries)
1376  MPI_Unpack(p_in, bufferSize, position, &index, 1,
1377  type, comm_);
1378  }
1379 
1380 
1381  template<typename T, typename A>
1382  inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& send,
1383  RemoteIndexList& receive,
1384  int remoteEntries,
1385  PairType** localSource,
1386  int localSourceEntries,
1387  PairType** localDest,
1388  int localDestEntries,
1389  char* p_in,
1390  MPI_Datatype type,
1391  int* position,
1392  int bufferSize)
1393  {
1394  int n_in=0, sourceIndex=0, destIndex=0;
1395 
1396  //Check if we know the global index
1397  while(n_in<remoteEntries && (sourceIndex<localSourceEntries || destIndex<localDestEntries)) {
1398  // Unpack next index
1399  PairType index;
1400  MPI_Unpack(p_in, bufferSize, position, &index, 1,
1401  type, comm_);
1402  n_in++;
1403 
1404  // Advance until global index in localSource and localDest are >= than the one in the unpacked index
1405  while(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()<index.global())
1406  sourceIndex++;
1407 
1408  while(destIndex<localDestEntries && localDest[destIndex]->global()<index.global())
1409  destIndex++;
1410 
1411  // Add a remote index if we found the global index.
1412  if(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()==index.global())
1413  send.push_back(RemoteIndex(index.local().attribute(),
1414  localSource[sourceIndex]));
1415 
1416  if(destIndex < localDestEntries && localDest[destIndex]->global() == index.global())
1417  receive.push_back(RemoteIndex(index.local().attribute(),
1418  localDest[sourceIndex]));
1419  }
1420 
1421  }
1422 
1423  template<typename T, typename A>
1425  {
1426  typedef typename RemoteIndexMap::iterator Iterator;
1427  Iterator lend = remoteIndices_.end();
1428  for(Iterator lists=remoteIndices_.begin(); lists != lend; ++lists) {
1429  if(lists->second.first==lists->second.second) {
1430  // there is only one remote index list.
1431  delete lists->second.first;
1432  }else{
1433  delete lists->second.first;
1434  delete lists->second.second;
1435  }
1436  }
1437  remoteIndices_.clear();
1438  firstBuild=true;
1439  }
1440 
1441  template<typename T, typename A>
1443  {
1444  return remoteIndices_.size();
1445  }
1446 
1447  template<typename T, typename A>
1448  template<bool ignorePublic>
1450  {
1451  // Test wether a rebuild is Needed.
1452  if(firstBuild ||
1453  ignorePublic!=publicIgnored || !
1454  isSynced()) {
1455  free();
1456 
1457  buildRemote<ignorePublic>(includeSelf);
1458 
1459  sourceSeqNo_ = source_->seqNo();
1460  destSeqNo_ = target_->seqNo();
1461  firstBuild=false;
1462  publicIgnored=ignorePublic;
1463  }
1464 
1465 
1466  }
1467 
1468  template<typename T, typename A>
1469  inline bool RemoteIndices<T,A>::isSynced() const
1470  {
1471  return sourceSeqNo_==source_->seqNo() && destSeqNo_ ==target_->seqNo();
1472  }
1473 
1474  template<typename T, typename A>
1475  template<bool mode, bool send>
1477  {
1478 
1479  // The user are on their own now!
1480  // We assume they know what they are doing and just set the
1481  // remote indices to synced status.
1482  sourceSeqNo_ = source_->seqNo();
1483  destSeqNo_ = target_->seqNo();
1484 
1485  typename RemoteIndexMap::iterator found = remoteIndices_.find(process);
1486 
1487  if(found == remoteIndices_.end())
1488  {
1489  if(source_ != target_)
1490  found = remoteIndices_.insert(found, std::make_pair(process,
1491  std::make_pair(new RemoteIndexList(),
1492  new RemoteIndexList())));
1493  else{
1494  RemoteIndexList* rlist = new RemoteIndexList();
1495  found = remoteIndices_.insert(found,
1496  std::make_pair(process,
1497  std::make_pair(rlist, rlist)));
1498  }
1499  }
1500 
1501  firstBuild = false;
1502 
1503  if(send)
1504  return RemoteIndexListModifier<T,A,mode>(*source_, *(found->second.first));
1505  else
1506  return RemoteIndexListModifier<T,A,mode>(*target_, *(found->second.second));
1507  }
1508 
1509  template<typename T, typename A>
1510  inline typename RemoteIndices<T,A>::const_iterator
1512  {
1513  return remoteIndices_.find(proc);
1514  }
1515 
1516  template<typename T, typename A>
1517  inline typename RemoteIndices<T,A>::const_iterator
1519  {
1520  return remoteIndices_.begin();
1521  }
1522 
1523  template<typename T, typename A>
1524  inline typename RemoteIndices<T,A>::const_iterator
1526  {
1527  return remoteIndices_.end();
1528  }
1529 
1530 
1531  template<typename T, typename A>
1533  {
1534  if(neighbours()!=ri.neighbours())
1535  return false;
1536 
1537  typedef RemoteIndexList RList;
1538  typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
1539 
1540  const const_iterator rend = remoteIndices_.end();
1541 
1542  for(const_iterator rindex = remoteIndices_.begin(), rindex1=ri.remoteIndices_.begin(); rindex!=rend; ++rindex, ++rindex1) {
1543  if(rindex->first != rindex1->first)
1544  return false;
1545  if(*(rindex->second.first) != *(rindex1->second.first))
1546  return false;
1547  if(*(rindex->second.second) != *(rindex1->second.second))
1548  return false;
1549  }
1550  return true;
1551  }
1552 
1553  template<class T, class A, bool mode>
1555  RemoteIndexList& rList)
1556  : rList_(&rList), indexSet_(&indexSet), iter_(rList.beginModify()), end_(rList.end()), first_(true)
1557  {
1558  if(MODIFYINDEXSET) {
1559  assert(indexSet_);
1560  for(ConstIterator iter=iter_; iter != end_; ++iter)
1561  glist_.push_back(iter->localIndexPair().global());
1562  giter_ = glist_.beginModify();
1563  }
1564  }
1565 
1566  template<typename T, typename A, bool mode>
1568  : rList_(other.rList_), indexSet_(other.indexSet_),
1569  glist_(other.glist_), iter_(other.iter_), giter_(other.giter_), end_(other.end_),
1570  first_(other.first_), last_(other.last_)
1571  {}
1572 
1573  template<typename T, typename A, bool mode>
1575  {
1576  if(MODIFYINDEXSET) {
1577  // repair pointers to local index set.
1578 #ifdef DUNE_ISTL_WITH_CHECKING
1579  if(indexSet_->state()!=GROUND)
1580  DUNE_THROW(InvalidIndexSetState, "Index has to be in ground mode for repairing pointers to indices");
1581 #endif
1582  typedef typename ParallelIndexSet::const_iterator IndexIterator;
1583  typedef typename GlobalList::const_iterator GlobalIterator;
1584  typedef typename RemoteIndexList::iterator Iterator;
1585  GlobalIterator giter = glist_.begin();
1586  IndexIterator index = indexSet_->begin();
1587 
1588  for(Iterator iter=rList_->begin(); iter != end_; ++iter) {
1589  while(index->global()<*giter) {
1590  ++index;
1591 #ifdef DUNE_ISTL_WITH_CHECKING
1592  if(index == indexSet_->end())
1593  DUNE_THROW(InvalidPosition, "No such global index in set!");
1594 #endif
1595  }
1596 
1597 #ifdef DUNE_ISTL_WITH_CHECKING
1598  if(index->global() != *giter)
1599  DUNE_THROW(InvalidPosition, "No such global index in set!");
1600 #endif
1601  iter->localIndex_ = &(*index);
1602  }
1603  }
1604  }
1605 
1606  template<typename T, typename A, bool mode>
1608  {
1609  dune_static_assert(!mode,"Not allowed if the mode indicates that new indices"
1610  "might be added to the underlying index set. Use "
1611  "insert(const RemoteIndex&, const GlobalIndex&) instead");
1612 
1613 #ifdef DUNE_ISTL_WITH_CHECKING
1614  if(!first_ && index.localIndexPair().global()<last_)
1615  DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
1616 #endif
1617  // Move to the correct position
1618  while(iter_ != end_ && iter_->localIndexPair().global() < index.localIndexPair().global()) {
1619  ++iter_;
1620  }
1621 
1622  // No duplicate entries allowed
1623  assert(iter_==end_ || iter_->localIndexPair().global() != index.localIndexPair().global());
1624  iter_.insert(index);
1625  last_ = index.localIndexPair().global();
1626  first_ = false;
1627  }
1628 
1629  template<typename T, typename A, bool mode>
1631  {
1632  dune_static_assert(mode,"Not allowed if the mode indicates that no new indices"
1633  "might be added to the underlying index set. Use "
1634  "insert(const RemoteIndex&) instead");
1635 #ifdef DUNE_ISTL_WITH_CHECKING
1636  if(!first_ && global<last_)
1637  DUNE_THROW(InvalidPosition, "Modification of remote indices have to occur with ascending global index.");
1638 #endif
1639  // Move to the correct position
1640  while(iter_ != end_ && *giter_ < global) {
1641  ++giter_;
1642  ++iter_;
1643  }
1644 
1645  // No duplicate entries allowed
1646  assert(iter_->localIndexPair().global() != global);
1647  iter_.insert(index);
1648  giter_.insert(global);
1649 
1650  last_ = global;
1651  first_ = false;
1652  }
1653 
1654  template<typename T, typename A, bool mode>
1656  {
1657 #ifdef DUNE_ISTL_WITH_CHECKING
1658  if(!first_ && global<last_)
1659  DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
1660 #endif
1661 
1662  bool found= false;
1663 
1664  if(MODIFYINDEXSET) {
1665  // Move to the correct position
1666  while(iter_!=end_ && *giter_< global) {
1667  ++giter_;
1668  ++iter_;
1669  }
1670  if(*giter_ == global) {
1671  giter_.remove();
1672  iter_.remove();
1673  found=true;
1674  }
1675  }else{
1676  while(iter_!=end_ && iter_->localIndexPair().global() < global)
1677  ++iter_;
1678 
1679  if(iter_->localIndexPair().global()==global) {
1680  iter_.remove();
1681  found = true;
1682  }
1683  }
1684 
1685  last_ = global;
1686  first_ = false;
1687  return found;
1688  }
1689 
1690  template<typename T, typename A>
1691  template<bool send>
1693  {
1694  return CollectiveIterator<T,A>(remoteIndices_, send);
1695  }
1696 
1697  template<typename T, typename A>
1698  inline MPI_Comm RemoteIndices<T,A>::communicator() const
1699  {
1700  return comm_;
1701 
1702  }
1703 
1704  template<typename T, typename A>
1706  {
1707  typedef typename RemoteIndexMap::const_iterator const_iterator;
1708 
1709  const const_iterator end=pmap.end();
1710  for(const_iterator process=pmap.begin(); process != end; ++process) {
1711  const RemoteIndexList* list = send ? process->second.first : process->second.second;
1712  typedef typename RemoteIndexList::const_iterator iterator;
1713  map_.insert(std::make_pair(process->first,
1714  std::pair<iterator, const iterator>(list->begin(), list->end())));
1715  }
1716  }
1717 
1718  template<typename T, typename A>
1719  inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index)
1720  {
1721  typedef typename Map::iterator iterator;
1722  typedef typename Map::const_iterator const_iterator;
1723  const const_iterator end = map_.end();
1724 
1725  for(iterator iter = map_.begin(); iter != end;) {
1726  // Step the iterator until we are >= index
1727  typename RemoteIndexList::const_iterator current = iter->second.first;
1728  typename RemoteIndexList::const_iterator rend = iter->second.second;
1729  RemoteIndex remoteIndex;
1730  if(current != rend)
1731  remoteIndex = *current;
1732 
1733  while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
1734  ++(iter->second.first);
1735 
1736  // erase from the map if there are no more entries.
1737  if(iter->second.first == iter->second.second)
1738  map_.erase(iter++);
1739  else{
1740  ++iter;
1741  }
1742  }
1743  index_=index;
1744  noattribute=true;
1745  }
1746 
1747  template<typename T, typename A>
1748  inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index,
1749  const Attribute& attribute)
1750  {
1751  typedef typename Map::iterator iterator;
1752  typedef typename Map::const_iterator const_iterator;
1753  const const_iterator end = map_.end();
1754 
1755  for(iterator iter = map_.begin(); iter != end;) {
1756  // Step the iterator until we are >= index
1757  typename RemoteIndexList::const_iterator current = iter->second.first;
1758  typename RemoteIndexList::const_iterator rend = iter->second.second;
1759  RemoteIndex remoteIndex;
1760  if(current != rend)
1761  remoteIndex = *current;
1762 
1763  // Move to global index or bigger
1764  while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
1765  ++(iter->second.first);
1766 
1767  // move to attribute or bigger
1768  while(iter->second.first!=iter->second.second
1769  && iter->second.first->localIndexPair().global()==index
1770  && iter->second.first->localIndexPair().local().attribute()<attribute)
1771  ++(iter->second.first);
1772 
1773  // erase from the map if there are no more entries.
1774  if(iter->second.first == iter->second.second)
1775  map_.erase(iter++);
1776  else{
1777  ++iter;
1778  }
1779  }
1780  index_=index;
1781  attribute_=attribute;
1782  noattribute=false;
1783  }
1784 
1785  template<typename T, typename A>
1787  {
1788  typedef typename Map::iterator iterator;
1789  typedef typename Map::const_iterator const_iterator;
1790  const const_iterator end = map_.end();
1791 
1792  for(iterator iter = map_.begin(); iter != end;) {
1793  // Step the iterator until we are >= index
1794  typename RemoteIndexList::const_iterator current = iter->second.first;
1795  typename RemoteIndexList::const_iterator rend = iter->second.second;
1796 
1797  // move all iterators pointing to the current global index to next value
1798  if(iter->second.first->localIndexPair().global()==index_ &&
1799  (noattribute || iter->second.first->localIndexPair().local().attribute() == attribute_))
1800  ++(iter->second.first);
1801 
1802  // erase from the map if there are no more entries.
1803  if(iter->second.first == iter->second.second)
1804  map_.erase(iter++);
1805  else{
1806  ++iter;
1807  }
1808  }
1809  return *this;
1810  }
1811 
1812  template<typename T, typename A>
1814  {
1815  return map_.empty();
1816  }
1817 
1818  template<typename T, typename A>
1819  inline typename CollectiveIterator<T,A>::iterator
1821  {
1822  if(noattribute)
1823  return iterator(map_.begin(), map_.end(), index_);
1824  else
1825  return iterator(map_.begin(), map_.end(), index_,
1826  attribute_);
1827  }
1828 
1829  template<typename T, typename A>
1830  inline typename CollectiveIterator<T,A>::iterator
1832  {
1833  return iterator(map_.end(), map_.end(), index_);
1834  }
1835 
1836  template<typename TG, typename TA>
1837  inline std::ostream& operator<<(std::ostream& os, const RemoteIndex<TG,TA>& index)
1838  {
1839  os<<"[global="<<index.localIndexPair().global()<<", remote attribute="<<index.attribute()<<" local attribute="<<index.localIndexPair().local().attribute()<<"]";
1840  return os;
1841  }
1842 
1843  template<typename T, typename A>
1844  inline std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices)
1845  {
1846  int rank;
1847  MPI_Comm_rank(indices.comm_, &rank);
1848 
1849  typedef typename RemoteIndices<T,A>::RemoteIndexList RList;
1850  typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
1851 
1852  const const_iterator rend = indices.remoteIndices_.end();
1853 
1854  for(const_iterator rindex = indices.remoteIndices_.begin(); rindex!=rend; ++rindex) {
1855  os<<rank<<": Prozess "<<rindex->first<<":";
1856 
1857  if(!rindex->second.first->empty()) {
1858  os<<" send:";
1859 
1860  const typename RList::const_iterator send= rindex->second.first->end();
1861 
1862  for(typename RList::const_iterator index = rindex->second.first->begin();
1863  index != send; ++index)
1864  os<<*index<<" ";
1865  os<<std::endl;
1866  }
1867  if(!rindex->second.second->empty()) {
1868  os<<rank<<": Prozess "<<rindex->first<<": "<<"receive: ";
1869 
1870  const typename RList::const_iterator rend= rindex->second.second->end();
1871 
1872  for(typename RList::const_iterator index = rindex->second.second->begin();
1873  index != rend; ++index)
1874  os<<*index<<" ";
1875  }
1876  os<<std::endl<<std::flush;
1877  }
1878  return os;
1879  }
1881 }
1882 
1883 #endif
1884 #endif