DUNE PDELab (2.7)

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
26namespace 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
205
209
214
215
220
224 typedef typename LocalIndex::Attribute Attribute;
225
230
231
235 using Allocator = typename std::allocator_traits<A>::template rebind_alloc<RemoteIndex>;
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
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
402
403 private:
405 RemoteIndices(const RemoteIndices&) = delete;
406
408 const ParallelIndexSet* source_;
409
411 const ParallelIndexSet* target_;
412
414 MPI_Comm comm_;
415
418 std::set<int> neighbourIds;
419
421 const static int commTag_=333;
422
427 int sourceSeqNo_;
428
433 int destSeqNo_;
434
438 bool publicIgnored;
439
443 bool firstBuild;
444
445 /*
446 * @brief If true, sending from indices of the processor to other
447 * indices on the same processor is enabled even if the same indexset is used
448 * on both the
449 * sending and receiving side.
450 */
451 bool includeSelf;
452
455 PairType;
456
463 RemoteIndexMap remoteIndices_;
464
475 template<bool ignorePublic>
476 inline void buildRemote(bool includeSelf);
477
483 inline int noPublic(const ParallelIndexSet& indexSet);
484
496 template<bool ignorePublic>
497 inline void packEntries(PairType** myPairs, const ParallelIndexSet& indexSet,
498 char* p_out, MPI_Datatype type, int bufferSize,
499 int* position, int n);
500
514 inline void unpackIndices(RemoteIndexList& remote, int remoteEntries,
515 PairType** local, int localEntries, char* p_in,
516 MPI_Datatype type, int* position, int bufferSize,
517 bool fromOurself);
518
519 inline void unpackIndices(RemoteIndexList& send, RemoteIndexList& receive,
520 int remoteEntries, PairType** localSource,
521 int localSourceEntries, PairType** localDest,
522 int localDestEntries, char* p_in,
523 MPI_Datatype type, int* position, int bufferSize);
524
525 void unpackCreateRemote(char* p_in, PairType** sourcePairs, PairType** DestPairs,
526 int remoteProc, int sourcePublish, int destPublish,
527 int bufferSize, bool sendTwo, bool fromOurSelf=false);
528 };
529
547 template<class T, class A, bool mode>
549 {
550
551 template<typename T1, typename A1>
552 friend class RemoteIndices;
553
554 public:
555 class InvalidPosition : public RangeError
556 {};
557
558 enum {
567 MODIFYINDEXSET=mode
568 };
569
574
579
584
588 typedef typename LocalIndex::Attribute Attribute;
589
594
598 typedef A Allocator;
599
603
608
613
627 void insert(const RemoteIndex& index);
628
629
644 void insert(const RemoteIndex& index, const GlobalIndex& global);
645
653 bool remove(const GlobalIndex& global);
654
668
669
671
677 : glist_()
678 {}
679
680 private:
681
688 RemoteIndexList& rList);
689
690 typedef SLList<GlobalIndex,Allocator> GlobalList;
691 typedef typename GlobalList::ModifyIterator GlobalModifyIterator;
692 RemoteIndexList* rList_;
693 const ParallelIndexSet* indexSet_;
694 GlobalList glist_;
695 ModifyIterator iter_;
696 GlobalModifyIterator giter_;
697 ConstIterator end_;
698 bool first_;
699 GlobalIndex last_;
700 };
701
706 template<class T, class A>
708 {
709
713 typedef T ParallelIndexSet;
714
718 typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
719
723 typedef typename ParallelIndexSet::LocalIndex LocalIndex;
724
728 typedef typename LocalIndex::Attribute Attribute;
729
732
734 using Allocator = typename std::allocator_traits<A>::template rebind_alloc<RemoteIndex>;
735
738
740 typedef std::map<int,std::pair<typename RemoteIndexList::const_iterator,
741 const typename RemoteIndexList::const_iterator> >
742 Map;
743
744 public:
745
747 typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
749
755 inline CollectiveIterator(const RemoteIndexMap& map_, bool send);
756
765 inline void advance(const GlobalIndex& global);
766
776 inline void advance(const GlobalIndex& global, const Attribute& attribute);
777
778 CollectiveIterator& operator++();
779
783 inline bool empty();
784
792 {
793 public:
794 typedef typename Map::iterator RealIterator;
795 typedef typename Map::iterator ConstRealIterator;
796
797
799 iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex& index)
800 : iter_(iter), end_(end), index_(index), hasAttribute(false)
801 {
802 // Move to the first valid entry
803 while(iter_!=end_ && iter_->second.first->localIndexPair().global()!=index_)
804 ++iter_;
805 }
806
807 iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex index,
808 Attribute attribute)
809 : iter_(iter), end_(end), index_(index), attribute_(attribute), hasAttribute(true)
810 {
811 // Move to the first valid entry or the end
812 while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_
813 || iter_->second.first->localIndexPair().local().attribute()!=attribute))
814 ++iter_;
815 }
817 iterator(const iterator& other)
818 : iter_(other.iter_), end_(other.end_), index_(other.index_)
819 { }
820
823 {
824 ++iter_;
825 // If entry is not valid move on
826 while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_ ||
827 (hasAttribute &&
828 iter_->second.first->localIndexPair().local().attribute()!=attribute_)))
829 ++iter_;
830 assert(iter_==end_ ||
831 (iter_->second.first->localIndexPair().global()==index_));
832 assert(iter_==end_ || !hasAttribute ||
833 (iter_->second.first->localIndexPair().local().attribute()==attribute_));
834 return *this;
835 }
836
838 const RemoteIndex& operator*() const
839 {
840 return *(iter_->second.first);
841 }
842
844 int process() const
845 {
846 return iter_->first;
847 }
848
850 const RemoteIndex* operator->() const
851 {
852 return iter_->second.first.operator->();
853 }
854
856 bool operator==(const iterator& other) const
857 {
858 return other.iter_==iter_;
859 }
860
862 bool operator!=(const iterator& other) const
863 {
864 return other.iter_!=iter_;
865 }
866
867 private:
868 iterator();
869
870 RealIterator iter_;
871 RealIterator end_;
872 GlobalIndex index_;
873 Attribute attribute_;
874 bool hasAttribute;
875 };
876
877 iterator begin();
878
879 iterator end();
880
881 private:
882
883 Map map_;
884 GlobalIndex index_;
885 Attribute attribute_;
886 bool noattribute;
887 };
888
889 template<typename TG, typename TA>
890 MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::getType()
891 {
892 if(type==MPI_DATATYPE_NULL) {
893 int length[2] = {1, 1};
894 MPI_Aint base;
895 MPI_Aint disp[2];
896 MPI_Datatype types[2] = {MPITraits<TG>::getType(),
897 MPITraits<ParallelLocalIndex<TA> >::getType()};
898 IndexPair<TG,ParallelLocalIndex<TA> > rep;
899 MPI_Get_address(&rep, &base); // lower bound of the datatype
900 MPI_Get_address(&(rep.global_), &disp[0]);
901 MPI_Get_address(&(rep.local_), &disp[1]);
902 for (MPI_Aint& d : disp)
903 d -= base;
904
905 MPI_Datatype tmp;
906 MPI_Type_create_struct(2, length, disp, types, &tmp);
907
908 MPI_Type_create_resized(tmp, 0, sizeof(IndexPair<TG,ParallelLocalIndex<TA> >), &type);
909 MPI_Type_commit(&type);
910
911 MPI_Type_free(&tmp);
912 }
913 return type;
914 }
915
916 template<typename TG, typename TA>
917 MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::type=MPI_DATATYPE_NULL;
918
919 template<typename T1, typename T2>
920 RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute, const PairType* local)
921 : localIndex_(local), attribute_(attribute)
922 {}
923
924 template<typename T1, typename T2>
926 : localIndex_(0), attribute_(attribute)
927 {}
928
929 template<typename T1, typename T2>
931 : localIndex_(0), attribute_()
932 {}
933 template<typename T1, typename T2>
934 inline bool RemoteIndex<T1,T2>::operator==(const RemoteIndex& ri) const
935 {
936 return localIndex_==ri.localIndex_ && attribute_==ri.attribute;
937 }
938
939 template<typename T1, typename T2>
940 inline bool RemoteIndex<T1,T2>::operator!=(const RemoteIndex& ri) const
941 {
942 return localIndex_!=ri.localIndex_ || attribute_!=ri.attribute_;
943 }
944
945 template<typename T1, typename T2>
946 inline const T2 RemoteIndex<T1,T2>::attribute() const
947 {
948 return T2(attribute_);
949 }
950
951 template<typename T1, typename T2>
953 {
954 return *localIndex_;
955 }
956
957 template<typename T, typename A>
959 const ParallelIndexSet& destination,
960 const MPI_Comm& comm,
961 const std::vector<int>& neighbours,
962 bool includeSelf_)
963 : source_(&source), target_(&destination), comm_(comm),
964 sourceSeqNo_(-1), destSeqNo_(-1), publicIgnored(false), firstBuild(true),
965 includeSelf(includeSelf_)
966 {
967 setNeighbours(neighbours);
968 }
969
970 template<typename T, typename A>
972 {
973 includeSelf=b;
974 }
975
976 template<typename T, typename A>
978 : source_(0), target_(0), sourceSeqNo_(-1),
979 destSeqNo_(-1), publicIgnored(false), firstBuild(true),
980 includeSelf(false)
981 {}
982
983 template<class T, typename A>
985 const ParallelIndexSet& destination,
986 const MPI_Comm& comm,
987 const std::vector<int>& neighbours)
988 {
989 free();
990 source_ = &source;
991 target_ = &destination;
992 comm_ = comm;
993 firstBuild = true;
994 setNeighbours(neighbours);
995 }
996
997 template<typename T, typename A>
1000 {
1001 return *source_;
1002 }
1003
1004
1005 template<typename T, typename A>
1008 {
1009 return *target_;
1010 }
1011
1012
1013 template<typename T, typename A>
1015 {
1016 free();
1017 }
1018
1019 template<typename T, typename A>
1020 template<bool ignorePublic>
1022 const ParallelIndexSet& indexSet,
1023 char* p_out, MPI_Datatype type,
1024 int bufferSize,
1025 int *position, int n)
1026 {
1028 // fill with own indices
1029 typedef typename ParallelIndexSet::const_iterator const_iterator;
1030 typedef IndexPair<GlobalIndex,LocalIndex> PairType;
1031 const const_iterator end = indexSet.end();
1032
1033 //Now pack the source indices
1034 int i=0;
1035 for(const_iterator index = indexSet.begin(); index != end; ++index)
1036 if(ignorePublic || index->local().isPublic()) {
1037
1038 MPI_Pack(const_cast<PairType*>(&(*index)), 1,
1039 type,
1040 p_out, bufferSize, position, comm_);
1041 pairs[i++] = const_cast<PairType*>(&(*index));
1042
1043 }
1044 assert(i==n);
1045 }
1046
1047 template<typename T, typename A>
1048 inline int RemoteIndices<T,A>::noPublic(const ParallelIndexSet& indexSet)
1049 {
1050 typedef typename ParallelIndexSet::const_iterator const_iterator;
1051
1052 int noPublic=0;
1053
1054 const const_iterator end=indexSet.end();
1055 for(const_iterator index=indexSet.begin(); index!=end; ++index)
1056 if(index->local().isPublic())
1057 noPublic++;
1058
1059 return noPublic;
1060
1061 }
1062
1063
1064 template<typename T, typename A>
1065 inline void RemoteIndices<T,A>::unpackCreateRemote(char* p_in, PairType** sourcePairs,
1066 PairType** destPairs, int remoteProc,
1067 int sourcePublish, int destPublish,
1068 int bufferSize, bool sendTwo,
1069 bool fromOurSelf)
1070 {
1071
1072 // unpack the number of indices we received
1073 int noRemoteSource=-1, noRemoteDest=-1;
1074 char twoIndexSets=0;
1075 int position=0;
1076 // Did we receive two index sets?
1077 MPI_Unpack(p_in, bufferSize, &position, &twoIndexSets, 1, MPI_CHAR, comm_);
1078 // The number of source indices received
1079 MPI_Unpack(p_in, bufferSize, &position, &noRemoteSource, 1, MPI_INT, comm_);
1080 // The number of destination indices received
1081 MPI_Unpack(p_in, bufferSize, &position, &noRemoteDest, 1, MPI_INT, comm_);
1082
1083
1084 // Indices for which we receive
1085 RemoteIndexList* receive= new RemoteIndexList();
1086 // Indices for which we send
1087 RemoteIndexList* send=0;
1088
1089 MPI_Datatype type= MPITraits<PairType>::getType();
1090
1091 if(!twoIndexSets) {
1092 if(sendTwo) {
1093 send = new RemoteIndexList();
1094 // Create both remote index sets simultaneously
1095 unpackIndices(*send, *receive, noRemoteSource, sourcePairs, sourcePublish,
1096 destPairs, destPublish, p_in, type, &position, bufferSize);
1097 }else{
1098 // we only need one list
1099 unpackIndices(*receive, noRemoteSource, sourcePairs, sourcePublish,
1100 p_in, type, &position, bufferSize, fromOurSelf);
1101 send=receive;
1102 }
1103 }else{
1104
1105 int oldPos=position;
1106 // Two index sets received
1107 unpackIndices(*receive, noRemoteSource, destPairs, destPublish,
1108 p_in, type, &position, bufferSize, fromOurSelf);
1109 if(!sendTwo)
1110 //unpack source entries again as destination entries
1111 position=oldPos;
1112
1113 send = new RemoteIndexList();
1114 unpackIndices(*send, noRemoteDest, sourcePairs, sourcePublish,
1115 p_in, type, &position, bufferSize, fromOurSelf);
1116 }
1117
1118 if(receive->empty() && send->empty()) {
1119 if(send==receive) {
1120 delete send;
1121 }else{
1122 delete send;
1123 delete receive;
1124 }
1125 }else{
1126 remoteIndices_.insert(std::make_pair(remoteProc,
1127 std::make_pair(send,receive)));
1128 }
1129 }
1130
1131
1132 template<typename T, typename A>
1133 template<bool ignorePublic>
1134 inline void RemoteIndices<T,A>::buildRemote(bool includeSelf_)
1135 {
1136 // Processor configuration
1137 int rank, procs;
1138 MPI_Comm_rank(comm_, &rank);
1139 MPI_Comm_size(comm_, &procs);
1140
1141 // number of local indices to publish
1142 // The indices of the destination will be send.
1143 int sourcePublish, destPublish;
1144
1145 // Do we need to send two index sets?
1146 char sendTwo = (source_ != target_);
1147
1148 if(procs==1 && !(sendTwo || includeSelf_))
1149 // Nothing to communicate
1150 return;
1151
1152 sourcePublish = (ignorePublic) ? source_->size() : noPublic(*source_);
1153
1154 if(sendTwo)
1155 destPublish = (ignorePublic) ? target_->size() : noPublic(*target_);
1156 else
1157 // we only need to send one set of indices
1158 destPublish = 0;
1159
1160 int maxPublish, publish=sourcePublish+destPublish;
1161
1162 // Calucate maximum number of indices send
1163 MPI_Allreduce(&publish, &maxPublish, 1, MPI_INT, MPI_MAX, comm_);
1164
1165 // allocate buffers
1166 typedef IndexPair<GlobalIndex,LocalIndex> PairType;
1167
1168 PairType** destPairs;
1169 PairType** sourcePairs = new PairType*[sourcePublish>0 ? sourcePublish : 1];
1170
1171 if(sendTwo)
1172 destPairs = new PairType*[destPublish>0 ? destPublish : 1];
1173 else
1174 destPairs=sourcePairs;
1175
1176 char** buffer = new char*[2];
1177 int bufferSize;
1178 int position=0;
1179 int intSize;
1180 int charSize;
1181
1182 // calculate buffer size
1183 MPI_Datatype type = MPITraits<PairType>::getType();
1184
1185 MPI_Pack_size(maxPublish, type, comm_,
1186 &bufferSize);
1187 MPI_Pack_size(1, MPI_INT, comm_,
1188 &intSize);
1189 MPI_Pack_size(1, MPI_CHAR, comm_,
1190 &charSize);
1191 // Our message will contain the following:
1192 // a bool whether two index sets where sent
1193 // the size of the source and the dest indexset,
1194 // then the source and destination indices
1195 bufferSize += 2 * intSize + charSize;
1196
1197 if(bufferSize<=0) bufferSize=1;
1198
1199 buffer[0] = new char[bufferSize];
1200 buffer[1] = new char[bufferSize];
1201
1202
1203 // pack entries into buffer[0], p_out below!
1204 MPI_Pack(&sendTwo, 1, MPI_CHAR, buffer[0], bufferSize, &position,
1205 comm_);
1206
1207 // The number of indices we send for each index set
1208 MPI_Pack(&sourcePublish, 1, MPI_INT, buffer[0], bufferSize, &position,
1209 comm_);
1210 MPI_Pack(&destPublish, 1, MPI_INT, buffer[0], bufferSize, &position,
1211 comm_);
1212
1213 // Now pack the source indices and setup the destination pairs
1214 packEntries<ignorePublic>(sourcePairs, *source_, buffer[0], type,
1215 bufferSize, &position, sourcePublish);
1216 // If necessary send the dest indices and setup the source pairs
1217 if(sendTwo)
1218 packEntries<ignorePublic>(destPairs, *target_, buffer[0], type,
1219 bufferSize, &position, destPublish);
1220
1221
1222 // Update remote indices for ourself
1223 if(sendTwo|| includeSelf_)
1224 unpackCreateRemote(buffer[0], sourcePairs, destPairs, rank, sourcePublish,
1225 destPublish, bufferSize, sendTwo, includeSelf_);
1226
1227 neighbourIds.erase(rank);
1228
1229 if(neighbourIds.size()==0)
1230 {
1231 Dune::dvverb<<rank<<": Sending messages in a ring"<<std::endl;
1232 // send messages in ring
1233 for(int proc=1; proc<procs; proc++) {
1234 // pointers to the current input and output buffers
1235 char* p_out = buffer[1-(proc%2)];
1236 char* p_in = buffer[proc%2];
1237
1238 MPI_Status status;
1239 if(rank%2==0) {
1240 MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
1241 commTag_, comm_);
1242 MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
1243 commTag_, comm_, &status);
1244 }else{
1245 MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
1246 commTag_, comm_, &status);
1247 MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
1248 commTag_, comm_);
1249 }
1250
1251
1252 // The process these indices are from
1253 int remoteProc = (rank+procs-proc)%procs;
1254
1255 unpackCreateRemote(p_in, sourcePairs, destPairs, remoteProc, sourcePublish,
1256 destPublish, bufferSize, sendTwo);
1257
1258 }
1259
1260 }
1261 else
1262 {
1263 MPI_Request* requests=new MPI_Request[neighbourIds.size()];
1264 MPI_Request* req=requests;
1265
1266 typedef typename std::set<int>::size_type size_type;
1267 size_type noNeighbours=neighbourIds.size();
1268
1269 // setup sends
1270 for(std::set<int>::iterator neighbour=neighbourIds.begin();
1271 neighbour!= neighbourIds.end(); ++neighbour) {
1272 // Only send the information to the neighbouring processors
1273 MPI_Issend(buffer[0], position , MPI_PACKED, *neighbour, commTag_, comm_, req++);
1274 }
1275
1276 //Test for received messages
1277
1278 for(size_type received=0; received <noNeighbours; ++received)
1279 {
1280 MPI_Status status;
1281 // probe for next message
1282 MPI_Probe(MPI_ANY_SOURCE, commTag_, comm_, &status);
1283 int remoteProc=status.MPI_SOURCE;
1284 int size;
1285 MPI_Get_count(&status, MPI_PACKED, &size);
1286 // receive message
1287 MPI_Recv(buffer[1], size, MPI_PACKED, remoteProc,
1288 commTag_, comm_, &status);
1289
1290 unpackCreateRemote(buffer[1], sourcePairs, destPairs, remoteProc, sourcePublish,
1291 destPublish, bufferSize, sendTwo);
1292 }
1293 // wait for completion of pending requests
1294 MPI_Status* statuses = new MPI_Status[neighbourIds.size()];
1295
1296 if(MPI_ERR_IN_STATUS==MPI_Waitall(neighbourIds.size(), requests, statuses)) {
1297 for(size_type i=0; i < neighbourIds.size(); ++i)
1298 if(statuses[i].MPI_ERROR!=MPI_SUCCESS) {
1299 std::cerr<<rank<<": MPI_Error occurred while receiving message."<<std::endl;
1300 MPI_Abort(comm_, 999);
1301 }
1302 }
1303 delete[] requests;
1304 delete[] statuses;
1305 }
1306
1307
1308 // delete allocated memory
1309 if(destPairs!=sourcePairs)
1310 delete[] destPairs;
1311
1312 delete[] sourcePairs;
1313 delete[] buffer[0];
1314 delete[] buffer[1];
1315 delete[] buffer;
1316 }
1317
1318 template<typename T, typename A>
1319 inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& remote,
1320 int remoteEntries,
1321 PairType** local,
1322 int localEntries,
1323 char* p_in,
1324 MPI_Datatype type,
1325 int* position,
1326 int bufferSize,
1327 bool fromOurSelf)
1328 {
1329 if(remoteEntries==0)
1330 return;
1331
1332 PairType index;
1333 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1334 type, comm_);
1335 GlobalIndex oldGlobal=index.global();
1336 int n_in=0, localIndex=0;
1337
1338 //Check if we know the global index
1339 while(localIndex<localEntries) {
1340 if(local[localIndex]->global()==index.global()) {
1341 int oldLocalIndex=localIndex;
1342
1343 while(localIndex<localEntries &&
1344 local[localIndex]->global()==index.global()) {
1345 if(!fromOurSelf || index.local().attribute() !=
1346 local[localIndex]->local().attribute())
1347 // if index is from us it has to have a different attribute
1348 remote.push_back(RemoteIndex(index.local().attribute(),
1349 local[localIndex]));
1350 localIndex++;
1351 }
1352
1353 // unpack next remote index
1354 if((++n_in) < remoteEntries) {
1355 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1356 type, comm_);
1357 if(index.global()==oldGlobal)
1358 // Restart comparison for the same global indices
1359 localIndex=oldLocalIndex;
1360 else
1361 oldGlobal=index.global();
1362 }else{
1363 // No more received indices
1364 break;
1365 }
1366 continue;
1367 }
1368
1369 if (local[localIndex]->global()<index.global()) {
1370 // compare with next entry in our list
1371 ++localIndex;
1372 }else{
1373 // We do not know the index, unpack next
1374 if((++n_in) < remoteEntries) {
1375 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1376 type, comm_);
1377 oldGlobal=index.global();
1378 }else
1379 // No more received indices
1380 break;
1381 }
1382 }
1383
1384 // Unpack the other received indices without doing anything
1385 while(++n_in < remoteEntries)
1386 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1387 type, comm_);
1388 }
1389
1390
1391 template<typename T, typename A>
1392 inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& send,
1393 RemoteIndexList& receive,
1394 int remoteEntries,
1395 PairType** localSource,
1396 int localSourceEntries,
1397 PairType** localDest,
1398 int localDestEntries,
1399 char* p_in,
1400 MPI_Datatype type,
1401 int* position,
1402 int bufferSize)
1403 {
1404 int n_in=0, sourceIndex=0, destIndex=0;
1405
1406 //Check if we know the global index
1407 while(n_in<remoteEntries && (sourceIndex<localSourceEntries || destIndex<localDestEntries)) {
1408 // Unpack next index
1409 PairType index;
1410 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1411 type, comm_);
1412 n_in++;
1413
1414 // Advance until global index in localSource and localDest are >= than the one in the unpacked index
1415 while(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()<index.global())
1416 sourceIndex++;
1417
1418 while(destIndex<localDestEntries && localDest[destIndex]->global()<index.global())
1419 destIndex++;
1420
1421 // Add a remote index if we found the global index.
1422 if(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()==index.global())
1423 send.push_back(RemoteIndex(index.local().attribute(),
1424 localSource[sourceIndex]));
1425
1426 if(destIndex < localDestEntries && localDest[destIndex]->global() == index.global())
1427 receive.push_back(RemoteIndex(index.local().attribute(),
1428 localDest[sourceIndex]));
1429 }
1430
1431 }
1432
1433 template<typename T, typename A>
1435 {
1436 typedef typename RemoteIndexMap::iterator Iterator;
1437 Iterator lend = remoteIndices_.end();
1438 for(Iterator lists=remoteIndices_.begin(); lists != lend; ++lists) {
1439 if(lists->second.first==lists->second.second) {
1440 // there is only one remote index list.
1441 delete lists->second.first;
1442 }else{
1443 delete lists->second.first;
1444 delete lists->second.second;
1445 }
1446 }
1447 remoteIndices_.clear();
1448 firstBuild=true;
1449 }
1450
1451 template<typename T, typename A>
1453 {
1454 return remoteIndices_.size();
1455 }
1456
1457 template<typename T, typename A>
1458 template<bool ignorePublic>
1460 {
1461 // Test whether a rebuild is Needed.
1462 if(firstBuild ||
1463 ignorePublic!=publicIgnored || !
1464 isSynced()) {
1465 free();
1466
1467 buildRemote<ignorePublic>(includeSelf);
1468
1469 sourceSeqNo_ = source_->seqNo();
1470 destSeqNo_ = target_->seqNo();
1471 firstBuild=false;
1472 publicIgnored=ignorePublic;
1473 }
1474
1475
1476 }
1477
1478 template<typename T, typename A>
1480 {
1481 return sourceSeqNo_==source_->seqNo() && destSeqNo_ ==target_->seqNo();
1482 }
1483
1484 template<typename T, typename A>
1485 template<bool mode, bool send>
1487 {
1488
1489 // The user are on their own now!
1490 // We assume they know what they are doing and just set the
1491 // remote indices to synced status.
1492 sourceSeqNo_ = source_->seqNo();
1493 destSeqNo_ = target_->seqNo();
1494
1495 typename RemoteIndexMap::iterator found = remoteIndices_.find(process);
1496
1497 if(found == remoteIndices_.end())
1498 {
1499 if(source_ != target_)
1500 found = remoteIndices_.insert(found, std::make_pair(process,
1501 std::make_pair(new RemoteIndexList(),
1502 new RemoteIndexList())));
1503 else{
1504 RemoteIndexList* rlist = new RemoteIndexList();
1505 found = remoteIndices_.insert(found,
1506 std::make_pair(process,
1507 std::make_pair(rlist, rlist)));
1508 }
1509 }
1510
1511 firstBuild = false;
1512
1513 if(send)
1514 return RemoteIndexListModifier<T,A,mode>(*source_, *(found->second.first));
1515 else
1516 return RemoteIndexListModifier<T,A,mode>(*target_, *(found->second.second));
1517 }
1518
1519 template<typename T, typename A>
1520 inline typename RemoteIndices<T,A>::const_iterator
1522 {
1523 return remoteIndices_.find(proc);
1524 }
1525
1526 template<typename T, typename A>
1527 inline typename RemoteIndices<T,A>::const_iterator
1529 {
1530 return remoteIndices_.begin();
1531 }
1532
1533 template<typename T, typename A>
1534 inline typename RemoteIndices<T,A>::const_iterator
1536 {
1537 return remoteIndices_.end();
1538 }
1539
1540
1541 template<typename T, typename A>
1543 {
1544 if(neighbours()!=ri.neighbours())
1545 return false;
1546
1547 typedef RemoteIndexList RList;
1548 typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
1549
1550 const const_iterator rend = remoteIndices_.end();
1551
1552 for(const_iterator rindex = remoteIndices_.begin(), rindex1=ri.remoteIndices_.begin(); rindex!=rend; ++rindex, ++rindex1) {
1553 if(rindex->first != rindex1->first)
1554 return false;
1555 if(*(rindex->second.first) != *(rindex1->second.first))
1556 return false;
1557 if(*(rindex->second.second) != *(rindex1->second.second))
1558 return false;
1559 }
1560 return true;
1561 }
1562
1563 template<class T, class A, bool mode>
1564 RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(const ParallelIndexSet& indexSet,
1565 RemoteIndexList& rList)
1566 : rList_(&rList), indexSet_(&indexSet), iter_(rList.beginModify()), end_(rList.end()), first_(true)
1567 {
1568 if(MODIFYINDEXSET) {
1569 assert(indexSet_);
1570 for(ConstIterator iter=iter_; iter != end_; ++iter)
1571 glist_.push_back(iter->localIndexPair().global());
1572 giter_ = glist_.beginModify();
1573 }
1574 }
1575
1576 template<typename T, typename A, bool mode>
1577 RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(const RemoteIndexListModifier<T,A,mode>& other)
1578 : rList_(other.rList_), indexSet_(other.indexSet_),
1579 glist_(other.glist_), iter_(other.iter_), giter_(other.giter_), end_(other.end_),
1580 first_(other.first_), last_(other.last_)
1581 {}
1582
1583 template<typename T, typename A, bool mode>
1585 {
1586 if(MODIFYINDEXSET) {
1587 // repair pointers to local index set.
1588#ifdef DUNE_ISTL_WITH_CHECKING
1589 if(indexSet_->state()!=GROUND)
1590 DUNE_THROW(InvalidIndexSetState, "Index has to be in ground mode for repairing pointers to indices");
1591#endif
1592 typedef typename ParallelIndexSet::const_iterator IndexIterator;
1593 typedef typename GlobalList::const_iterator GlobalIterator;
1594 typedef typename RemoteIndexList::iterator Iterator;
1595 GlobalIterator giter = glist_.begin();
1596 IndexIterator index = indexSet_->begin();
1597
1598 for(Iterator iter=rList_->begin(); iter != end_; ++iter) {
1599 while(index->global()<*giter) {
1600 ++index;
1601#ifdef DUNE_ISTL_WITH_CHECKING
1602 if(index == indexSet_->end())
1603 DUNE_THROW(InvalidPosition, "No such global index in set!");
1604#endif
1605 }
1606
1607#ifdef DUNE_ISTL_WITH_CHECKING
1608 if(index->global() != *giter)
1609 DUNE_THROW(InvalidPosition, "No such global index in set!");
1610#endif
1611 iter->localIndex_ = &(*index);
1612 }
1613 }
1614 }
1615
1616 template<typename T, typename A, bool mode>
1618 {
1619 static_assert(!mode,"Not allowed if the mode indicates that new indices"
1620 "might be added to the underlying index set. Use "
1621 "insert(const RemoteIndex&, const GlobalIndex&) instead");
1622
1623#ifdef DUNE_ISTL_WITH_CHECKING
1624 if(!first_ && index.localIndexPair().global()<last_)
1625 DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
1626#endif
1627 // Move to the correct position
1628 while(iter_ != end_ && iter_->localIndexPair().global() < index.localIndexPair().global()) {
1629 ++iter_;
1630 }
1631
1632 // No duplicate entries allowed
1633 assert(iter_==end_ || iter_->localIndexPair().global() != index.localIndexPair().global());
1634 iter_.insert(index);
1635 last_ = index.localIndexPair().global();
1636 first_ = false;
1637 }
1638
1639 template<typename T, typename A, bool mode>
1641 {
1642 static_assert(mode,"Not allowed if the mode indicates that no new indices"
1643 "might be added to the underlying index set. Use "
1644 "insert(const RemoteIndex&) instead");
1645#ifdef DUNE_ISTL_WITH_CHECKING
1646 if(!first_ && global<last_)
1647 DUNE_THROW(InvalidPosition, "Modification of remote indices have to occur with ascending global index.");
1648#endif
1649 // Move to the correct position
1650 while(iter_ != end_ && *giter_ < global) {
1651 ++giter_;
1652 ++iter_;
1653 }
1654
1655 // No duplicate entries allowed
1656 assert(iter_->localIndexPair().global() != global);
1657 iter_.insert(index);
1658 giter_.insert(global);
1659
1660 last_ = global;
1661 first_ = false;
1662 }
1663
1664 template<typename T, typename A, bool mode>
1666 {
1667#ifdef DUNE_ISTL_WITH_CHECKING
1668 if(!first_ && global<last_)
1669 DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
1670#endif
1671
1672 bool found= false;
1673
1674 if(MODIFYINDEXSET) {
1675 // Move to the correct position
1676 while(iter_!=end_ && *giter_< global) {
1677 ++giter_;
1678 ++iter_;
1679 }
1680 if(*giter_ == global) {
1681 giter_.remove();
1682 iter_.remove();
1683 found=true;
1684 }
1685 }else{
1686 while(iter_!=end_ && iter_->localIndexPair().global() < global)
1687 ++iter_;
1688
1689 if(iter_->localIndexPair().global()==global) {
1690 iter_.remove();
1691 found = true;
1692 }
1693 }
1694
1695 last_ = global;
1696 first_ = false;
1697 return found;
1698 }
1699
1700 template<typename T, typename A>
1701 template<bool send>
1703 {
1704 return CollectiveIterator<T,A>(remoteIndices_, send);
1705 }
1706
1707 template<typename T, typename A>
1708 inline MPI_Comm RemoteIndices<T,A>::communicator() const
1709 {
1710 return comm_;
1711
1712 }
1713
1714 template<typename T, typename A>
1716 {
1717 typedef typename RemoteIndexMap::const_iterator const_iterator;
1718
1719 const const_iterator end=pmap.end();
1720 for(const_iterator process=pmap.begin(); process != end; ++process) {
1721 const RemoteIndexList* list = send ? process->second.first : process->second.second;
1723 map_.insert(std::make_pair(process->first,
1724 std::pair<iterator, const iterator>(list->begin(), list->end())));
1725 }
1726 }
1727
1728 template<typename T, typename A>
1729 inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index)
1730 {
1731 typedef typename Map::iterator iterator;
1732 typedef typename Map::const_iterator const_iterator;
1733 const const_iterator end = map_.end();
1734
1735 for(iterator iter = map_.begin(); iter != end;) {
1736 // Step the iterator until we are >= index
1737 typename RemoteIndexList::const_iterator current = iter->second.first;
1738 typename RemoteIndexList::const_iterator rend = iter->second.second;
1739 RemoteIndex remoteIndex;
1740 if(current != rend)
1741 remoteIndex = *current;
1742
1743 while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
1744 ++(iter->second.first);
1745
1746 // erase from the map if there are no more entries.
1747 if(iter->second.first == iter->second.second)
1748 map_.erase(iter++);
1749 else{
1750 ++iter;
1751 }
1752 }
1753 index_=index;
1754 noattribute=true;
1755 }
1756
1757 template<typename T, typename A>
1758 inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index,
1759 const Attribute& attribute)
1760 {
1761 typedef typename Map::iterator iterator;
1762 typedef typename Map::const_iterator const_iterator;
1763 const const_iterator end = map_.end();
1764
1765 for(iterator iter = map_.begin(); iter != end;) {
1766 // Step the iterator until we are >= index
1767 typename RemoteIndexList::const_iterator current = iter->second.first;
1768 typename RemoteIndexList::const_iterator rend = iter->second.second;
1769 RemoteIndex remoteIndex;
1770 if(current != rend)
1771 remoteIndex = *current;
1772
1773 // Move to global index or bigger
1774 while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
1775 ++(iter->second.first);
1776
1777 // move to attribute or bigger
1778 while(iter->second.first!=iter->second.second
1779 && iter->second.first->localIndexPair().global()==index
1780 && iter->second.first->localIndexPair().local().attribute()<attribute)
1781 ++(iter->second.first);
1782
1783 // erase from the map if there are no more entries.
1784 if(iter->second.first == iter->second.second)
1785 map_.erase(iter++);
1786 else{
1787 ++iter;
1788 }
1789 }
1790 index_=index;
1791 attribute_=attribute;
1792 noattribute=false;
1793 }
1794
1795 template<typename T, typename A>
1797 {
1798 typedef typename Map::iterator iterator;
1799 typedef typename Map::const_iterator const_iterator;
1800 const const_iterator end = map_.end();
1801
1802 for(iterator iter = map_.begin(); iter != end;) {
1803 // Step the iterator until we are >= index
1804 typename RemoteIndexList::const_iterator current = iter->second.first;
1805 typename RemoteIndexList::const_iterator rend = iter->second.second;
1806
1807 // move all iterators pointing to the current global index to next value
1808 if(iter->second.first->localIndexPair().global()==index_ &&
1809 (noattribute || iter->second.first->localIndexPair().local().attribute() == attribute_))
1810 ++(iter->second.first);
1811
1812 // erase from the map if there are no more entries.
1813 if(iter->second.first == iter->second.second)
1814 map_.erase(iter++);
1815 else{
1816 ++iter;
1817 }
1818 }
1819 return *this;
1820 }
1821
1822 template<typename T, typename A>
1824 {
1825 return map_.empty();
1826 }
1827
1828 template<typename T, typename A>
1829 inline typename CollectiveIterator<T,A>::iterator
1831 {
1832 if(noattribute)
1833 return iterator(map_.begin(), map_.end(), index_);
1834 else
1835 return iterator(map_.begin(), map_.end(), index_,
1836 attribute_);
1837 }
1838
1839 template<typename T, typename A>
1840 inline typename CollectiveIterator<T,A>::iterator
1841 CollectiveIterator<T,A>::end()
1842 {
1843 return iterator(map_.end(), map_.end(), index_);
1844 }
1845
1846 template<typename TG, typename TA>
1847 inline std::ostream& operator<<(std::ostream& os, const RemoteIndex<TG,TA>& index)
1848 {
1849 os<<"[global="<<index.localIndexPair().global()<<", remote attribute="<<index.attribute()<<" local attribute="<<index.localIndexPair().local().attribute()<<"]";
1850 return os;
1851 }
1852
1853 template<typename T, typename A>
1854 inline std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices)
1855 {
1856 int rank;
1857 MPI_Comm_rank(indices.comm_, &rank);
1858
1859 typedef typename RemoteIndices<T,A>::RemoteIndexList RList;
1860 typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
1861
1862 const const_iterator rend = indices.remoteIndices_.end();
1863
1864 for(const_iterator rindex = indices.remoteIndices_.begin(); rindex!=rend; ++rindex) {
1865 os<<rank<<": Prozess "<<rindex->first<<":";
1866
1867 if(!rindex->second.first->empty()) {
1868 os<<" send:";
1869
1870 const typename RList::const_iterator send= rindex->second.first->end();
1871
1872 for(typename RList::const_iterator index = rindex->second.first->begin();
1873 index != send; ++index)
1874 os<<*index<<" ";
1875 os<<std::endl;
1876 }
1877 if(!rindex->second.second->empty()) {
1878 os<<rank<<": Prozess "<<rindex->first<<": "<<"receive: ";
1879
1880 for(const auto& index : *(rindex->second.second))
1881 os << index << " ";
1882 }
1883 os<<std::endl<<std::flush;
1884 }
1885 return os;
1886 }
1888}
1889
1890#endif // HAVE_MPI
1891
1892#endif
Iterator over the valid underlying iterators.
Definition: remoteindices.hh:792
iterator(const RealIterator &iter, const ConstRealIterator &end, GlobalIndex &index)
Definition: remoteindices.hh:799
iterator(const iterator &other)
Definition: remoteindices.hh:817
const RemoteIndex & operator*() const
Definition: remoteindices.hh:838
iterator & operator++()
Definition: remoteindices.hh:822
const RemoteIndex * operator->() const
Definition: remoteindices.hh:850
bool operator==(const iterator &other) const
Definition: remoteindices.hh:856
int process() const
Definition: remoteindices.hh:844
bool operator!=(const iterator &other) const
Definition: remoteindices.hh:862
A collective iterator for moving over the remote indices for all processes collectively.
Definition: remoteindices.hh:708
CollectiveIterator(const RemoteIndexMap &map_, bool send)
Constructor.
Definition: remoteindices.hh:1715
bool empty()
Checks whether there are still iterators in the map.
Definition: remoteindices.hh:1823
void advance(const GlobalIndex &global)
Advances all underlying iterators.
Definition: remoteindices.hh:1729
std::map< int, std::pair< RemoteIndexList *, RemoteIndexList * > > RemoteIndexMap
The type of the map from rank to remote index list.
Definition: remoteindices.hh:748
A constant random access iterator for the Dune::ArrayList class.
Definition: arraylist.hh:377
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:549
void repairLocalIndexPointers()
Repair the pointers to the local index pairs.
Definition: remoteindices.hh:1584
@ MODIFYINDEXSET
If true the index set corresponding to the remote indices might get modified.
Definition: remoteindices.hh:567
Dune::SLList< RemoteIndex, Allocator > RemoteIndexList
The type of the remote index list.
Definition: remoteindices.hh:602
A Allocator
The type of the allocator for the remote index list.
Definition: remoteindices.hh:598
void insert(const RemoteIndex &index)
Insert an index to the list.
Definition: remoteindices.hh:1617
ParallelIndexSet::GlobalIndex GlobalIndex
The type of the global index.
Definition: remoteindices.hh:578
ParallelIndexSet::LocalIndex LocalIndex
The type of the local index.
Definition: remoteindices.hh:583
RemoteIndexList::const_iterator ConstIterator
The type of the remote index list iterator.
Definition: remoteindices.hh:612
SLListModifyIterator< RemoteIndex, Allocator > ModifyIterator
The type of the modifying iterator of the remote index list.
Definition: remoteindices.hh:607
bool remove(const GlobalIndex &global)
Remove a remote index.
Definition: remoteindices.hh:1665
T ParallelIndexSet
Type of the index set we use.
Definition: remoteindices.hh:573
RemoteIndexListModifier()
Default constructor.
Definition: remoteindices.hh:676
LocalIndex::Attribute Attribute
The type of the attribute.
Definition: remoteindices.hh:588
Dune::RemoteIndex< GlobalIndex, Attribute > RemoteIndex
Type of the remote indices we manage.
Definition: remoteindices.hh:593
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:946
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:952
RemoteIndex()
Parameterless Constructor.
Definition: remoteindices.hh:930
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:82
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:984
void free()
Free the index lists.
Definition: remoteindices.hh:1434
ParallelIndexSet::GlobalIndex GlobalIndex
The type of the global index.
Definition: remoteindices.hh:213
void rebuild()
Rebuilds the set of remote indices.
Definition: remoteindices.hh:1459
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:1708
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:1702
void setIncludeSelf(bool includeSelf)
Tell whether sending from indices of the processor to other indices on the same processor is enabled ...
Definition: remoteindices.hh:971
const_iterator end() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1535
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:1486
const ParallelIndexSet & sourceIndexSet() const
Get the index set at the source.
Definition: remoteindices.hh:999
~RemoteIndices()
Destructor.
Definition: remoteindices.hh:1014
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:1452
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:1007
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:958
typename std::allocator_traits< A >::template rebind_alloc< RemoteIndex > Allocator
The type of the allocator for the remote index list.
Definition: remoteindices.hh:235
const_iterator find(int proc) const
Find an iterator over the remote index lists of a specific process.
Definition: remoteindices.hh:1521
bool isSynced() const
Checks whether the remote indices are synced with the indexsets.
Definition: remoteindices.hh:1479
const_iterator begin() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1528
ParallelIndexSet::LocalIndex LocalIndex
The type of the local index.
Definition: remoteindices.hh:219
A constant iterator for the SLList.
Definition: sllist.hh:369
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:297
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.
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:780
SLListConstIterator< RemoteIndex, Allocator > 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:768
#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 equality.
Definition: iteratorfacades.hh:235
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:257
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:14
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.111.3 (Jul 15, 22:36, 2024)