Dune Core Modules (2.6.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 #ifndef DUNE_REMOTEINDICES_HH
4 #define DUNE_REMOTEINDICES_HH
5 
6 #if HAVE_MPI
7 
8 #include <cassert>
9 #include <iostream>
10 #include <ostream>
11 #include <map>
12 #include <memory>
13 #include <set>
14 #include <utility>
15 #include <vector>
16 
17 #include <mpi.h>
18 
23 #include <dune/common/sllist.hh>
25 
26 namespace Dune {
38  template<typename TG, typename TA>
40  {
41  public:
42  inline static MPI_Datatype getType();
43  private:
44  static MPI_Datatype type;
45  };
46 
47 
48  template<typename T, typename A>
49  class RemoteIndices;
50 
51  template<typename T1, typename T2>
52  class RemoteIndex;
53 
54  // forward declaration needed for friend declaration.
55  template<typename T>
56  class IndicesSyncer;
57 
58  template<typename T1, typename T2>
59  std::ostream& operator<<(std::ostream& os, const RemoteIndex<T1,T2>& index);
60 
61 
62  template<typename T, typename A, bool mode>
64 
65 
69  template<typename T1, typename T2>
71  {
72  template<typename T>
73  friend class IndicesSyncer;
74 
75  template<typename T, typename A, typename A1>
76  friend void repairLocalIndexPointers(std::map<int,SLList<std::pair<typename T::GlobalIndex, typename T::LocalIndex::Attribute>,A> >&,
78  const T&);
79 
80  template<typename T, typename A, bool mode>
81  friend class RemoteIndexListModifier;
82 
83  public:
88  typedef T1 GlobalIndex;
97  typedef T2 Attribute;
98 
104 
109  const Attribute attribute() const;
110 
116  const PairType& localIndexPair() const;
117 
121  RemoteIndex();
122 
123 
129  RemoteIndex(const T2& attribute,
130  const PairType* local);
131 
132 
138  RemoteIndex(const T2& attribute);
139 
140  bool operator==(const RemoteIndex& ri) const;
141 
142  bool operator!=(const RemoteIndex& ri) const;
143  private:
145  const PairType* localIndex_;
146 
148  char attribute_;
149  };
150 
151  template<class T, class A>
152  std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices);
153 
154  class InterfaceBuilder;
155 
156  template<class T, class A>
157  class CollectiveIterator;
158 
159  // forward declaration needed for friend declaration.
160  template<class T>
161  class IndicesSyncer;
162 
163  // forward declaration needed for friend declaration.
164  template<typename T1, typename T2>
166 
167 
184  template<class T, class A=std::allocator<RemoteIndex<typename T::GlobalIndex,
185  typename T::LocalIndex::Attribute> > >
187  {
188  friend class InterfaceBuilder;
189  friend class IndicesSyncer<T>;
190  template<typename T1, typename A2, typename A1>
191  friend void repairLocalIndexPointers(std::map<int,SLList<std::pair<typename T1::GlobalIndex, typename T1::LocalIndex::Attribute>,A2> >&,
193  const T1&);
194 
195  template<class G, class T1, class T2>
196  friend void fillIndexSetHoles(const G& graph, Dune::OwnerOverlapCopyCommunication<T1,T2>& oocomm);
197  friend std::ostream& operator<<<>(std::ostream&, const RemoteIndices<T>&);
198 
199  public:
200 
204  typedef T ParallelIndexSet;
205 
209 
214 
215 
220 
224  typedef typename LocalIndex::Attribute Attribute;
225 
230 
231 
235  typedef typename A::template rebind<RemoteIndex>::other Allocator;
236 
240 
242  typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
244 
245  typedef typename RemoteIndexMap::const_iterator const_iterator;
246 
264  inline RemoteIndices(const ParallelIndexSet& source, const ParallelIndexSet& destination,
265  const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>(), bool includeSelf=false);
266 
267  RemoteIndices();
268 
276  void setIncludeSelf(bool includeSelf);
277 
294  void setIndexSets(const ParallelIndexSet& source, const ParallelIndexSet& destination,
295  const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>());
296 
297  template<typename C>
298  void setNeighbours(const C& neighbours)
299  {
300  neighbourIds.clear();
301  neighbourIds.insert(neighbours.begin(), neighbours.end());
302 
303  }
304 
305  const std::set<int>& getNeighbours() const
306  {
307  return neighbourIds;
308  }
309 
314 
324  template<bool ignorePublic>
325  void rebuild();
326 
327  bool operator==(const RemoteIndices& ri);
328 
336  inline bool isSynced() const;
337 
341  inline MPI_Comm communicator() const;
342 
357  template<bool mode, bool send>
359 
366  inline const_iterator find(int proc) const;
367 
372  inline const_iterator begin() const;
373 
378  inline const_iterator end() const;
379 
383  template<bool send>
385 
389  inline void free();
390 
395  inline int neighbours() const;
396 
398  inline const ParallelIndexSet& sourceIndexSet() const;
399 
401  inline const ParallelIndexSet& destinationIndexSet() const;
402 
403  private:
406  {}
407 
409  const ParallelIndexSet* source_;
410 
412  const ParallelIndexSet* target_;
413 
415  MPI_Comm comm_;
416 
419  std::set<int> neighbourIds;
420 
422  const static int commTag_=333;
423 
428  int sourceSeqNo_;
429 
434  int destSeqNo_;
435 
439  bool publicIgnored;
440 
444  bool firstBuild;
445 
446  /*
447  * @brief If true, sending from indices of the processor to other
448  * indices on the same processor is enabled even if the same indexset is used
449  * on both the
450  * sending and receiving side.
451  */
452  bool includeSelf;
453 
456  PairType;
457 
464  RemoteIndexMap remoteIndices_;
465 
476  template<bool ignorePublic>
477  inline void buildRemote(bool includeSelf);
478 
484  inline int noPublic(const ParallelIndexSet& indexSet);
485 
497  template<bool ignorePublic>
498  inline void packEntries(PairType** myPairs, const ParallelIndexSet& indexSet,
499  char* p_out, MPI_Datatype type, int bufferSize,
500  int* position, int n);
501 
515  inline void unpackIndices(RemoteIndexList& remote, int remoteEntries,
516  PairType** local, int localEntries, char* p_in,
517  MPI_Datatype type, int* position, int bufferSize,
518  bool fromOurself);
519 
520  inline void unpackIndices(RemoteIndexList& send, RemoteIndexList& receive,
521  int remoteEntries, PairType** localSource,
522  int localSourceEntries, PairType** localDest,
523  int localDestEntries, char* p_in,
524  MPI_Datatype type, int* position, int bufferSize);
525 
526  void unpackCreateRemote(char* p_in, PairType** sourcePairs, PairType** DestPairs,
527  int remoteProc, int sourcePublish, int destPublish,
528  int bufferSize, bool sendTwo, bool fromOurSelf=false);
529  };
530 
548  template<class T, class A, bool mode>
550  {
551 
552  template<typename T1, typename A1>
553  friend class RemoteIndices;
554 
555  public:
556  class InvalidPosition : public RangeError
557  {};
558 
559  enum {
568  MODIFYINDEXSET=mode
569  };
570 
574  typedef T ParallelIndexSet;
575 
580 
585 
589  typedef typename LocalIndex::Attribute Attribute;
590 
595 
599  typedef A Allocator;
600 
604 
609 
614 
628  void insert(const RemoteIndex& index);
629 
630 
645  void insert(const RemoteIndex& index, const GlobalIndex& global);
646 
654  bool remove(const GlobalIndex& global);
655 
669 
670 
672 
678  : glist_()
679  {}
680 
681  private:
682 
689  RemoteIndexList& rList);
690 
691  typedef SLList<GlobalIndex,Allocator> GlobalList;
692  typedef typename GlobalList::ModifyIterator GlobalModifyIterator;
693  RemoteIndexList* rList_;
694  const ParallelIndexSet* indexSet_;
695  GlobalList glist_;
696  ModifyIterator iter_;
697  GlobalModifyIterator giter_;
698  ConstIterator end_;
699  bool first_;
700  GlobalIndex last_;
701  };
702 
707  template<class T, class A>
709  {
710 
714  typedef T ParallelIndexSet;
715 
719  typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
720 
724  typedef typename ParallelIndexSet::LocalIndex LocalIndex;
725 
729  typedef typename LocalIndex::Attribute Attribute;
730 
733 
735  typedef typename A::template rebind<RemoteIndex>::other Allocator;
736 
739 
741  typedef std::map<int,std::pair<typename RemoteIndexList::const_iterator,
742  const typename RemoteIndexList::const_iterator> >
743  Map;
744 
745  public:
746 
748  typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
750 
756  inline CollectiveIterator(const RemoteIndexMap& map_, bool send);
757 
766  inline void advance(const GlobalIndex& global);
767 
777  inline void advance(const GlobalIndex& global, const Attribute& attribute);
778 
779  CollectiveIterator& operator++();
780 
784  inline bool empty();
785 
792  class iterator
793  {
794  public:
795  typedef typename Map::iterator RealIterator;
796  typedef typename Map::iterator ConstRealIterator;
797 
798 
800  iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex& index)
801  : iter_(iter), end_(end), index_(index), hasAttribute(false)
802  {
803  // Move to the first valid entry
804  while(iter_!=end_ && iter_->second.first->localIndexPair().global()!=index_)
805  ++iter_;
806  }
807 
808  iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex index,
809  Attribute attribute)
810  : iter_(iter), end_(end), index_(index), attribute_(attribute), hasAttribute(true)
811  {
812  // Move to the first valid entry or the end
813  while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_
814  || iter_->second.first->localIndexPair().local().attribute()!=attribute))
815  ++iter_;
816  }
818  iterator(const iterator& other)
819  : iter_(other.iter_), end_(other.end_), index_(other.index_)
820  { }
821 
824  {
825  ++iter_;
826  // If entry is not valid move on
827  while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_ ||
828  (hasAttribute &&
829  iter_->second.first->localIndexPair().local().attribute()!=attribute_)))
830  ++iter_;
831  assert(iter_==end_ ||
832  (iter_->second.first->localIndexPair().global()==index_));
833  assert(iter_==end_ || !hasAttribute ||
834  (iter_->second.first->localIndexPair().local().attribute()==attribute_));
835  return *this;
836  }
837 
839  const RemoteIndex& operator*() const
840  {
841  return *(iter_->second.first);
842  }
843 
845  int process() const
846  {
847  return iter_->first;
848  }
849 
851  const RemoteIndex* operator->() const
852  {
853  return iter_->second.first.operator->();
854  }
855 
857  bool operator==(const iterator& other)
858  {
859  return other.iter_==iter_;
860  }
861 
863  bool operator!=(const iterator& other)
864  {
865  return other.iter_!=iter_;
866  }
867 
868  private:
869  iterator();
870 
871  RealIterator iter_;
872  RealIterator end_;
873  GlobalIndex index_;
874  Attribute attribute_;
875  bool hasAttribute;
876  };
877 
878  iterator begin();
879 
880  iterator end();
881 
882  private:
883 
884  Map map_;
885  GlobalIndex index_;
886  Attribute attribute_;
887  bool noattribute;
888  };
889 
890  template<typename TG, typename TA>
891  MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::getType()
892  {
893  if(type==MPI_DATATYPE_NULL) {
894  int length[2] = {1, 1};
895  MPI_Aint base;
896  MPI_Aint disp[2];
897  MPI_Datatype types[2] = {MPITraits<TG>::getType(),
898  MPITraits<ParallelLocalIndex<TA> >::getType()};
899  IndexPair<TG,ParallelLocalIndex<TA> > rep;
900  MPI_Get_address(&rep, &base); // lower bound of the datatype
901  MPI_Get_address(&(rep.global_), &disp[0]);
902  MPI_Get_address(&(rep.local_), &disp[1]);
903  for (MPI_Aint& d : disp)
904  d -= base;
905 
906  MPI_Datatype tmp;
907  MPI_Type_create_struct(2, length, disp, types, &tmp);
908 
909  MPI_Type_create_resized(tmp, 0, sizeof(IndexPair<TG,ParallelLocalIndex<TA> >), &type);
910  MPI_Type_commit(&type);
911 
912  MPI_Type_free(&tmp);
913  }
914  return type;
915  }
916 
917  template<typename TG, typename TA>
918  MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::type=MPI_DATATYPE_NULL;
919 
920  template<typename T1, typename T2>
921  RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute, const PairType* local)
922  : localIndex_(local), attribute_(attribute)
923  {}
924 
925  template<typename T1, typename T2>
926  RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute)
927  : localIndex_(0), attribute_(attribute)
928  {}
929 
930  template<typename T1, typename T2>
932  : localIndex_(0), attribute_()
933  {}
934  template<typename T1, typename T2>
935  inline bool RemoteIndex<T1,T2>::operator==(const RemoteIndex& ri) const
936  {
937  return localIndex_==ri.localIndex_ && attribute_==ri.attribute;
938  }
939 
940  template<typename T1, typename T2>
941  inline bool RemoteIndex<T1,T2>::operator!=(const RemoteIndex& ri) const
942  {
943  return localIndex_!=ri.localIndex_ || attribute_!=ri.attribute_;
944  }
945 
946  template<typename T1, typename T2>
947  inline const T2 RemoteIndex<T1,T2>::attribute() const
948  {
949  return T2(attribute_);
950  }
951 
952  template<typename T1, typename T2>
954  {
955  return *localIndex_;
956  }
957 
958  template<typename T, typename A>
960  const ParallelIndexSet& destination,
961  const MPI_Comm& comm,
962  const std::vector<int>& neighbours,
963  bool includeSelf_)
964  : source_(&source), target_(&destination), comm_(comm),
965  sourceSeqNo_(-1), destSeqNo_(-1), publicIgnored(false), firstBuild(true),
966  includeSelf(includeSelf_)
967  {
968  setNeighbours(neighbours);
969  }
970 
971  template<typename T, typename A>
973  {
974  includeSelf=b;
975  }
976 
977  template<typename T, typename A>
979  : source_(0), target_(0), sourceSeqNo_(-1),
980  destSeqNo_(-1), publicIgnored(false), firstBuild(true),
981  includeSelf(false)
982  {}
983 
984  template<class T, typename A>
986  const ParallelIndexSet& destination,
987  const MPI_Comm& comm,
988  const std::vector<int>& neighbours)
989  {
990  free();
991  source_ = &source;
992  target_ = &destination;
993  comm_ = comm;
994  firstBuild = true;
995  setNeighbours(neighbours);
996  }
997 
998  template<typename T, typename A>
1001  {
1002  return *source_;
1003  }
1004 
1005 
1006  template<typename T, typename A>
1007  const typename RemoteIndices<T,A>::ParallelIndexSet&
1009  {
1010  return *target_;
1011  }
1012 
1013 
1014  template<typename T, typename A>
1016  {
1017  free();
1018  }
1019 
1020  template<typename T, typename A>
1021  template<bool ignorePublic>
1023  const ParallelIndexSet& indexSet,
1024  char* p_out, MPI_Datatype type,
1025  int bufferSize,
1026  int *position, int n)
1027  {
1029  // fill with own indices
1030  typedef typename ParallelIndexSet::const_iterator const_iterator;
1031  typedef IndexPair<GlobalIndex,LocalIndex> PairType;
1032  const const_iterator end = indexSet.end();
1033 
1034  //Now pack the source indices
1035  int i=0;
1036  for(const_iterator index = indexSet.begin(); index != end; ++index)
1037  if(ignorePublic || index->local().isPublic()) {
1038 
1039  MPI_Pack(const_cast<PairType*>(&(*index)), 1,
1040  type,
1041  p_out, bufferSize, position, comm_);
1042  pairs[i++] = const_cast<PairType*>(&(*index));
1043 
1044  }
1045  assert(i==n);
1046  }
1047 
1048  template<typename T, typename A>
1049  inline int RemoteIndices<T,A>::noPublic(const ParallelIndexSet& indexSet)
1050  {
1051  typedef typename ParallelIndexSet::const_iterator const_iterator;
1052 
1053  int noPublic=0;
1054 
1055  const const_iterator end=indexSet.end();
1056  for(const_iterator index=indexSet.begin(); index!=end; ++index)
1057  if(index->local().isPublic())
1058  noPublic++;
1059 
1060  return noPublic;
1061 
1062  }
1063 
1064 
1065  template<typename T, typename A>
1066  inline void RemoteIndices<T,A>::unpackCreateRemote(char* p_in, PairType** sourcePairs,
1067  PairType** destPairs, int remoteProc,
1068  int sourcePublish, int destPublish,
1069  int bufferSize, bool sendTwo,
1070  bool fromOurSelf)
1071  {
1072 
1073  // unpack the number of indices we received
1074  int noRemoteSource=-1, noRemoteDest=-1;
1075  char twoIndexSets=0;
1076  int position=0;
1077  // Did we receive two index sets?
1078  MPI_Unpack(p_in, bufferSize, &position, &twoIndexSets, 1, MPI_CHAR, comm_);
1079  // The number of source indices received
1080  MPI_Unpack(p_in, bufferSize, &position, &noRemoteSource, 1, MPI_INT, comm_);
1081  // The number of destination indices received
1082  MPI_Unpack(p_in, bufferSize, &position, &noRemoteDest, 1, MPI_INT, comm_);
1083 
1084 
1085  // Indices for which we receive
1086  RemoteIndexList* receive= new RemoteIndexList();
1087  // Indices for which we send
1088  RemoteIndexList* send=0;
1089 
1090  MPI_Datatype type= MPITraits<PairType>::getType();
1091 
1092  if(!twoIndexSets) {
1093  if(sendTwo) {
1094  send = new RemoteIndexList();
1095  // Create both remote index sets simultaneously
1096  unpackIndices(*send, *receive, noRemoteSource, sourcePairs, sourcePublish,
1097  destPairs, destPublish, p_in, type, &position, bufferSize);
1098  }else{
1099  // we only need one list
1100  unpackIndices(*receive, noRemoteSource, sourcePairs, sourcePublish,
1101  p_in, type, &position, bufferSize, fromOurSelf);
1102  send=receive;
1103  }
1104  }else{
1105 
1106  int oldPos=position;
1107  // Two index sets received
1108  unpackIndices(*receive, noRemoteSource, destPairs, destPublish,
1109  p_in, type, &position, bufferSize, fromOurSelf);
1110  if(!sendTwo)
1111  //unpack source entries again as destination entries
1112  position=oldPos;
1113 
1114  send = new RemoteIndexList();
1115  unpackIndices(*send, noRemoteDest, sourcePairs, sourcePublish,
1116  p_in, type, &position, bufferSize, fromOurSelf);
1117  }
1118 
1119  if(receive->empty() && send->empty()) {
1120  if(send==receive) {
1121  delete send;
1122  }else{
1123  delete send;
1124  delete receive;
1125  }
1126  }else{
1127  remoteIndices_.insert(std::make_pair(remoteProc,
1128  std::make_pair(send,receive)));
1129  }
1130  }
1131 
1132 
1133  template<typename T, typename A>
1134  template<bool ignorePublic>
1135  inline void RemoteIndices<T,A>::buildRemote(bool includeSelf_)
1136  {
1137  // Processor configuration
1138  int rank, procs;
1139  MPI_Comm_rank(comm_, &rank);
1140  MPI_Comm_size(comm_, &procs);
1141 
1142  // number of local indices to publish
1143  // The indices of the destination will be send.
1144  int sourcePublish, destPublish;
1145 
1146  // Do we need to send two index sets?
1147  char sendTwo = (source_ != target_);
1148 
1149  if(procs==1 && !(sendTwo || includeSelf_))
1150  // Nothing to communicate
1151  return;
1152 
1153  sourcePublish = (ignorePublic) ? source_->size() : noPublic(*source_);
1154 
1155  if(sendTwo)
1156  destPublish = (ignorePublic) ? target_->size() : noPublic(*target_);
1157  else
1158  // we only need to send one set of indices
1159  destPublish = 0;
1160 
1161  int maxPublish, publish=sourcePublish+destPublish;
1162 
1163  // Calucate maximum number of indices send
1164  MPI_Allreduce(&publish, &maxPublish, 1, MPI_INT, MPI_MAX, comm_);
1165 
1166  // allocate buffers
1167  typedef IndexPair<GlobalIndex,LocalIndex> PairType;
1168 
1169  PairType** destPairs;
1170  PairType** sourcePairs = new PairType*[sourcePublish>0 ? sourcePublish : 1];
1171 
1172  if(sendTwo)
1173  destPairs = new PairType*[destPublish>0 ? destPublish : 1];
1174  else
1175  destPairs=sourcePairs;
1176 
1177  char** buffer = new char*[2];
1178  int bufferSize;
1179  int position=0;
1180  int intSize;
1181  int charSize;
1182 
1183  // calculate buffer size
1184  MPI_Datatype type = MPITraits<PairType>::getType();
1185 
1186  MPI_Pack_size(maxPublish, type, comm_,
1187  &bufferSize);
1188  MPI_Pack_size(1, MPI_INT, comm_,
1189  &intSize);
1190  MPI_Pack_size(1, MPI_CHAR, comm_,
1191  &charSize);
1192  // Our message will contain the following:
1193  // a bool whether two index sets where sent
1194  // the size of the source and the dest indexset,
1195  // then the source and destination indices
1196  bufferSize += 2 * intSize + charSize;
1197 
1198  if(bufferSize<=0) bufferSize=1;
1199 
1200  buffer[0] = new char[bufferSize];
1201  buffer[1] = new char[bufferSize];
1202 
1203 
1204  // pack entries into buffer[0], p_out below!
1205  MPI_Pack(&sendTwo, 1, MPI_CHAR, buffer[0], bufferSize, &position,
1206  comm_);
1207 
1208  // The number of indices we send for each index set
1209  MPI_Pack(&sourcePublish, 1, MPI_INT, buffer[0], bufferSize, &position,
1210  comm_);
1211  MPI_Pack(&destPublish, 1, MPI_INT, buffer[0], bufferSize, &position,
1212  comm_);
1213 
1214  // Now pack the source indices and setup the destination pairs
1215  packEntries<ignorePublic>(sourcePairs, *source_, buffer[0], type,
1216  bufferSize, &position, sourcePublish);
1217  // If necessary send the dest indices and setup the source pairs
1218  if(sendTwo)
1219  packEntries<ignorePublic>(destPairs, *target_, buffer[0], type,
1220  bufferSize, &position, destPublish);
1221 
1222 
1223  // Update remote indices for ourself
1224  if(sendTwo|| includeSelf_)
1225  unpackCreateRemote(buffer[0], sourcePairs, destPairs, rank, sourcePublish,
1226  destPublish, bufferSize, sendTwo, includeSelf_);
1227 
1228  neighbourIds.erase(rank);
1229 
1230  if(neighbourIds.size()==0)
1231  {
1232  Dune::dvverb<<rank<<": Sending messages in a ring"<<std::endl;
1233  // send messages in ring
1234  for(int proc=1; proc<procs; proc++) {
1235  // pointers to the current input and output buffers
1236  char* p_out = buffer[1-(proc%2)];
1237  char* p_in = buffer[proc%2];
1238 
1239  MPI_Status status;
1240  if(rank%2==0) {
1241  MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
1242  commTag_, comm_);
1243  MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
1244  commTag_, comm_, &status);
1245  }else{
1246  MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
1247  commTag_, comm_, &status);
1248  MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
1249  commTag_, comm_);
1250  }
1251 
1252 
1253  // The process these indices are from
1254  int remoteProc = (rank+procs-proc)%procs;
1255 
1256  unpackCreateRemote(p_in, sourcePairs, destPairs, remoteProc, sourcePublish,
1257  destPublish, bufferSize, sendTwo);
1258 
1259  }
1260 
1261  }
1262  else
1263  {
1264  MPI_Request* requests=new MPI_Request[neighbourIds.size()];
1265  MPI_Request* req=requests;
1266 
1267  typedef typename std::set<int>::size_type size_type;
1268  size_type noNeighbours=neighbourIds.size();
1269 
1270  // setup sends
1271  for(std::set<int>::iterator neighbour=neighbourIds.begin();
1272  neighbour!= neighbourIds.end(); ++neighbour) {
1273  // Only send the information to the neighbouring processors
1274  MPI_Issend(buffer[0], position , MPI_PACKED, *neighbour, commTag_, comm_, req++);
1275  }
1276 
1277  //Test for received messages
1278 
1279  for(size_type received=0; received <noNeighbours; ++received)
1280  {
1281  MPI_Status status;
1282  // probe for next message
1283  MPI_Probe(MPI_ANY_SOURCE, commTag_, comm_, &status);
1284  int remoteProc=status.MPI_SOURCE;
1285  int size;
1286  MPI_Get_count(&status, MPI_PACKED, &size);
1287  // receive message
1288  MPI_Recv(buffer[1], size, MPI_PACKED, remoteProc,
1289  commTag_, comm_, &status);
1290 
1291  unpackCreateRemote(buffer[1], sourcePairs, destPairs, remoteProc, sourcePublish,
1292  destPublish, bufferSize, sendTwo);
1293  }
1294  // wait for completion of pending requests
1295  MPI_Status* statuses = new MPI_Status[neighbourIds.size()];
1296 
1297  if(MPI_ERR_IN_STATUS==MPI_Waitall(neighbourIds.size(), requests, statuses)) {
1298  for(size_type i=0; i < neighbourIds.size(); ++i)
1299  if(statuses[i].MPI_ERROR!=MPI_SUCCESS) {
1300  std::cerr<<rank<<": MPI_Error occurred while receiving message."<<std::endl;
1301  MPI_Abort(comm_, 999);
1302  }
1303  }
1304  delete[] requests;
1305  delete[] statuses;
1306  }
1307 
1308 
1309  // delete allocated memory
1310  if(destPairs!=sourcePairs)
1311  delete[] destPairs;
1312 
1313  delete[] sourcePairs;
1314  delete[] buffer[0];
1315  delete[] buffer[1];
1316  delete[] buffer;
1317  }
1318 
1319  template<typename T, typename A>
1320  inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& remote,
1321  int remoteEntries,
1322  PairType** local,
1323  int localEntries,
1324  char* p_in,
1325  MPI_Datatype type,
1326  int* position,
1327  int bufferSize,
1328  bool fromOurSelf)
1329  {
1330  if(remoteEntries==0)
1331  return;
1332 
1333  PairType index(1);
1334  MPI_Unpack(p_in, bufferSize, position, &index, 1,
1335  type, comm_);
1336  GlobalIndex oldGlobal=index.global();
1337  int n_in=0, localIndex=0;
1338 
1339  //Check if we know the global index
1340  while(localIndex<localEntries) {
1341  if(local[localIndex]->global()==index.global()) {
1342  int oldLocalIndex=localIndex;
1343 
1344  while(localIndex<localEntries &&
1345  local[localIndex]->global()==index.global()) {
1346  if(!fromOurSelf || index.local().attribute() !=
1347  local[localIndex]->local().attribute())
1348  // if index is from us it has to have a different attribute
1349  remote.push_back(RemoteIndex(index.local().attribute(),
1350  local[localIndex]));
1351  localIndex++;
1352  }
1353 
1354  // unpack next remote index
1355  if((++n_in) < remoteEntries) {
1356  MPI_Unpack(p_in, bufferSize, position, &index, 1,
1357  type, comm_);
1358  if(index.global()==oldGlobal)
1359  // Restart comparison for the same global indices
1360  localIndex=oldLocalIndex;
1361  else
1362  oldGlobal=index.global();
1363  }else{
1364  // No more received indices
1365  break;
1366  }
1367  continue;
1368  }
1369 
1370  if (local[localIndex]->global()<index.global()) {
1371  // compare with next entry in our list
1372  ++localIndex;
1373  }else{
1374  // We do not know the index, unpack next
1375  if((++n_in) < remoteEntries) {
1376  MPI_Unpack(p_in, bufferSize, position, &index, 1,
1377  type, comm_);
1378  oldGlobal=index.global();
1379  }else
1380  // No more received indices
1381  break;
1382  }
1383  }
1384 
1385  // Unpack the other received indices without doing anything
1386  while(++n_in < remoteEntries)
1387  MPI_Unpack(p_in, bufferSize, position, &index, 1,
1388  type, comm_);
1389  }
1390 
1391 
1392  template<typename T, typename A>
1393  inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& send,
1394  RemoteIndexList& receive,
1395  int remoteEntries,
1396  PairType** localSource,
1397  int localSourceEntries,
1398  PairType** localDest,
1399  int localDestEntries,
1400  char* p_in,
1401  MPI_Datatype type,
1402  int* position,
1403  int bufferSize)
1404  {
1405  int n_in=0, sourceIndex=0, destIndex=0;
1406 
1407  //Check if we know the global index
1408  while(n_in<remoteEntries && (sourceIndex<localSourceEntries || destIndex<localDestEntries)) {
1409  // Unpack next index
1410  PairType index;
1411  MPI_Unpack(p_in, bufferSize, position, &index, 1,
1412  type, comm_);
1413  n_in++;
1414 
1415  // Advance until global index in localSource and localDest are >= than the one in the unpacked index
1416  while(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()<index.global())
1417  sourceIndex++;
1418 
1419  while(destIndex<localDestEntries && localDest[destIndex]->global()<index.global())
1420  destIndex++;
1421 
1422  // Add a remote index if we found the global index.
1423  if(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()==index.global())
1424  send.push_back(RemoteIndex(index.local().attribute(),
1425  localSource[sourceIndex]));
1426 
1427  if(destIndex < localDestEntries && localDest[destIndex]->global() == index.global())
1428  receive.push_back(RemoteIndex(index.local().attribute(),
1429  localDest[sourceIndex]));
1430  }
1431 
1432  }
1433 
1434  template<typename T, typename A>
1436  {
1437  typedef typename RemoteIndexMap::iterator Iterator;
1438  Iterator lend = remoteIndices_.end();
1439  for(Iterator lists=remoteIndices_.begin(); lists != lend; ++lists) {
1440  if(lists->second.first==lists->second.second) {
1441  // there is only one remote index list.
1442  delete lists->second.first;
1443  }else{
1444  delete lists->second.first;
1445  delete lists->second.second;
1446  }
1447  }
1448  remoteIndices_.clear();
1449  firstBuild=true;
1450  }
1451 
1452  template<typename T, typename A>
1454  {
1455  return remoteIndices_.size();
1456  }
1457 
1458  template<typename T, typename A>
1459  template<bool ignorePublic>
1461  {
1462  // Test whether a rebuild is Needed.
1463  if(firstBuild ||
1464  ignorePublic!=publicIgnored || !
1465  isSynced()) {
1466  free();
1467 
1468  buildRemote<ignorePublic>(includeSelf);
1469 
1470  sourceSeqNo_ = source_->seqNo();
1471  destSeqNo_ = target_->seqNo();
1472  firstBuild=false;
1473  publicIgnored=ignorePublic;
1474  }
1475 
1476 
1477  }
1478 
1479  template<typename T, typename A>
1480  inline bool RemoteIndices<T,A>::isSynced() const
1481  {
1482  return sourceSeqNo_==source_->seqNo() && destSeqNo_ ==target_->seqNo();
1483  }
1484 
1485  template<typename T, typename A>
1486  template<bool mode, bool send>
1488  {
1489 
1490  // The user are on their own now!
1491  // We assume they know what they are doing and just set the
1492  // remote indices to synced status.
1493  sourceSeqNo_ = source_->seqNo();
1494  destSeqNo_ = target_->seqNo();
1495 
1496  typename RemoteIndexMap::iterator found = remoteIndices_.find(process);
1497 
1498  if(found == remoteIndices_.end())
1499  {
1500  if(source_ != target_)
1501  found = remoteIndices_.insert(found, std::make_pair(process,
1502  std::make_pair(new RemoteIndexList(),
1503  new RemoteIndexList())));
1504  else{
1505  RemoteIndexList* rlist = new RemoteIndexList();
1506  found = remoteIndices_.insert(found,
1507  std::make_pair(process,
1508  std::make_pair(rlist, rlist)));
1509  }
1510  }
1511 
1512  firstBuild = false;
1513 
1514  if(send)
1515  return RemoteIndexListModifier<T,A,mode>(*source_, *(found->second.first));
1516  else
1517  return RemoteIndexListModifier<T,A,mode>(*target_, *(found->second.second));
1518  }
1519 
1520  template<typename T, typename A>
1521  inline typename RemoteIndices<T,A>::const_iterator
1523  {
1524  return remoteIndices_.find(proc);
1525  }
1526 
1527  template<typename T, typename A>
1528  inline typename RemoteIndices<T,A>::const_iterator
1530  {
1531  return remoteIndices_.begin();
1532  }
1533 
1534  template<typename T, typename A>
1535  inline typename RemoteIndices<T,A>::const_iterator
1537  {
1538  return remoteIndices_.end();
1539  }
1540 
1541 
1542  template<typename T, typename A>
1544  {
1545  if(neighbours()!=ri.neighbours())
1546  return false;
1547 
1548  typedef RemoteIndexList RList;
1549  typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
1550 
1551  const const_iterator rend = remoteIndices_.end();
1552 
1553  for(const_iterator rindex = remoteIndices_.begin(), rindex1=ri.remoteIndices_.begin(); rindex!=rend; ++rindex, ++rindex1) {
1554  if(rindex->first != rindex1->first)
1555  return false;
1556  if(*(rindex->second.first) != *(rindex1->second.first))
1557  return false;
1558  if(*(rindex->second.second) != *(rindex1->second.second))
1559  return false;
1560  }
1561  return true;
1562  }
1563 
1564  template<class T, class A, bool mode>
1565  RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(const ParallelIndexSet& indexSet,
1566  RemoteIndexList& rList)
1567  : rList_(&rList), indexSet_(&indexSet), iter_(rList.beginModify()), end_(rList.end()), first_(true)
1568  {
1569  if(MODIFYINDEXSET) {
1570  assert(indexSet_);
1571  for(ConstIterator iter=iter_; iter != end_; ++iter)
1572  glist_.push_back(iter->localIndexPair().global());
1573  giter_ = glist_.beginModify();
1574  }
1575  }
1576 
1577  template<typename T, typename A, bool mode>
1578  RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(const RemoteIndexListModifier<T,A,mode>& other)
1579  : rList_(other.rList_), indexSet_(other.indexSet_),
1580  glist_(other.glist_), iter_(other.iter_), giter_(other.giter_), end_(other.end_),
1581  first_(other.first_), last_(other.last_)
1582  {}
1583 
1584  template<typename T, typename A, bool mode>
1586  {
1587  if(MODIFYINDEXSET) {
1588  // repair pointers to local index set.
1589 #ifdef DUNE_ISTL_WITH_CHECKING
1590  if(indexSet_->state()!=GROUND)
1591  DUNE_THROW(InvalidIndexSetState, "Index has to be in ground mode for repairing pointers to indices");
1592 #endif
1593  typedef typename ParallelIndexSet::const_iterator IndexIterator;
1594  typedef typename GlobalList::const_iterator GlobalIterator;
1595  typedef typename RemoteIndexList::iterator Iterator;
1596  GlobalIterator giter = glist_.begin();
1597  IndexIterator index = indexSet_->begin();
1598 
1599  for(Iterator iter=rList_->begin(); iter != end_; ++iter) {
1600  while(index->global()<*giter) {
1601  ++index;
1602 #ifdef DUNE_ISTL_WITH_CHECKING
1603  if(index == indexSet_->end())
1604  DUNE_THROW(InvalidPosition, "No such global index in set!");
1605 #endif
1606  }
1607 
1608 #ifdef DUNE_ISTL_WITH_CHECKING
1609  if(index->global() != *giter)
1610  DUNE_THROW(InvalidPosition, "No such global index in set!");
1611 #endif
1612  iter->localIndex_ = &(*index);
1613  }
1614  }
1615  }
1616 
1617  template<typename T, typename A, bool mode>
1619  {
1620  static_assert(!mode,"Not allowed if the mode indicates that new indices"
1621  "might be added to the underlying index set. Use "
1622  "insert(const RemoteIndex&, const GlobalIndex&) instead");
1623 
1624 #ifdef DUNE_ISTL_WITH_CHECKING
1625  if(!first_ && index.localIndexPair().global()<last_)
1626  DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
1627 #endif
1628  // Move to the correct position
1629  while(iter_ != end_ && iter_->localIndexPair().global() < index.localIndexPair().global()) {
1630  ++iter_;
1631  }
1632 
1633  // No duplicate entries allowed
1634  assert(iter_==end_ || iter_->localIndexPair().global() != index.localIndexPair().global());
1635  iter_.insert(index);
1636  last_ = index.localIndexPair().global();
1637  first_ = false;
1638  }
1639 
1640  template<typename T, typename A, bool mode>
1641  inline void RemoteIndexListModifier<T,A,mode>::insert(const RemoteIndex& index, const GlobalIndex& global)
1642  {
1643  static_assert(mode,"Not allowed if the mode indicates that no new indices"
1644  "might be added to the underlying index set. Use "
1645  "insert(const RemoteIndex&) instead");
1646 #ifdef DUNE_ISTL_WITH_CHECKING
1647  if(!first_ && global<last_)
1648  DUNE_THROW(InvalidPosition, "Modification of remote indices have to occur with ascending global index.");
1649 #endif
1650  // Move to the correct position
1651  while(iter_ != end_ && *giter_ < global) {
1652  ++giter_;
1653  ++iter_;
1654  }
1655 
1656  // No duplicate entries allowed
1657  assert(iter_->localIndexPair().global() != global);
1658  iter_.insert(index);
1659  giter_.insert(global);
1660 
1661  last_ = global;
1662  first_ = false;
1663  }
1664 
1665  template<typename T, typename A, bool mode>
1667  {
1668 #ifdef DUNE_ISTL_WITH_CHECKING
1669  if(!first_ && global<last_)
1670  DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
1671 #endif
1672 
1673  bool found= false;
1674 
1675  if(MODIFYINDEXSET) {
1676  // Move to the correct position
1677  while(iter_!=end_ && *giter_< global) {
1678  ++giter_;
1679  ++iter_;
1680  }
1681  if(*giter_ == global) {
1682  giter_.remove();
1683  iter_.remove();
1684  found=true;
1685  }
1686  }else{
1687  while(iter_!=end_ && iter_->localIndexPair().global() < global)
1688  ++iter_;
1689 
1690  if(iter_->localIndexPair().global()==global) {
1691  iter_.remove();
1692  found = true;
1693  }
1694  }
1695 
1696  last_ = global;
1697  first_ = false;
1698  return found;
1699  }
1700 
1701  template<typename T, typename A>
1702  template<bool send>
1704  {
1705  return CollectiveIterator<T,A>(remoteIndices_, send);
1706  }
1707 
1708  template<typename T, typename A>
1709  inline MPI_Comm RemoteIndices<T,A>::communicator() const
1710  {
1711  return comm_;
1712 
1713  }
1714 
1715  template<typename T, typename A>
1717  {
1718  typedef typename RemoteIndexMap::const_iterator const_iterator;
1719 
1720  const const_iterator end=pmap.end();
1721  for(const_iterator process=pmap.begin(); process != end; ++process) {
1722  const RemoteIndexList* list = send ? process->second.first : process->second.second;
1723  typedef typename RemoteIndexList::const_iterator iterator;
1724  map_.insert(std::make_pair(process->first,
1725  std::pair<iterator, const iterator>(list->begin(), list->end())));
1726  }
1727  }
1728 
1729  template<typename T, typename A>
1730  inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index)
1731  {
1732  typedef typename Map::iterator iterator;
1733  typedef typename Map::const_iterator const_iterator;
1734  const const_iterator end = map_.end();
1735 
1736  for(iterator iter = map_.begin(); iter != end;) {
1737  // Step the iterator until we are >= index
1738  typename RemoteIndexList::const_iterator current = iter->second.first;
1739  typename RemoteIndexList::const_iterator rend = iter->second.second;
1740  RemoteIndex remoteIndex;
1741  if(current != rend)
1742  remoteIndex = *current;
1743 
1744  while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
1745  ++(iter->second.first);
1746 
1747  // erase from the map if there are no more entries.
1748  if(iter->second.first == iter->second.second)
1749  map_.erase(iter++);
1750  else{
1751  ++iter;
1752  }
1753  }
1754  index_=index;
1755  noattribute=true;
1756  }
1757 
1758  template<typename T, typename A>
1759  inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index,
1760  const Attribute& attribute)
1761  {
1762  typedef typename Map::iterator iterator;
1763  typedef typename Map::const_iterator const_iterator;
1764  const const_iterator end = map_.end();
1765 
1766  for(iterator iter = map_.begin(); iter != end;) {
1767  // Step the iterator until we are >= index
1768  typename RemoteIndexList::const_iterator current = iter->second.first;
1769  typename RemoteIndexList::const_iterator rend = iter->second.second;
1770  RemoteIndex remoteIndex;
1771  if(current != rend)
1772  remoteIndex = *current;
1773 
1774  // Move to global index or bigger
1775  while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
1776  ++(iter->second.first);
1777 
1778  // move to attribute or bigger
1779  while(iter->second.first!=iter->second.second
1780  && iter->second.first->localIndexPair().global()==index
1781  && iter->second.first->localIndexPair().local().attribute()<attribute)
1782  ++(iter->second.first);
1783 
1784  // erase from the map if there are no more entries.
1785  if(iter->second.first == iter->second.second)
1786  map_.erase(iter++);
1787  else{
1788  ++iter;
1789  }
1790  }
1791  index_=index;
1792  attribute_=attribute;
1793  noattribute=false;
1794  }
1795 
1796  template<typename T, typename A>
1798  {
1799  typedef typename Map::iterator iterator;
1800  typedef typename Map::const_iterator const_iterator;
1801  const const_iterator end = map_.end();
1802 
1803  for(iterator iter = map_.begin(); iter != end;) {
1804  // Step the iterator until we are >= index
1805  typename RemoteIndexList::const_iterator current = iter->second.first;
1806  typename RemoteIndexList::const_iterator rend = iter->second.second;
1807 
1808  // move all iterators pointing to the current global index to next value
1809  if(iter->second.first->localIndexPair().global()==index_ &&
1810  (noattribute || iter->second.first->localIndexPair().local().attribute() == attribute_))
1811  ++(iter->second.first);
1812 
1813  // erase from the map if there are no more entries.
1814  if(iter->second.first == iter->second.second)
1815  map_.erase(iter++);
1816  else{
1817  ++iter;
1818  }
1819  }
1820  return *this;
1821  }
1822 
1823  template<typename T, typename A>
1825  {
1826  return map_.empty();
1827  }
1828 
1829  template<typename T, typename A>
1830  inline typename CollectiveIterator<T,A>::iterator
1832  {
1833  if(noattribute)
1834  return iterator(map_.begin(), map_.end(), index_);
1835  else
1836  return iterator(map_.begin(), map_.end(), index_,
1837  attribute_);
1838  }
1839 
1840  template<typename T, typename A>
1841  inline typename CollectiveIterator<T,A>::iterator
1842  CollectiveIterator<T,A>::end()
1843  {
1844  return iterator(map_.end(), map_.end(), index_);
1845  }
1846 
1847  template<typename TG, typename TA>
1848  inline std::ostream& operator<<(std::ostream& os, const RemoteIndex<TG,TA>& index)
1849  {
1850  os<<"[global="<<index.localIndexPair().global()<<", remote attribute="<<index.attribute()<<" local attribute="<<index.localIndexPair().local().attribute()<<"]";
1851  return os;
1852  }
1853 
1854  template<typename T, typename A>
1855  inline std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices)
1856  {
1857  int rank;
1858  MPI_Comm_rank(indices.comm_, &rank);
1859 
1860  typedef typename RemoteIndices<T,A>::RemoteIndexList RList;
1861  typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
1862 
1863  const const_iterator rend = indices.remoteIndices_.end();
1864 
1865  for(const_iterator rindex = indices.remoteIndices_.begin(); rindex!=rend; ++rindex) {
1866  os<<rank<<": Prozess "<<rindex->first<<":";
1867 
1868  if(!rindex->second.first->empty()) {
1869  os<<" send:";
1870 
1871  const typename RList::const_iterator send= rindex->second.first->end();
1872 
1873  for(typename RList::const_iterator index = rindex->second.first->begin();
1874  index != send; ++index)
1875  os<<*index<<" ";
1876  os<<std::endl;
1877  }
1878  if(!rindex->second.second->empty()) {
1879  os<<rank<<": Prozess "<<rindex->first<<": "<<"receive: ";
1880 
1881  for(const auto& index : *(rindex->second.second))
1882  os << index << " ";
1883  }
1884  os<<std::endl<<std::flush;
1885  }
1886  return os;
1887  }
1889 }
1890 
1891 #endif // HAVE_MPI
1892 
1893 #endif
Iterator over the valid underlying iterators.
Definition: remoteindices.hh:793
const RemoteIndex * operator->() const
Definition: remoteindices.hh:851
iterator & operator++()
Definition: remoteindices.hh:823
bool operator==(const iterator &other)
Definition: remoteindices.hh:857
iterator(const RealIterator &iter, const ConstRealIterator &end, GlobalIndex &index)
Definition: remoteindices.hh:800
bool operator!=(const iterator &other)
Definition: remoteindices.hh:863
iterator(const iterator &other)
Definition: remoteindices.hh:818
const RemoteIndex & operator*() const
Definition: remoteindices.hh:839
int process() const
Definition: remoteindices.hh:845
A collective iterator for moving over the remote indices for all processes collectively.
Definition: remoteindices.hh:709
CollectiveIterator(const RemoteIndexMap &map_, bool send)
Constructor.
Definition: remoteindices.hh:1716
bool empty()
Checks whether there are still iterators in the map.
Definition: remoteindices.hh:1824
void advance(const GlobalIndex &global)
Advances all underlying iterators.
Definition: remoteindices.hh:1730
std::map< int, std::pair< RemoteIndexList *, RemoteIndexList * > > RemoteIndexMap
The type of the map from rank to remote index list.
Definition: remoteindices.hh:749
A constant random access iterator for the Dune::ArrayList class.
Definition: arraylist.hh:379
A pair consisting of a global and local index.
Definition: indexset.hh:84
Class for recomputing missing indices of a distributed index set.
Definition: indicessyncer.hh:40
Base class of all classes representing a communication interface.
Definition: interface.hh:33
Exception indicating that the index set is not in the expected state.
Definition: indexset.hh:204
A class setting up standard communication for a two-valued attribute set with owner/overlap/copy sema...
Definition: owneroverlapcopy.hh:172
Manager class for the mapping between local indices and globally unique indices.
Definition: indexset.hh:217
An index present on the local process with an additional attribute flag.
Definition: plocalindex.hh:47
Default exception class for range errors.
Definition: exceptions.hh:252
Modifier for adding and/or deleting remote indices from the remote index list.
Definition: remoteindices.hh:550
void repairLocalIndexPointers()
Repair the pointers to the local index pairs.
Definition: remoteindices.hh:1585
Dune::SLList< RemoteIndex, Allocator > RemoteIndexList
The type of the remote index list.
Definition: remoteindices.hh:603
@ MODIFYINDEXSET
If true the index set corresponding to the remote indices might get modified.
Definition: remoteindices.hh:568
A Allocator
The type of the allocator for the remote index list.
Definition: remoteindices.hh:599
void insert(const RemoteIndex &index)
Insert an index to the list.
Definition: remoteindices.hh:1618
ParallelIndexSet::GlobalIndex GlobalIndex
The type of the global index.
Definition: remoteindices.hh:579
ParallelIndexSet::LocalIndex LocalIndex
The type of the local index.
Definition: remoteindices.hh:584
RemoteIndexList::const_iterator ConstIterator
The type of the remote index list iterator.
Definition: remoteindices.hh:613
SLListModifyIterator< RemoteIndex, Allocator > ModifyIterator
The type of the modifying iterator of the remote index list.
Definition: remoteindices.hh:608
bool remove(const GlobalIndex &global)
Remove a remote index.
Definition: remoteindices.hh:1666
T ParallelIndexSet
Type of the index set we use.
Definition: remoteindices.hh:574
RemoteIndexListModifier()
Default constructor.
Definition: remoteindices.hh:677
LocalIndex::Attribute Attribute
The type of the attribute.
Definition: remoteindices.hh:589
Dune::RemoteIndex< GlobalIndex, Attribute > RemoteIndex
Type of the remote indices we manage.
Definition: remoteindices.hh:594
Information about an index residing on another processor.
Definition: remoteindices.hh:71
const Attribute attribute() const
Get the attribute of the index on the remote process.
Definition: remoteindices.hh:947
T1 GlobalIndex
the type of the global index. This type has to provide at least a operator< for sorting.
Definition: remoteindices.hh:88
T2 Attribute
The type of the attributes. Normally this will be an enumeration like.
Definition: remoteindices.hh:97
IndexPair< GlobalIndex, ParallelLocalIndex< Attribute > > PairType
The type of the index pair.
Definition: remoteindices.hh:103
const PairType & localIndexPair() const
Get the corresponding local index pair.
Definition: remoteindices.hh:953
RemoteIndex()
Parameterless Constructor.
Definition: remoteindices.hh:931
The indices present on remote processes.
Definition: remoteindices.hh:187
Dune::RemoteIndex< GlobalIndex, Attribute > RemoteIndex
Type of the remote indices we manage.
Definition: remoteindices.hh:229
friend void fillIndexSetHoles(const G &graph, Dune::OwnerOverlapCopyCommunication< T1, T2 > &oocomm)
Fills the holes in an index set.
Definition: repartition.hh:81
void setIndexSets(const ParallelIndexSet &source, const ParallelIndexSet &destination, const MPI_Comm &comm, const std::vector< int > &neighbours=std::vector< int >())
Set the index sets and communicator we work with.
Definition: remoteindices.hh:985
void free()
Free the index lists.
Definition: remoteindices.hh:1435
ParallelIndexSet::GlobalIndex GlobalIndex
The type of the global index.
Definition: remoteindices.hh:213
void rebuild()
Rebuilds the set of remote indices.
Definition: remoteindices.hh:1460
T ParallelIndexSet
Type of the index set we use, e.g. ParallelLocalIndexSet.
Definition: remoteindices.hh:204
MPI_Comm communicator() const
Get the mpi communicator used.
Definition: remoteindices.hh:1709
LocalIndex::Attribute Attribute
The type of the attribute.
Definition: remoteindices.hh:224
CollectiveIteratorT iterator() const
Get an iterator for colletively iterating over the remote indices of all remote processes.
Definition: remoteindices.hh:1703
void setIncludeSelf(bool includeSelf)
Tell whether sending from indices of the processor to other indices on the same processor is enabled ...
Definition: remoteindices.hh:972
const_iterator end() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1536
std::map< int, std::pair< RemoteIndexList *, RemoteIndexList * > > RemoteIndexMap
The type of the map from rank to remote index list.
Definition: remoteindices.hh:243
RemoteIndexListModifier< T, A, mode > getModifier(int process)
Get a modifier for a remote index list.
Definition: remoteindices.hh:1487
const ParallelIndexSet & sourceIndexSet() const
Get the index set at the source.
Definition: remoteindices.hh:1000
~RemoteIndices()
Destructor.
Definition: remoteindices.hh:1015
Dune::SLList< RemoteIndex, Allocator > RemoteIndexList
The type of the remote index list.
Definition: remoteindices.hh:239
int neighbours() const
Get the number of processors we share indices with.
Definition: remoteindices.hh:1453
CollectiveIterator< T, A > CollectiveIteratorT
The type of the collective iterator over all remote indices.
Definition: remoteindices.hh:208
const ParallelIndexSet & destinationIndexSet() const
Get the index set at destination.
Definition: remoteindices.hh:1008
RemoteIndices(const ParallelIndexSet &source, const ParallelIndexSet &destination, const MPI_Comm &comm, const std::vector< int > &neighbours=std::vector< int >(), bool includeSelf=false)
Constructor.
Definition: remoteindices.hh:959
const_iterator find(int proc) const
Find an iterator over the remote index lists of a specific process.
Definition: remoteindices.hh:1522
bool isSynced() const
Checks whether the remote indices are synced with the indexsets.
Definition: remoteindices.hh:1480
const_iterator begin() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1529
ParallelIndexSet::LocalIndex LocalIndex
The type of the local index.
Definition: remoteindices.hh:219
A::template rebind< RemoteIndex >::other Allocator
The type of the allocator for the remote index list.
Definition: remoteindices.hh:235
A constant iterator for the SLList.
Definition: sllist.hh:369
A mutable iterator for the SLList.
Definition: sllist.hh:269
A single linked list.
Definition: sllist.hh:42
A few common exception classes.
ArrayList< IndexPair, N >::const_iterator const_iterator
The constant iterator over the pairs.
Definition: indexset.hh:307
void repairLocalIndexPointers(std::map< int, SLList< std::pair< typename T::GlobalIndex, typename T::LocalIndex::Attribute >, A > > &globalMap, RemoteIndices< T, A1 > &remoteIndices, const T &indexSet)
Repair the pointers to the local indices in the remote indices.
Definition: indicessyncer.hh:490
iterator begin()
Get an iterator over the indices positioned at the first index.
iterator end()
Get an iterator over the indices positioned after the last index.
const GlobalIndex & global() const
Get the global index.
TL LocalIndex
The type of the local index, e.g. ParallelLocalIndex.
Definition: indexset.hh:238
TG GlobalIndex
the type of the global index. This type has to provide at least a operator< for sorting.
Definition: indexset.hh:225
@ GROUND
The default mode. Indicates that the index set is ready to be used.
Definition: indexset.hh:185
iterator end()
Get an iterator pointing to the end of the list.
Definition: sllist.hh:788
SLListConstIterator< T, A > const_iterator
The constant iterator of the list.
Definition: sllist.hh:72
SLListModifyIterator< GlobalIndex, Allocator > ModifyIterator
The type of the iterator capable of deletion and insertion.
Definition: sllist.hh:101
iterator begin()
Get an iterator pointing to the first element in the list.
Definition: sllist.hh:776
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
EnableIfInterOperable< T1, T2, bool >::type operator!=(const ForwardIteratorFacade< T1, V1, R1, D > &lhs, const ForwardIteratorFacade< T2, V2, R2, D > &rhs)
Checks for inequality.
Definition: iteratorfacades.hh:255
EnableIfInterOperable< T1, T2, bool >::type operator==(const ForwardIteratorFacade< T1, V1, R1, D > &lhs, const ForwardIteratorFacade< T2, V2, R2, D > &rhs)
Checks for equality.
Definition: iteratorfacades.hh:233
DVVerbType dvverb(std::cout)
stream for very verbose output.
Definition: stdstreams.hh:93
Provides a map between global and local indices.
Traits classes for mapping types onto MPI_Datatype.
Dune namespace.
Definition: alignedallocator.hh:10
Provides classes for use as the local index in ParallelIndexSet for distributed computing.
Implements a singly linked list together with the necessary iterators.
Standard Dune debug streams.
A traits class describing the mapping of types onto MPI_Datatypes.
Definition: mpitraits.hh:38
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 24, 22:30, 2024)