Dune Core Modules (2.4.2)

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#include "indexset.hh"
7#include "plocalindex.hh"
10#include <dune/common/sllist.hh>
12#include <map>
13#include <set>
14#include <utility>
15#include <iostream>
16#include <algorithm>
17#include <iterator>
18#if HAVE_MPI
19#include "mpitraits.hh"
20#include <mpi.h>
21
22namespace Dune {
34 template<typename TG, typename TA>
36 {
37 public:
38 inline static MPI_Datatype getType();
39 private:
40 static MPI_Datatype type;
41 };
42
43
44 template<typename T, typename A>
45 class RemoteIndices;
46
47 template<typename T1, typename T2>
48 class RemoteIndex;
49
50 template<typename T>
51 class IndicesSyncer;
52
53 template<typename T1, typename T2>
54 std::ostream& operator<<(std::ostream& os, const RemoteIndex<T1,T2>& index);
55
56
57 template<typename T, typename A, bool mode>
59
60
64 template<typename T1, typename T2>
66 {
67 template<typename T>
68 friend class IndicesSyncer;
69
70 template<typename T, typename A, typename A1>
71 friend void repairLocalIndexPointers(std::map<int,SLList<std::pair<typename T::GlobalIndex, typename T::LocalIndex::Attribute>,A> >&,
73 const T&);
74
75 template<typename T, typename A, bool mode>
76 friend class RemoteIndexListModifier;
77
78 public:
83 typedef T1 GlobalIndex;
92 typedef T2 Attribute;
93
99
104 const Attribute attribute() const;
105
111 const PairType& localIndexPair() const;
112
116 RemoteIndex();
117
118
124 RemoteIndex(const T2& attribute,
125 const PairType* local);
126
127
133 RemoteIndex(const T2& attribute);
134
135 bool operator==(const RemoteIndex& ri) const;
136
137 bool operator!=(const RemoteIndex& ri) const;
138 private:
140 const PairType* localIndex_;
141
143 char attribute_;
144 };
145
146 template<class T, class A>
147 std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices);
148
149 class InterfaceBuilder;
150
151 template<class T, class A>
152 class CollectiveIterator;
153
154 template<class T>
155 class IndicesSyncer;
156
157 // forward declaration needed for friend declaration.
158 template<typename T1, typename T2>
160
161
178 template<class T, class A=std::allocator<RemoteIndex<typename T::GlobalIndex,
179 typename T::LocalIndex::Attribute> > >
181 {
182 friend class InterfaceBuilder;
183 friend class IndicesSyncer<T>;
184 template<typename T1, typename A2, typename A1>
185 friend void repairLocalIndexPointers(std::map<int,SLList<std::pair<typename T1::GlobalIndex, typename T1::LocalIndex::Attribute>,A2> >&,
187 const T1&);
188
189 template<class G, class T1, class T2>
190 friend void fillIndexSetHoles(const G& graph, Dune::OwnerOverlapCopyCommunication<T1,T2>& oocomm);
191 friend std::ostream& operator<<<>(std::ostream&, const RemoteIndices<T>&);
192
193 public:
194
199
203
208
209
214
218 typedef typename LocalIndex::Attribute Attribute;
219
224
225
229 typedef typename A::template rebind<RemoteIndex>::other Allocator;
230
234
236 typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
238
239 typedef typename RemoteIndexMap::const_iterator const_iterator;
240
258 inline RemoteIndices(const ParallelIndexSet& source, const ParallelIndexSet& destination,
259 const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>(), bool includeSelf=false);
260
262
270 void setIncludeSelf(bool includeSelf);
271
288 void setIndexSets(const ParallelIndexSet& source, const ParallelIndexSet& destination,
289 const MPI_Comm& comm, const std::vector<int>& neighbours=std::vector<int>());
290
291 template<typename C>
292 void setNeighbours(const C& neighbours)
293 {
294 neighbourIds.clear();
295 neighbourIds.insert(neighbours.begin(), neighbours.end());
296
297 }
298
299 const std::set<int>& getNeighbours() const
300 {
301 return neighbourIds;
302 }
303
308
318 template<bool ignorePublic>
319 void rebuild();
320
321 bool operator==(const RemoteIndices& ri);
322
330 inline bool isSynced() const;
331
335 inline MPI_Comm communicator() const;
336
351 template<bool mode, bool send>
353
360 inline const_iterator find(int proc) const;
361
366 inline const_iterator begin() const;
367
372 inline const_iterator end() const;
373
377 template<bool send>
379
383 inline void free();
384
389 inline int neighbours() const;
390
392 inline const ParallelIndexSet& sourceIndexSet() const;
393
396
397 private:
400 {}
401
403 const ParallelIndexSet* source_;
404
406 const ParallelIndexSet* target_;
407
409 MPI_Comm comm_;
410
413 std::set<int> neighbourIds;
414
416 const static int commTag_=333;
417
422 int sourceSeqNo_;
423
428 int destSeqNo_;
429
433 bool publicIgnored;
434
438 bool firstBuild;
439
440 /*
441 * @brief If true, sending from indices of the processor to other
442 * indices on the same processor is enabled even if the same indexset is used
443 * on both the
444 * sending and receiving side.
445 */
446 bool includeSelf;
447
450 PairType;
451
458 RemoteIndexMap remoteIndices_;
459
470 template<bool ignorePublic>
471 inline void buildRemote(bool includeSelf);
472
478 inline int noPublic(const ParallelIndexSet& indexSet);
479
491 template<bool ignorePublic>
492 inline void packEntries(PairType** myPairs, const ParallelIndexSet& indexSet,
493 char* p_out, MPI_Datatype type, int bufferSize,
494 int* position, int n);
495
509 inline void unpackIndices(RemoteIndexList& remote, int remoteEntries,
510 PairType** local, int localEntries, char* p_in,
511 MPI_Datatype type, int* positon, int bufferSize,
512 bool fromOurself);
513
514 inline void unpackIndices(RemoteIndexList& send, RemoteIndexList& receive,
515 int remoteEntries, PairType** localSource,
516 int localSourceEntries, PairType** localDest,
517 int localDestEntries, char* p_in,
518 MPI_Datatype type, int* position, int bufferSize);
519
520 void unpackCreateRemote(char* p_in, PairType** sourcePairs, PairType** DestPairs,
521 int remoteProc, int sourcePublish, int destPublish,
522 int bufferSize, bool sendTwo, bool fromOurSelf=false);
523 };
524
542 template<class T, class A, bool mode>
544 {
545
546 template<typename T1, typename A1>
547 friend class RemoteIndices;
548
549 public:
550 class InvalidPosition : public RangeError
551 {};
552
553 enum {
562 MODIFYINDEXSET=mode
563 };
564
569
574
579
583 typedef typename LocalIndex::Attribute Attribute;
584
589
593 typedef A Allocator;
594
598
603
608
622 void insert(const RemoteIndex& index) throw(InvalidPosition);
623
624
639 void insert(const RemoteIndex& index, const GlobalIndex& global) throw(InvalidPosition);
640
648 bool remove(const GlobalIndex& global) throw(InvalidPosition);
649
663
664
666
672 : glist_()
673 {}
674
675 private:
676
683 RemoteIndexList& rList);
684
685 typedef SLList<GlobalIndex,Allocator> GlobalList;
686 typedef typename GlobalList::ModifyIterator GlobalModifyIterator;
687 RemoteIndexList* rList_;
688 const ParallelIndexSet* indexSet_;
689 GlobalList glist_;
690 ModifyIterator iter_;
691 GlobalModifyIterator giter_;
692 ConstIterator end_;
693 bool first_;
694 GlobalIndex last_;
695 };
696
701 template<class T, class A>
703 {
704
708 typedef T ParallelIndexSet;
709
713 typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
714
718 typedef typename ParallelIndexSet::LocalIndex LocalIndex;
719
723 typedef typename LocalIndex::Attribute Attribute;
724
727
729 typedef typename A::template rebind<RemoteIndex>::other Allocator;
730
733
735 typedef std::map<int,std::pair<typename RemoteIndexList::const_iterator,
736 const typename RemoteIndexList::const_iterator> >
737 Map;
738
739 public:
740
742 typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
744
750 inline CollectiveIterator(const RemoteIndexMap& map_, bool send);
751
760 inline void advance(const GlobalIndex& global);
761
771 inline void advance(const GlobalIndex& global, const Attribute& attribute);
772
773 CollectiveIterator& operator++();
774
778 inline bool empty();
779
787 {
788 public:
789 typedef typename Map::iterator RealIterator;
790 typedef typename Map::iterator ConstRealIterator;
791
792
794 iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex& index)
795 : iter_(iter), end_(end), index_(index), hasAttribute(false)
796 {
797 // Move to the first valid entry
798 while(iter_!=end_ && iter_->second.first->localIndexPair().global()!=index_)
799 ++iter_;
800 }
801
802 iterator(const RealIterator& iter, const ConstRealIterator& end, GlobalIndex index,
803 Attribute attribute)
804 : iter_(iter), end_(end), index_(index), attribute_(attribute), hasAttribute(true)
805 {
806 // Move to the first valid entry or the end
807 while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_
808 || iter_->second.first->localIndexPair().local().attribute()!=attribute))
809 ++iter_;
810 }
812 iterator(const iterator& other)
813 : iter_(other.iter_), end_(other.end_), index_(other.index_)
814 { }
815
818 {
819 ++iter_;
820 // If entry is not valid move on
821 while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_ ||
822 (hasAttribute &&
823 iter_->second.first->localIndexPair().local().attribute()!=attribute_)))
824 ++iter_;
825 assert(iter_==end_ ||
826 (iter_->second.first->localIndexPair().global()==index_));
827 assert(iter_==end_ || !hasAttribute ||
828 (iter_->second.first->localIndexPair().local().attribute()==attribute_));
829 return *this;
830 }
831
833 const RemoteIndex& operator*() const
834 {
835 return *(iter_->second.first);
836 }
837
839 int process() const
840 {
841 return iter_->first;
842 }
843
845 const RemoteIndex* operator->() const
846 {
847 return iter_->second.first.operator->();
848 }
849
851 bool operator==(const iterator& other)
852 {
853 return other.iter_==iter_;
854 }
855
857 bool operator!=(const iterator& other)
858 {
859 return other.iter_!=iter_;
860 }
861
862 private:
863 iterator();
864
865 RealIterator iter_;
866 RealIterator end_;
867 GlobalIndex index_;
868 Attribute attribute_;
869 bool hasAttribute;
870 };
871
872 iterator begin();
873
874 iterator end();
875
876 private:
877
878 Map map_;
879 GlobalIndex index_;
880 Attribute attribute_;
881 bool noattribute;
882 };
883
884 template<typename TG, typename TA>
885 MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::getType()
886 {
887 if(type==MPI_DATATYPE_NULL) {
888 int length[4];
889 MPI_Aint disp[4];
890 MPI_Datatype types[4] = {MPI_LB, MPITraits<TG>::getType(),
891 MPITraits<ParallelLocalIndex<TA> >::getType(), MPI_UB};
892 IndexPair<TG,ParallelLocalIndex<TA> > rep[2];
893 length[0]=length[1]=length[2]=length[3]=1;
894 MPI_Get_address(rep, disp); // lower bound of the datatype
895 MPI_Get_address(&(rep[0].global_), disp+1);
896 MPI_Get_address(&(rep[0].local_), disp+2);
897 MPI_Get_address(rep+1, disp+3); // upper bound of the datatype
898 for(int i=3; i >= 0; --i)
899 disp[i] -= disp[0];
900 MPI_Type_create_struct(4, length, disp, types, &type);
901 MPI_Type_commit(&type);
902 }
903 return type;
904 }
905
906 template<typename TG, typename TA>
907 MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::type=MPI_DATATYPE_NULL;
908
909 template<typename T1, typename T2>
910 RemoteIndex<T1,T2>::RemoteIndex(const T2& attribute, const PairType* local)
911 : localIndex_(local), attribute_(attribute)
912 {}
913
914 template<typename T1, typename T2>
916 : localIndex_(0), attribute_(attribute)
917 {}
918
919 template<typename T1, typename T2>
921 : localIndex_(0), attribute_()
922 {}
923 template<typename T1, typename T2>
924 inline bool RemoteIndex<T1,T2>::operator==(const RemoteIndex& ri) const
925 {
926 return localIndex_==ri.localIndex_ && attribute_==ri.attribute;
927 }
928
929 template<typename T1, typename T2>
930 inline bool RemoteIndex<T1,T2>::operator!=(const RemoteIndex& ri) const
931 {
932 return localIndex_!=ri.localIndex_ || attribute_!=ri.attribute_;
933 }
934
935 template<typename T1, typename T2>
936 inline const T2 RemoteIndex<T1,T2>::attribute() const
937 {
938 return T2(attribute_);
939 }
940
941 template<typename T1, typename T2>
943 {
944 return *localIndex_;
945 }
946
947 template<typename T, typename A>
949 const ParallelIndexSet& destination,
950 const MPI_Comm& comm,
951 const std::vector<int>& neighbours,
952 bool includeSelf_)
953 : source_(&source), target_(&destination), comm_(comm),
954 sourceSeqNo_(-1), destSeqNo_(-1), publicIgnored(false), firstBuild(true),
955 includeSelf(includeSelf_)
956 {
957 setNeighbours(neighbours);
958 }
959
960 template<typename T, typename A>
962 {
963 includeSelf=b;
964 }
965
966 template<typename T, typename A>
968 : source_(0), target_(0), sourceSeqNo_(-1),
969 destSeqNo_(-1), publicIgnored(false), firstBuild(true)
970 {}
971
972 template<class T, typename A>
974 const ParallelIndexSet& destination,
975 const MPI_Comm& comm,
976 const std::vector<int>& neighbours)
977 {
978 free();
979 source_ = &source;
980 target_ = &destination;
981 comm_ = comm;
982 firstBuild = true;
983 setNeighbours(neighbours);
984 }
985
986 template<typename T, typename A>
989 {
990 return *source_;
991 }
992
993
994 template<typename T, typename A>
997 {
998 return *target_;
999 }
1000
1001
1002 template<typename T, typename A>
1004 {
1005 free();
1006 }
1007
1008 template<typename T, typename A>
1009 template<bool ignorePublic>
1011 const ParallelIndexSet& indexSet,
1012 char* p_out, MPI_Datatype type,
1013 int bufferSize,
1014 int *position, int n)
1015 {
1017 // fill with own indices
1018 typedef typename ParallelIndexSet::const_iterator const_iterator;
1019 typedef IndexPair<GlobalIndex,LocalIndex> PairType;
1020 const const_iterator end = indexSet.end();
1021
1022 //Now pack the source indices
1023 int i=0;
1024 for(const_iterator index = indexSet.begin(); index != end; ++index)
1025 if(ignorePublic || index->local().isPublic()) {
1026
1027 MPI_Pack(const_cast<PairType*>(&(*index)), 1,
1028 type,
1029 p_out, bufferSize, position, comm_);
1030 pairs[i++] = const_cast<PairType*>(&(*index));
1031
1032 }
1033 assert(i==n);
1034 }
1035
1036 template<typename T, typename A>
1037 inline int RemoteIndices<T,A>::noPublic(const ParallelIndexSet& indexSet)
1038 {
1039 typedef typename ParallelIndexSet::const_iterator const_iterator;
1040
1041 int noPublic=0;
1042
1043 const const_iterator end=indexSet.end();
1044 for(const_iterator index=indexSet.begin(); index!=end; ++index)
1045 if(index->local().isPublic())
1046 noPublic++;
1047
1048 return noPublic;
1049
1050 }
1051
1052
1053 template<typename T, typename A>
1054 inline void RemoteIndices<T,A>::unpackCreateRemote(char* p_in, PairType** sourcePairs,
1055 PairType** destPairs, int remoteProc,
1056 int sourcePublish, int destPublish,
1057 int bufferSize, bool sendTwo,
1058 bool fromOurSelf)
1059 {
1060
1061 // unpack the number of indices we received
1062 int noRemoteSource=-1, noRemoteDest=-1;
1063 char twoIndexSets=0;
1064 int position=0;
1065 // Did we receive two index sets?
1066 MPI_Unpack(p_in, bufferSize, &position, &twoIndexSets, 1, MPI_CHAR, comm_);
1067 // The number of source indices received
1068 MPI_Unpack(p_in, bufferSize, &position, &noRemoteSource, 1, MPI_INT, comm_);
1069 // The number of destination indices received
1070 MPI_Unpack(p_in, bufferSize, &position, &noRemoteDest, 1, MPI_INT, comm_);
1071
1072
1073 // Indices for which we receive
1074 RemoteIndexList* receive= new RemoteIndexList();
1075 // Indices for which we send
1076 RemoteIndexList* send=0;
1077
1078 MPI_Datatype type= MPITraits<PairType>::getType();
1079
1080 if(!twoIndexSets) {
1081 if(sendTwo) {
1082 send = new RemoteIndexList();
1083 // Create both remote index sets simultaneously
1084 unpackIndices(*send, *receive, noRemoteSource, sourcePairs, sourcePublish,
1085 destPairs, destPublish, p_in, type, &position, bufferSize);
1086 }else{
1087 // we only need one list
1088 unpackIndices(*receive, noRemoteSource, sourcePairs, sourcePublish,
1089 p_in, type, &position, bufferSize, fromOurSelf);
1090 send=receive;
1091 }
1092 }else{
1093
1094 int oldPos=position;
1095 // Two index sets received
1096 unpackIndices(*receive, noRemoteSource, destPairs, destPublish,
1097 p_in, type, &position, bufferSize, fromOurSelf);
1098 if(!sendTwo)
1099 //unpack source entries again as destination entries
1100 position=oldPos;
1101
1102 send = new RemoteIndexList();
1103 unpackIndices(*send, noRemoteDest, sourcePairs, sourcePublish,
1104 p_in, type, &position, bufferSize, fromOurSelf);
1105 }
1106
1107 if(receive->empty() && send->empty()) {
1108 if(send==receive) {
1109 delete send;
1110 }else{
1111 delete send;
1112 delete receive;
1113 }
1114 }else{
1115 remoteIndices_.insert(std::make_pair(remoteProc,
1116 std::make_pair(send,receive)));
1117 }
1118 }
1119
1120
1121 template<typename T, typename A>
1122 template<bool ignorePublic>
1123 inline void RemoteIndices<T,A>::buildRemote(bool includeSelf_)
1124 {
1125 // Processor configuration
1126 int rank, procs;
1127 MPI_Comm_rank(comm_, &rank);
1128 MPI_Comm_size(comm_, &procs);
1129
1130 // number of local indices to publish
1131 // The indices of the destination will be send.
1132 int sourcePublish, destPublish;
1133
1134 // Do we need to send two index sets?
1135 char sendTwo = (source_ != target_);
1136
1137 if(procs==1 && !(sendTwo || includeSelf_))
1138 // Nothing to communicate
1139 return;
1140
1141 sourcePublish = (ignorePublic) ? source_->size() : noPublic(*source_);
1142
1143 if(sendTwo)
1144 destPublish = (ignorePublic) ? target_->size() : noPublic(*target_);
1145 else
1146 // we only need to send one set of indices
1147 destPublish = 0;
1148
1149 int maxPublish, publish=sourcePublish+destPublish;
1150
1151 // Calucate maximum number of indices send
1152 MPI_Allreduce(&publish, &maxPublish, 1, MPI_INT, MPI_MAX, comm_);
1153
1154 // allocate buffers
1155 typedef IndexPair<GlobalIndex,LocalIndex> PairType;
1156
1157 PairType** destPairs;
1158 PairType** sourcePairs = new PairType*[sourcePublish>0 ? sourcePublish : 1];
1159
1160 if(sendTwo)
1161 destPairs = new PairType*[destPublish>0 ? destPublish : 1];
1162 else
1163 destPairs=sourcePairs;
1164
1165 char** buffer = new char*[2];
1166 int bufferSize;
1167 int position=0;
1168 int intSize;
1169 int charSize;
1170
1171 // calculate buffer size
1172 MPI_Datatype type = MPITraits<PairType>::getType();
1173
1174 MPI_Pack_size(maxPublish, type, comm_,
1175 &bufferSize);
1176 MPI_Pack_size(1, MPI_INT, comm_,
1177 &intSize);
1178 MPI_Pack_size(1, MPI_CHAR, comm_,
1179 &charSize);
1180 // Our message will contain the following:
1181 // a bool wether two index sets where sent
1182 // the size of the source and the dest indexset,
1183 // then the source and destination indices
1184 bufferSize += 2 * intSize + charSize;
1185
1186 if(bufferSize<=0) bufferSize=1;
1187
1188 buffer[0] = new char[bufferSize];
1189 buffer[1] = new char[bufferSize];
1190
1191
1192 // pack entries into buffer[0], p_out below!
1193 MPI_Pack(&sendTwo, 1, MPI_CHAR, buffer[0], bufferSize, &position,
1194 comm_);
1195
1196 // The number of indices we send for each index set
1197 MPI_Pack(&sourcePublish, 1, MPI_INT, buffer[0], bufferSize, &position,
1198 comm_);
1199 MPI_Pack(&destPublish, 1, MPI_INT, buffer[0], bufferSize, &position,
1200 comm_);
1201
1202 // Now pack the source indices and setup the destination pairs
1203 packEntries<ignorePublic>(sourcePairs, *source_, buffer[0], type,
1204 bufferSize, &position, sourcePublish);
1205 // If necessary send the dest indices and setup the source pairs
1206 if(sendTwo)
1207 packEntries<ignorePublic>(destPairs, *target_, buffer[0], type,
1208 bufferSize, &position, destPublish);
1209
1210
1211 // Update remote indices for ourself
1212 if(sendTwo|| includeSelf_)
1213 unpackCreateRemote(buffer[0], sourcePairs, destPairs, rank, sourcePublish,
1214 destPublish, bufferSize, sendTwo, includeSelf_);
1215
1216 neighbourIds.erase(rank);
1217
1218 if(neighbourIds.size()==0)
1219 {
1220 Dune::dvverb<<rank<<": Sending messages in a ring"<<std::endl;
1221 // send messages in ring
1222 for(int proc=1; proc<procs; proc++) {
1223 // pointers to the current input and output buffers
1224 char* p_out = buffer[1-(proc%2)];
1225 char* p_in = buffer[proc%2];
1226
1227 MPI_Status status;
1228 if(rank%2==0) {
1229 MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
1230 commTag_, comm_);
1231 MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
1232 commTag_, comm_, &status);
1233 }else{
1234 MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
1235 commTag_, comm_, &status);
1236 MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
1237 commTag_, comm_);
1238 }
1239
1240
1241 // The process these indices are from
1242 int remoteProc = (rank+procs-proc)%procs;
1243
1244 unpackCreateRemote(p_in, sourcePairs, destPairs, remoteProc, sourcePublish,
1245 destPublish, bufferSize, sendTwo);
1246
1247 }
1248
1249 }
1250 else
1251 {
1252 MPI_Request* requests=new MPI_Request[neighbourIds.size()];
1253 MPI_Request* req=requests;
1254
1255 typedef typename std::set<int>::size_type size_type;
1256 size_type noNeighbours=neighbourIds.size();
1257
1258 // setup sends
1259 for(std::set<int>::iterator neighbour=neighbourIds.begin();
1260 neighbour!= neighbourIds.end(); ++neighbour) {
1261 // Only send the information to the neighbouring processors
1262 MPI_Issend(buffer[0], position , MPI_PACKED, *neighbour, commTag_, comm_, req++);
1263 }
1264
1265 //Test for received messages
1266
1267 for(size_type received=0; received <noNeighbours; ++received)
1268 {
1269 MPI_Status status;
1270 // probe for next message
1271 MPI_Probe(MPI_ANY_SOURCE, commTag_, comm_, &status);
1272 int remoteProc=status.MPI_SOURCE;
1273 int size;
1274 MPI_Get_count(&status, MPI_PACKED, &size);
1275 // receive message
1276 MPI_Recv(buffer[1], size, MPI_PACKED, remoteProc,
1277 commTag_, comm_, &status);
1278
1279 unpackCreateRemote(buffer[1], sourcePairs, destPairs, remoteProc, sourcePublish,
1280 destPublish, bufferSize, sendTwo);
1281 }
1282 // wait for completion of pending requests
1283 MPI_Status* statuses = new MPI_Status[neighbourIds.size()];
1284
1285 if(MPI_ERR_IN_STATUS==MPI_Waitall(neighbourIds.size(), requests, statuses)) {
1286 for(size_type i=0; i < neighbourIds.size(); ++i)
1287 if(statuses[i].MPI_ERROR!=MPI_SUCCESS) {
1288 std::cerr<<rank<<": MPI_Error occurred while receiving message."<<std::endl;
1289 MPI_Abort(comm_, 999);
1290 }
1291 }
1292 delete[] requests;
1293 delete[] statuses;
1294 }
1295
1296
1297 // delete allocated memory
1298 if(destPairs!=sourcePairs)
1299 delete[] destPairs;
1300
1301 delete[] sourcePairs;
1302 delete[] buffer[0];
1303 delete[] buffer[1];
1304 delete[] buffer;
1305 }
1306
1307 template<typename T, typename A>
1308 inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& remote,
1309 int remoteEntries,
1310 PairType** local,
1311 int localEntries,
1312 char* p_in,
1313 MPI_Datatype type,
1314 int* position,
1315 int bufferSize,
1316 bool fromOurSelf)
1317 {
1318 if(remoteEntries==0)
1319 return;
1320
1321 PairType index(1);
1322 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1323 type, comm_);
1324 GlobalIndex oldGlobal=index.global();
1325 int n_in=0, localIndex=0;
1326
1327 //Check if we know the global index
1328 while(localIndex<localEntries) {
1329 if(local[localIndex]->global()==index.global()) {
1330 int oldLocalIndex=localIndex;
1331
1332 while(localIndex<localEntries &&
1333 local[localIndex]->global()==index.global()) {
1334 if(!fromOurSelf || index.local().attribute() !=
1335 local[localIndex]->local().attribute())
1336 // if index is from us it has to have a different attribute
1337 remote.push_back(RemoteIndex(index.local().attribute(),
1338 local[localIndex]));
1339 localIndex++;
1340 }
1341
1342 // unpack next remote index
1343 if((++n_in) < remoteEntries) {
1344 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1345 type, comm_);
1346 if(index.global()==oldGlobal)
1347 // Restart comparison for the same global indices
1348 localIndex=oldLocalIndex;
1349 else
1350 oldGlobal=index.global();
1351 }else{
1352 // No more received indices
1353 break;
1354 }
1355 continue;
1356 }
1357
1358 if (local[localIndex]->global()<index.global()) {
1359 // compare with next entry in our list
1360 ++localIndex;
1361 }else{
1362 // We do not know the index, unpack next
1363 if((++n_in) < remoteEntries) {
1364 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1365 type, comm_);
1366 oldGlobal=index.global();
1367 }else
1368 // No more received indices
1369 break;
1370 }
1371 }
1372
1373 // Unpack the other received indices without doing anything
1374 while(++n_in < remoteEntries)
1375 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1376 type, comm_);
1377 }
1378
1379
1380 template<typename T, typename A>
1381 inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& send,
1382 RemoteIndexList& receive,
1383 int remoteEntries,
1384 PairType** localSource,
1385 int localSourceEntries,
1386 PairType** localDest,
1387 int localDestEntries,
1388 char* p_in,
1389 MPI_Datatype type,
1390 int* position,
1391 int bufferSize)
1392 {
1393 int n_in=0, sourceIndex=0, destIndex=0;
1394
1395 //Check if we know the global index
1396 while(n_in<remoteEntries && (sourceIndex<localSourceEntries || destIndex<localDestEntries)) {
1397 // Unpack next index
1398 PairType index;
1399 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1400 type, comm_);
1401 n_in++;
1402
1403 // Advance until global index in localSource and localDest are >= than the one in the unpacked index
1404 while(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()<index.global())
1405 sourceIndex++;
1406
1407 while(destIndex<localDestEntries && localDest[destIndex]->global()<index.global())
1408 destIndex++;
1409
1410 // Add a remote index if we found the global index.
1411 if(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()==index.global())
1412 send.push_back(RemoteIndex(index.local().attribute(),
1413 localSource[sourceIndex]));
1414
1415 if(destIndex < localDestEntries && localDest[destIndex]->global() == index.global())
1416 receive.push_back(RemoteIndex(index.local().attribute(),
1417 localDest[sourceIndex]));
1418 }
1419
1420 }
1421
1422 template<typename T, typename A>
1424 {
1425 typedef typename RemoteIndexMap::iterator Iterator;
1426 Iterator lend = remoteIndices_.end();
1427 for(Iterator lists=remoteIndices_.begin(); lists != lend; ++lists) {
1428 if(lists->second.first==lists->second.second) {
1429 // there is only one remote index list.
1430 delete lists->second.first;
1431 }else{
1432 delete lists->second.first;
1433 delete lists->second.second;
1434 }
1435 }
1436 remoteIndices_.clear();
1437 firstBuild=true;
1438 }
1439
1440 template<typename T, typename A>
1442 {
1443 return remoteIndices_.size();
1444 }
1445
1446 template<typename T, typename A>
1447 template<bool ignorePublic>
1449 {
1450 // Test wether a rebuild is Needed.
1451 if(firstBuild ||
1452 ignorePublic!=publicIgnored || !
1453 isSynced()) {
1454 free();
1455
1456 buildRemote<ignorePublic>(includeSelf);
1457
1458 sourceSeqNo_ = source_->seqNo();
1459 destSeqNo_ = target_->seqNo();
1460 firstBuild=false;
1461 publicIgnored=ignorePublic;
1462 }
1463
1464
1465 }
1466
1467 template<typename T, typename A>
1469 {
1470 return sourceSeqNo_==source_->seqNo() && destSeqNo_ ==target_->seqNo();
1471 }
1472
1473 template<typename T, typename A>
1474 template<bool mode, bool send>
1476 {
1477
1478 // The user are on their own now!
1479 // We assume they know what they are doing and just set the
1480 // remote indices to synced status.
1481 sourceSeqNo_ = source_->seqNo();
1482 destSeqNo_ = target_->seqNo();
1483
1484 typename RemoteIndexMap::iterator found = remoteIndices_.find(process);
1485
1486 if(found == remoteIndices_.end())
1487 {
1488 if(source_ != target_)
1489 found = remoteIndices_.insert(found, std::make_pair(process,
1490 std::make_pair(new RemoteIndexList(),
1491 new RemoteIndexList())));
1492 else{
1493 RemoteIndexList* rlist = new RemoteIndexList();
1494 found = remoteIndices_.insert(found,
1495 std::make_pair(process,
1496 std::make_pair(rlist, rlist)));
1497 }
1498 }
1499
1500 firstBuild = false;
1501
1502 if(send)
1503 return RemoteIndexListModifier<T,A,mode>(*source_, *(found->second.first));
1504 else
1505 return RemoteIndexListModifier<T,A,mode>(*target_, *(found->second.second));
1506 }
1507
1508 template<typename T, typename A>
1509 inline typename RemoteIndices<T,A>::const_iterator
1511 {
1512 return remoteIndices_.find(proc);
1513 }
1514
1515 template<typename T, typename A>
1516 inline typename RemoteIndices<T,A>::const_iterator
1518 {
1519 return remoteIndices_.begin();
1520 }
1521
1522 template<typename T, typename A>
1523 inline typename RemoteIndices<T,A>::const_iterator
1525 {
1526 return remoteIndices_.end();
1527 }
1528
1529
1530 template<typename T, typename A>
1532 {
1533 if(neighbours()!=ri.neighbours())
1534 return false;
1535
1536 typedef RemoteIndexList RList;
1537 typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
1538
1539 const const_iterator rend = remoteIndices_.end();
1540
1541 for(const_iterator rindex = remoteIndices_.begin(), rindex1=ri.remoteIndices_.begin(); rindex!=rend; ++rindex, ++rindex1) {
1542 if(rindex->first != rindex1->first)
1543 return false;
1544 if(*(rindex->second.first) != *(rindex1->second.first))
1545 return false;
1546 if(*(rindex->second.second) != *(rindex1->second.second))
1547 return false;
1548 }
1549 return true;
1550 }
1551
1552 template<class T, class A, bool mode>
1553 RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(const ParallelIndexSet& indexSet,
1554 RemoteIndexList& rList)
1555 : rList_(&rList), indexSet_(&indexSet), iter_(rList.beginModify()), end_(rList.end()), first_(true)
1556 {
1557 if(MODIFYINDEXSET) {
1558 assert(indexSet_);
1559 for(ConstIterator iter=iter_; iter != end_; ++iter)
1560 glist_.push_back(iter->localIndexPair().global());
1561 giter_ = glist_.beginModify();
1562 }
1563 }
1564
1565 template<typename T, typename A, bool mode>
1566 RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(const RemoteIndexListModifier<T,A,mode>& other)
1567 : rList_(other.rList_), indexSet_(other.indexSet_),
1568 glist_(other.glist_), iter_(other.iter_), giter_(other.giter_), end_(other.end_),
1569 first_(other.first_), last_(other.last_)
1570 {}
1571
1572 template<typename T, typename A, bool mode>
1574 {
1575 if(MODIFYINDEXSET) {
1576 // repair pointers to local index set.
1577#ifdef DUNE_ISTL_WITH_CHECKING
1578 if(indexSet_->state()!=GROUND)
1579 DUNE_THROW(InvalidIndexSetState, "Index has to be in ground mode for repairing pointers to indices");
1580#endif
1581 typedef typename ParallelIndexSet::const_iterator IndexIterator;
1582 typedef typename GlobalList::const_iterator GlobalIterator;
1583 typedef typename RemoteIndexList::iterator Iterator;
1584 GlobalIterator giter = glist_.begin();
1585 IndexIterator index = indexSet_->begin();
1586
1587 for(Iterator iter=rList_->begin(); iter != end_; ++iter) {
1588 while(index->global()<*giter) {
1589 ++index;
1590#ifdef DUNE_ISTL_WITH_CHECKING
1591 if(index == indexSet_->end())
1592 DUNE_THROW(InvalidPosition, "No such global index in set!");
1593#endif
1594 }
1595
1596#ifdef DUNE_ISTL_WITH_CHECKING
1597 if(index->global() != *giter)
1598 DUNE_THROW(InvalidPosition, "No such global index in set!");
1599#endif
1600 iter->localIndex_ = &(*index);
1601 }
1602 }
1603 }
1604
1605 template<typename T, typename A, bool mode>
1606 inline void RemoteIndexListModifier<T,A,mode>::insert(const RemoteIndex& index) throw(InvalidPosition)
1607 {
1608 static_assert(!mode,"Not allowed if the mode indicates that new indices"
1609 "might be added to the underlying index set. Use "
1610 "insert(const RemoteIndex&, const GlobalIndex&) instead");
1611
1612#ifdef DUNE_ISTL_WITH_CHECKING
1613 if(!first_ && index.localIndexPair().global()<last_)
1614 DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
1615#endif
1616 // Move to the correct position
1617 while(iter_ != end_ && iter_->localIndexPair().global() < index.localIndexPair().global()) {
1618 ++iter_;
1619 }
1620
1621 // No duplicate entries allowed
1622 assert(iter_==end_ || iter_->localIndexPair().global() != index.localIndexPair().global());
1623 iter_.insert(index);
1624 last_ = index.localIndexPair().global();
1625 first_ = false;
1626 }
1627
1628 template<typename T, typename A, bool mode>
1629 inline void RemoteIndexListModifier<T,A,mode>::insert(const RemoteIndex& index, const GlobalIndex& global) throw(InvalidPosition)
1630 {
1631 static_assert(mode,"Not allowed if the mode indicates that no new indices"
1632 "might be added to the underlying index set. Use "
1633 "insert(const RemoteIndex&) instead");
1634#ifdef DUNE_ISTL_WITH_CHECKING
1635 if(!first_ && global<last_)
1636 DUNE_THROW(InvalidPosition, "Modification of remote indices have to occur with ascending global index.");
1637#endif
1638 // Move to the correct position
1639 while(iter_ != end_ && *giter_ < global) {
1640 ++giter_;
1641 ++iter_;
1642 }
1643
1644 // No duplicate entries allowed
1645 assert(iter_->localIndexPair().global() != global);
1646 iter_.insert(index);
1647 giter_.insert(global);
1648
1649 last_ = global;
1650 first_ = false;
1651 }
1652
1653 template<typename T, typename A, bool mode>
1654 bool RemoteIndexListModifier<T,A,mode>::remove(const GlobalIndex& global) throw(InvalidPosition)
1655 {
1656#ifdef DUNE_ISTL_WITH_CHECKING
1657 if(!first_ && global<last_)
1658 DUNE_THROW(InvalidPosition, "Modifcation of remote indices have to occur with ascending global index.");
1659#endif
1660
1661 bool found= false;
1662
1663 if(MODIFYINDEXSET) {
1664 // Move to the correct position
1665 while(iter_!=end_ && *giter_< global) {
1666 ++giter_;
1667 ++iter_;
1668 }
1669 if(*giter_ == global) {
1670 giter_.remove();
1671 iter_.remove();
1672 found=true;
1673 }
1674 }else{
1675 while(iter_!=end_ && iter_->localIndexPair().global() < global)
1676 ++iter_;
1677
1678 if(iter_->localIndexPair().global()==global) {
1679 iter_.remove();
1680 found = true;
1681 }
1682 }
1683
1684 last_ = global;
1685 first_ = false;
1686 return found;
1687 }
1688
1689 template<typename T, typename A>
1690 template<bool send>
1692 {
1693 return CollectiveIterator<T,A>(remoteIndices_, send);
1694 }
1695
1696 template<typename T, typename A>
1697 inline MPI_Comm RemoteIndices<T,A>::communicator() const
1698 {
1699 return comm_;
1700
1701 }
1702
1703 template<typename T, typename A>
1705 {
1706 typedef typename RemoteIndexMap::const_iterator const_iterator;
1707
1708 const const_iterator end=pmap.end();
1709 for(const_iterator process=pmap.begin(); process != end; ++process) {
1710 const RemoteIndexList* list = send ? process->second.first : process->second.second;
1712 map_.insert(std::make_pair(process->first,
1713 std::pair<iterator, const iterator>(list->begin(), list->end())));
1714 }
1715 }
1716
1717 template<typename T, typename A>
1718 inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index)
1719 {
1720 typedef typename Map::iterator iterator;
1721 typedef typename Map::const_iterator const_iterator;
1722 const const_iterator end = map_.end();
1723
1724 for(iterator iter = map_.begin(); iter != end;) {
1725 // Step the iterator until we are >= index
1726 typename RemoteIndexList::const_iterator current = iter->second.first;
1727 typename RemoteIndexList::const_iterator rend = iter->second.second;
1728 RemoteIndex remoteIndex;
1729 if(current != rend)
1730 remoteIndex = *current;
1731
1732 while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
1733 ++(iter->second.first);
1734
1735 // erase from the map if there are no more entries.
1736 if(iter->second.first == iter->second.second)
1737 map_.erase(iter++);
1738 else{
1739 ++iter;
1740 }
1741 }
1742 index_=index;
1743 noattribute=true;
1744 }
1745
1746 template<typename T, typename A>
1747 inline void CollectiveIterator<T,A>::advance(const GlobalIndex& index,
1748 const Attribute& attribute)
1749 {
1750 typedef typename Map::iterator iterator;
1751 typedef typename Map::const_iterator const_iterator;
1752 const const_iterator end = map_.end();
1753
1754 for(iterator iter = map_.begin(); iter != end;) {
1755 // Step the iterator until we are >= index
1756 typename RemoteIndexList::const_iterator current = iter->second.first;
1757 typename RemoteIndexList::const_iterator rend = iter->second.second;
1758 RemoteIndex remoteIndex;
1759 if(current != rend)
1760 remoteIndex = *current;
1761
1762 // Move to global index or bigger
1763 while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
1764 ++(iter->second.first);
1765
1766 // move to attribute or bigger
1767 while(iter->second.first!=iter->second.second
1768 && iter->second.first->localIndexPair().global()==index
1769 && iter->second.first->localIndexPair().local().attribute()<attribute)
1770 ++(iter->second.first);
1771
1772 // erase from the map if there are no more entries.
1773 if(iter->second.first == iter->second.second)
1774 map_.erase(iter++);
1775 else{
1776 ++iter;
1777 }
1778 }
1779 index_=index;
1780 attribute_=attribute;
1781 noattribute=false;
1782 }
1783
1784 template<typename T, typename A>
1786 {
1787 typedef typename Map::iterator iterator;
1788 typedef typename Map::const_iterator const_iterator;
1789 const const_iterator end = map_.end();
1790
1791 for(iterator iter = map_.begin(); iter != end;) {
1792 // Step the iterator until we are >= index
1793 typename RemoteIndexList::const_iterator current = iter->second.first;
1794 typename RemoteIndexList::const_iterator rend = iter->second.second;
1795
1796 // move all iterators pointing to the current global index to next value
1797 if(iter->second.first->localIndexPair().global()==index_ &&
1798 (noattribute || iter->second.first->localIndexPair().local().attribute() == attribute_))
1799 ++(iter->second.first);
1800
1801 // erase from the map if there are no more entries.
1802 if(iter->second.first == iter->second.second)
1803 map_.erase(iter++);
1804 else{
1805 ++iter;
1806 }
1807 }
1808 return *this;
1809 }
1810
1811 template<typename T, typename A>
1813 {
1814 return map_.empty();
1815 }
1816
1817 template<typename T, typename A>
1818 inline typename CollectiveIterator<T,A>::iterator
1820 {
1821 if(noattribute)
1822 return iterator(map_.begin(), map_.end(), index_);
1823 else
1824 return iterator(map_.begin(), map_.end(), index_,
1825 attribute_);
1826 }
1827
1828 template<typename T, typename A>
1829 inline typename CollectiveIterator<T,A>::iterator
1830 CollectiveIterator<T,A>::end()
1831 {
1832 return iterator(map_.end(), map_.end(), index_);
1833 }
1834
1835 template<typename TG, typename TA>
1836 inline std::ostream& operator<<(std::ostream& os, const RemoteIndex<TG,TA>& index)
1837 {
1838 os<<"[global="<<index.localIndexPair().global()<<", remote attribute="<<index.attribute()<<" local attribute="<<index.localIndexPair().local().attribute()<<"]";
1839 return os;
1840 }
1841
1842 template<typename T, typename A>
1843 inline std::ostream& operator<<(std::ostream& os, const RemoteIndices<T,A>& indices)
1844 {
1845 int rank;
1846 MPI_Comm_rank(indices.comm_, &rank);
1847
1848 typedef typename RemoteIndices<T,A>::RemoteIndexList RList;
1849 typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
1850
1851 const const_iterator rend = indices.remoteIndices_.end();
1852
1853 for(const_iterator rindex = indices.remoteIndices_.begin(); rindex!=rend; ++rindex) {
1854 os<<rank<<": Prozess "<<rindex->first<<":";
1855
1856 if(!rindex->second.first->empty()) {
1857 os<<" send:";
1858
1859 const typename RList::const_iterator send= rindex->second.first->end();
1860
1861 for(typename RList::const_iterator index = rindex->second.first->begin();
1862 index != send; ++index)
1863 os<<*index<<" ";
1864 os<<std::endl;
1865 }
1866 if(!rindex->second.second->empty()) {
1867 os<<rank<<": Prozess "<<rindex->first<<": "<<"receive: ";
1868
1869 const typename RList::const_iterator rend= rindex->second.second->end();
1870
1871 for(typename RList::const_iterator index = rindex->second.second->begin();
1872 index != rend; ++index)
1873 os << *index << " ";
1874 }
1875 os<<std::endl<<std::flush;
1876 }
1877 return os;
1878 }
1880}
1881
1882#endif
1883#endif
Iterator over the valid underlying iterators.
Definition: remoteindices.hh:787
bool operator==(const iterator &other)
Definition: remoteindices.hh:851
iterator(const RealIterator &iter, const ConstRealIterator &end, GlobalIndex &index)
Definition: remoteindices.hh:794
bool operator!=(const iterator &other)
Definition: remoteindices.hh:857
iterator(const iterator &other)
Definition: remoteindices.hh:812
const RemoteIndex & operator*() const
Definition: remoteindices.hh:833
iterator & operator++()
Definition: remoteindices.hh:817
const RemoteIndex * operator->() const
Definition: remoteindices.hh:845
int process() const
Definition: remoteindices.hh:839
A collective iterator for moving over the remote indices for all processes collectively.
Definition: remoteindices.hh:703
CollectiveIterator(const RemoteIndexMap &map_, bool send)
Constructor.
Definition: remoteindices.hh:1704
bool empty()
Checks whether there are still iterators in the map.
Definition: remoteindices.hh:1812
void advance(const GlobalIndex &global)
Advances all underlying iterators.
Definition: remoteindices.hh:1718
std::map< int, std::pair< RemoteIndexList *, RemoteIndexList * > > RemoteIndexMap
The type of the map from rank to remote index list.
Definition: remoteindices.hh:743
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:173
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:279
Modifier for adding and/or deleting remote indices from the remote index list.
Definition: remoteindices.hh:544
void repairLocalIndexPointers()
Repair the pointers to the local index pairs.
Definition: remoteindices.hh:1573
Dune::SLList< RemoteIndex, Allocator > RemoteIndexList
The type of the remote index list.
Definition: remoteindices.hh:597
A Allocator
The type of the allocator for the remote index list.
Definition: remoteindices.hh:593
void insert(const RemoteIndex &index)
Insert an index to the list.
Definition: remoteindices.hh:1606
ParallelIndexSet::GlobalIndex GlobalIndex
The type of the global index.
Definition: remoteindices.hh:573
ParallelIndexSet::LocalIndex LocalIndex
The type of the local index.
Definition: remoteindices.hh:578
@ MODIFYINDEXSET
If true the index set corresponding to the remote indices might get modified.
Definition: remoteindices.hh:562
RemoteIndexList::const_iterator ConstIterator
The type of the remote index list iterator.
Definition: remoteindices.hh:607
SLListModifyIterator< RemoteIndex, Allocator > ModifyIterator
The type of the modifying iterator of the remote index list.
Definition: remoteindices.hh:602
bool remove(const GlobalIndex &global)
Remove a remote index.
Definition: remoteindices.hh:1654
T ParallelIndexSet
Type of the index set we use.
Definition: remoteindices.hh:568
RemoteIndexListModifier()
Default constructor.
Definition: remoteindices.hh:671
LocalIndex::Attribute Attribute
The type of the attribute.
Definition: remoteindices.hh:583
Dune::RemoteIndex< GlobalIndex, Attribute > RemoteIndex
Type of the remote indices we manage.
Definition: remoteindices.hh:588
Information about an index residing on another processor.
Definition: remoteindices.hh:66
const Attribute attribute() const
Get the attribute of the index on the remote process.
Definition: remoteindices.hh:936
T1 GlobalIndex
the type of the global index. This type has to provide at least a operator< for sorting.
Definition: remoteindices.hh:83
T2 Attribute
The type of the attributes. Normally this will be an enumeration like.
Definition: remoteindices.hh:92
IndexPair< GlobalIndex, ParallelLocalIndex< Attribute > > PairType
The type of the index pair.
Definition: remoteindices.hh:98
const PairType & localIndexPair() const
Get the corresponding local index pair.
Definition: remoteindices.hh:942
RemoteIndex()
Parameterless Constructor.
Definition: remoteindices.hh:920
The indices present on remote processes.
Definition: remoteindices.hh:181
Dune::RemoteIndex< GlobalIndex, Attribute > RemoteIndex
Type of the remote indices we manage.
Definition: remoteindices.hh:223
friend void fillIndexSetHoles(const G &graph, Dune::OwnerOverlapCopyCommunication< T1, T2 > &oocomm)
Fills the holes in an index set.
Definition: repartition.hh:58
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:973
void free()
Free the index lists.
Definition: remoteindices.hh:1423
ParallelIndexSet::GlobalIndex GlobalIndex
The type of the global index.
Definition: remoteindices.hh:207
void rebuild()
Rebuilds the set of remote indices.
Definition: remoteindices.hh:1448
T ParallelIndexSet
Type of the index set we use, e.g. ParallelLocalIndexSet.
Definition: remoteindices.hh:198
MPI_Comm communicator() const
Get the mpi communicator used.
Definition: remoteindices.hh:1697
LocalIndex::Attribute Attribute
The type of the attribute.
Definition: remoteindices.hh:218
CollectiveIteratorT iterator() const
Get an iterator for colletively iterating over the remote indices of all remote processes.
Definition: remoteindices.hh:1691
void setIncludeSelf(bool includeSelf)
Tell whether sending from indices of the processor to other indices on the same processor is enabled ...
Definition: remoteindices.hh:961
const_iterator end() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1524
std::map< int, std::pair< RemoteIndexList *, RemoteIndexList * > > RemoteIndexMap
The type of the map from rank to remote index list.
Definition: remoteindices.hh:237
RemoteIndexListModifier< T, A, mode > getModifier(int process)
Get a modifier for a remote index list.
Definition: remoteindices.hh:1475
const ParallelIndexSet & sourceIndexSet() const
Get the index set at the source.
Definition: remoteindices.hh:988
~RemoteIndices()
Destructor.
Definition: remoteindices.hh:1003
Dune::SLList< RemoteIndex, Allocator > RemoteIndexList
The type of the remote index list.
Definition: remoteindices.hh:233
int neighbours() const
Get the number of processors we share indices with.
Definition: remoteindices.hh:1441
CollectiveIterator< T, A > CollectiveIteratorT
The type of the collective iterator over all remote indices.
Definition: remoteindices.hh:202
const ParallelIndexSet & destinationIndexSet() const
Get the index set at destination.
Definition: remoteindices.hh:996
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:948
const_iterator find(int proc) const
Find an iterator over the remote index lists of a specific process.
Definition: remoteindices.hh:1510
bool isSynced() const
Checks whether the remote indices are synced with the indexsets.
Definition: remoteindices.hh:1468
const_iterator begin() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1517
ParallelIndexSet::LocalIndex LocalIndex
The type of the local index.
Definition: remoteindices.hh:213
A::template rebind< RemoteIndex >::other Allocator
The type of the allocator for the remote index list.
Definition: remoteindices.hh:229
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:305
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
std::ostream & operator<<(std::ostream &s, const array< T, N > &e)
Output operator for array.
Definition: array.hh:26
iterator end()
Get an iterator pointing to the end of the list.
Definition: sllist.hh:788
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:776
#define DUNE_THROW(E, m)
Definition: exceptions.hh:243
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:230
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:252
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: alignment.hh:10
Provides classes for use as the local index in ParallelIndexSet for distributed computing.
An stl-compliant pool allocator.
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:37
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 20, 23:31, 2024)