3 #ifndef DUNE_REMOTEINDICES_HH
4 #define DUNE_REMOTEINDICES_HH
37 template<
typename TG,
typename TA>
39 class MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >
42 inline static MPI_Datatype getType();
44 static MPI_Datatype type;
48 template<
typename T,
typename A>
51 template<
typename T1,
typename T2>
58 template<
typename T1,
typename T2>
59 std::ostream&
operator<<(std::ostream& os,
const RemoteIndex<T1,T2>& index);
62 template<
typename T,
typename A,
bool mode>
63 class RemoteIndexListModifier;
69 template<
typename T1,
typename T2>
73 friend class IndicesSyncer;
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> >&,
80 template<
typename T,
typename A,
bool mode>
81 friend class RemoteIndexListModifier;
88 typedef T1 GlobalIndex;
102 typedef IndexPair<GlobalIndex,ParallelLocalIndex<Attribute> >
109 const Attribute attribute()
const;
116 const PairType& localIndexPair()
const;
129 RemoteIndex(
const T2& attribute,
130 const PairType* local);
138 RemoteIndex(
const T2& attribute);
145 const PairType* localIndex_;
151 template<
class T,
class A>
152 std::ostream&
operator<<(std::ostream& os,
const RemoteIndices<T,A>& indices);
154 class InterfaceBuilder;
156 template<
class T,
class A>
157 class CollectiveIterator;
164 template<
typename T1,
typename T2>
165 class OwnerOverlapCopyCommunication;
184 template<
class T,
class A=std::allocator<RemoteIndex<
typename T::GlobalIndex,
185 typename T::LocalIndex::Attribute> > >
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> >&,
192 RemoteIndices<T1,A1>&,
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>&);
204 typedef T ParallelIndexSet;
208 typedef CollectiveIterator<T,A> CollectiveIteratorT;
224 typedef typename LocalIndex::Attribute Attribute;
229 typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
235 typedef typename A::template rebind<RemoteIndex>::other Allocator;
242 typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
245 typedef typename RemoteIndexMap::const_iterator const_iterator;
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);
276 void setIncludeSelf(
bool includeSelf);
294 void setIndexSets(
const ParallelIndexSet& source,
const ParallelIndexSet& destination,
295 const MPI_Comm& comm,
const std::vector<int>& neighbours=std::vector<int>());
298 void setNeighbours(
const C& neighbours)
300 neighbourIds.clear();
301 neighbourIds.insert(neighbours.begin(), neighbours.end());
305 const std::set<int>& getNeighbours()
const
324 template<
bool ignorePublic>
336 inline bool isSynced()
const;
341 inline MPI_Comm communicator()
const;
357 template<
bool mode,
bool send>
358 inline RemoteIndexListModifier<T,A,mode> getModifier(
int process);
366 inline const_iterator find(
int proc)
const;
372 inline const_iterator begin()
const;
378 inline const_iterator end()
const;
384 inline CollectiveIteratorT iterator()
const;
395 inline int neighbours()
const;
398 inline const ParallelIndexSet& sourceIndexSet()
const;
401 inline const ParallelIndexSet& destinationIndexSet()
const;
405 RemoteIndices(
const RemoteIndices&)
409 const ParallelIndexSet* source_;
412 const ParallelIndexSet* target_;
419 std::set<int> neighbourIds;
422 const static int commTag_=333;
455 typedef IndexPair<GlobalIndex, LocalIndex>
464 RemoteIndexMap remoteIndices_;
476 template<
bool ignorePublic>
477 inline void buildRemote(
bool includeSelf);
484 inline int noPublic(
const ParallelIndexSet& indexSet);
497 template<
bool ignorePublic>
498 inline void packEntries(PairType** myPairs,
const ParallelIndexSet& indexSet,
499 char* p_out, MPI_Datatype type,
int bufferSize,
500 int* position,
int n);
515 inline void unpackIndices(RemoteIndexList& remote,
int remoteEntries,
516 PairType** local,
int localEntries,
char* p_in,
517 MPI_Datatype type,
int* position,
int bufferSize,
520 inline void unpackIndices(RemoteIndexList& send, RemoteIndexList& receive,
521 int remoteEntries, PairType** localSource,
522 int localSourceEntries, PairType** localDest,
523 int localDestEntries,
char* p_in,
524 MPI_Datatype type,
int* position,
int bufferSize);
526 void unpackCreateRemote(
char* p_in, PairType** sourcePairs, PairType** DestPairs,
527 int remoteProc,
int sourcePublish,
int destPublish,
528 int bufferSize,
bool sendTwo,
bool fromOurSelf=
false);
548 template<
class T,
class A,
bool mode>
549 class RemoteIndexListModifier
552 template<
typename T1,
typename A1>
553 friend class RemoteIndices;
556 class InvalidPosition :
public RangeError
574 typedef T ParallelIndexSet;
589 typedef typename LocalIndex::Attribute Attribute;
594 typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
608 typedef SLListModifyIterator<RemoteIndex,Allocator> ModifyIterator;
613 typedef typename RemoteIndexList::const_iterator ConstIterator;
628 void insert(
const RemoteIndex& index);
645 void insert(
const RemoteIndex& index,
const GlobalIndex& global);
654 bool remove(
const GlobalIndex& global);
668 void repairLocalIndexPointers();
671 RemoteIndexListModifier(
const RemoteIndexListModifier&);
677 RemoteIndexListModifier()
688 RemoteIndexListModifier(
const ParallelIndexSet& indexSet,
689 RemoteIndexList& rList);
691 typedef SLList<GlobalIndex,Allocator> GlobalList;
692 typedef typename GlobalList::ModifyIterator GlobalModifyIterator;
693 RemoteIndexList* rList_;
694 const ParallelIndexSet* indexSet_;
696 ModifyIterator iter_;
697 GlobalModifyIterator giter_;
707 template<
class T,
class A>
708 class CollectiveIterator
714 typedef T ParallelIndexSet;
729 typedef typename LocalIndex::Attribute Attribute;
732 typedef Dune::RemoteIndex<GlobalIndex,Attribute> RemoteIndex;
735 typedef typename A::template rebind<RemoteIndex>::other Allocator;
741 typedef std::map<int,std::pair<
typename RemoteIndexList::const_iterator,
742 const typename RemoteIndexList::const_iterator> >
748 typedef std::map<int, std::pair<RemoteIndexList*,RemoteIndexList*> >
756 inline CollectiveIterator(
const RemoteIndexMap& map_,
bool send);
766 inline void advance(
const GlobalIndex& global);
777 inline void advance(
const GlobalIndex& global,
const Attribute& attribute);
779 CollectiveIterator& operator++();
795 typedef typename Map::iterator RealIterator;
796 typedef typename Map::iterator ConstRealIterator;
800 iterator(
const RealIterator& iter,
const ConstRealIterator& end, GlobalIndex& index)
801 : iter_(iter), end_(end), index_(index), hasAttribute(false)
804 while(iter_!=end_ && iter_->second.first->localIndexPair().global()!=index_)
808 iterator(
const RealIterator& iter,
const ConstRealIterator& end, GlobalIndex index,
810 : iter_(iter), end_(end), index_(index), attribute_(attribute), hasAttribute(true)
813 while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_
814 || iter_->second.first->localIndexPair().local().attribute()!=attribute))
818 iterator(
const iterator& other)
819 : iter_(other.iter_), end_(other.end_), index_(other.index_)
823 iterator& operator++()
827 while(iter_!=end_ && (iter_->second.first->localIndexPair().global()!=index_ ||
829 iter_->second.first->localIndexPair().local().attribute()!=attribute_)))
831 assert(iter_==end_ ||
832 (iter_->second.first->localIndexPair().global()==index_));
833 assert(iter_==end_ || !hasAttribute ||
834 (iter_->second.first->localIndexPair().local().attribute()==attribute_));
841 return *(iter_->second.first);
851 const RemoteIndex* operator->()
const
853 return iter_->second.first.operator->();
859 return other.iter_==iter_;
865 return other.iter_!=iter_;
874 Attribute attribute_;
886 Attribute attribute_;
890 template<
typename TG,
typename TA>
891 MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::getType()
893 if(type==MPI_DATATYPE_NULL) {
894 int length[2] = {1, 1};
897 MPI_Datatype types[2] = {MPITraits<TG>::getType(),
898 MPITraits<ParallelLocalIndex<TA> >::getType()};
899 IndexPair<TG,ParallelLocalIndex<TA> > rep;
900 MPI_Get_address(&rep, &base);
901 MPI_Get_address(&(rep.global_), &disp[0]);
902 MPI_Get_address(&(rep.local_), &disp[1]);
903 for (MPI_Aint& d : disp)
907 MPI_Type_create_struct(2, length, disp, types, &tmp);
909 MPI_Type_create_resized(tmp, 0,
sizeof(IndexPair<TG,ParallelLocalIndex<TA> >), &type);
910 MPI_Type_commit(&type);
917 template<
typename TG,
typename TA>
918 MPI_Datatype MPITraits<IndexPair<TG,ParallelLocalIndex<TA> > >::type=MPI_DATATYPE_NULL;
920 template<
typename T1,
typename T2>
921 RemoteIndex<T1,T2>::RemoteIndex(
const T2& attribute,
const PairType* local)
922 : localIndex_(local), attribute_(attribute)
925 template<
typename T1,
typename T2>
926 RemoteIndex<T1,T2>::RemoteIndex(
const T2& attribute)
927 : localIndex_(0), attribute_(attribute)
930 template<
typename T1,
typename T2>
931 RemoteIndex<T1,T2>::RemoteIndex()
932 : localIndex_(0), attribute_()
934 template<
typename T1,
typename T2>
937 return localIndex_==ri.localIndex_ && attribute_==ri.attribute;
940 template<
typename T1,
typename T2>
943 return localIndex_!=ri.localIndex_ || attribute_!=ri.attribute_;
946 template<
typename T1,
typename T2>
947 inline const T2 RemoteIndex<T1,T2>::attribute()
const
949 return T2(attribute_);
952 template<
typename T1,
typename T2>
953 inline const IndexPair<T1,ParallelLocalIndex<T2> >& RemoteIndex<T1,T2>::localIndexPair()
const
958 template<
typename T,
typename A>
959 inline RemoteIndices<T,A>::RemoteIndices(
const ParallelIndexSet& source,
960 const ParallelIndexSet& destination,
961 const MPI_Comm& comm,
962 const std::vector<int>& neighbours,
964 : source_(&source), target_(&destination), comm_(comm),
965 sourceSeqNo_(-1), destSeqNo_(-1), publicIgnored(false), firstBuild(true),
966 includeSelf(includeSelf_)
968 setNeighbours(neighbours);
971 template<
typename T,
typename A>
972 void RemoteIndices<T,A>::setIncludeSelf(
bool b)
977 template<
typename T,
typename A>
978 RemoteIndices<T,A>::RemoteIndices()
979 : source_(0), target_(0), sourceSeqNo_(-1),
980 destSeqNo_(-1), publicIgnored(false), firstBuild(true),
984 template<
class T,
typename A>
985 void RemoteIndices<T,A>::setIndexSets(
const ParallelIndexSet& source,
986 const ParallelIndexSet& destination,
987 const MPI_Comm& comm,
988 const std::vector<int>& neighbours)
992 target_ = &destination;
995 setNeighbours(neighbours);
998 template<
typename T,
typename A>
999 const typename RemoteIndices<T,A>::ParallelIndexSet&
1000 RemoteIndices<T,A>::sourceIndexSet()
const
1006 template<
typename T,
typename A>
1007 const typename RemoteIndices<T,A>::ParallelIndexSet&
1008 RemoteIndices<T,A>::destinationIndexSet()
const
1014 template<
typename T,
typename A>
1015 RemoteIndices<T,A>::~RemoteIndices()
1020 template<
typename T,
typename A>
1021 template<
bool ignorePublic>
1022 inline void RemoteIndices<T,A>::packEntries(IndexPair<GlobalIndex,LocalIndex>** pairs,
1023 const ParallelIndexSet& indexSet,
1024 char* p_out, MPI_Datatype type,
1026 int *position,
int n)
1030 typedef typename ParallelIndexSet::const_iterator const_iterator;
1031 typedef IndexPair<GlobalIndex,LocalIndex> PairType;
1032 const const_iterator end = indexSet.end();
1036 for(const_iterator index = indexSet.begin(); index != end; ++index)
1037 if(ignorePublic || index->local().isPublic()) {
1039 MPI_Pack(
const_cast<PairType*
>(&(*index)), 1,
1041 p_out, bufferSize, position, comm_);
1042 pairs[i++] =
const_cast<PairType*
>(&(*index));
1048 template<
typename T,
typename A>
1049 inline int RemoteIndices<T,A>::noPublic(
const ParallelIndexSet& indexSet)
1051 typedef typename ParallelIndexSet::const_iterator const_iterator;
1055 const const_iterator end=indexSet.end();
1056 for(const_iterator index=indexSet.begin(); index!=end; ++index)
1057 if(index->local().isPublic())
1065 template<
typename T,
typename A>
1066 inline void RemoteIndices<T,A>::unpackCreateRemote(
char* p_in, PairType** sourcePairs,
1067 PairType** destPairs,
int remoteProc,
1068 int sourcePublish,
int destPublish,
1069 int bufferSize,
bool sendTwo,
1074 int noRemoteSource=-1, noRemoteDest=-1;
1075 char twoIndexSets=0;
1078 MPI_Unpack(p_in, bufferSize, &position, &twoIndexSets, 1, MPI_CHAR, comm_);
1080 MPI_Unpack(p_in, bufferSize, &position, &noRemoteSource, 1, MPI_INT, comm_);
1082 MPI_Unpack(p_in, bufferSize, &position, &noRemoteDest, 1, MPI_INT, comm_);
1086 RemoteIndexList* receive=
new RemoteIndexList();
1088 RemoteIndexList* send=0;
1090 MPI_Datatype type= MPITraits<PairType>::getType();
1094 send =
new RemoteIndexList();
1096 unpackIndices(*send, *receive, noRemoteSource, sourcePairs, sourcePublish,
1097 destPairs, destPublish, p_in, type, &position, bufferSize);
1100 unpackIndices(*receive, noRemoteSource, sourcePairs, sourcePublish,
1101 p_in, type, &position, bufferSize, fromOurSelf);
1106 int oldPos=position;
1108 unpackIndices(*receive, noRemoteSource, destPairs, destPublish,
1109 p_in, type, &position, bufferSize, fromOurSelf);
1114 send =
new RemoteIndexList();
1115 unpackIndices(*send, noRemoteDest, sourcePairs, sourcePublish,
1116 p_in, type, &position, bufferSize, fromOurSelf);
1119 if(receive->empty() && send->empty()) {
1127 remoteIndices_.insert(std::make_pair(remoteProc,
1128 std::make_pair(send,receive)));
1133 template<
typename T,
typename A>
1134 template<
bool ignorePublic>
1135 inline void RemoteIndices<T,A>::buildRemote(
bool includeSelf_)
1139 MPI_Comm_rank(comm_, &rank);
1140 MPI_Comm_size(comm_, &procs);
1144 int sourcePublish, destPublish;
1147 char sendTwo = (source_ != target_);
1149 if(procs==1 && !(sendTwo || includeSelf_))
1153 sourcePublish = (ignorePublic) ? source_->size() : noPublic(*source_);
1156 destPublish = (ignorePublic) ? target_->size() : noPublic(*target_);
1161 int maxPublish, publish=sourcePublish+destPublish;
1164 MPI_Allreduce(&publish, &maxPublish, 1, MPI_INT, MPI_MAX, comm_);
1167 typedef IndexPair<GlobalIndex,LocalIndex> PairType;
1169 PairType** destPairs;
1170 PairType** sourcePairs =
new PairType*[sourcePublish>0 ? sourcePublish : 1];
1173 destPairs =
new PairType*[destPublish>0 ? destPublish : 1];
1175 destPairs=sourcePairs;
1177 char** buffer =
new char*[2];
1184 MPI_Datatype type = MPITraits<PairType>::getType();
1186 MPI_Pack_size(maxPublish, type, comm_,
1188 MPI_Pack_size(1, MPI_INT, comm_,
1190 MPI_Pack_size(1, MPI_CHAR, comm_,
1196 bufferSize += 2 * intSize + charSize;
1198 if(bufferSize<=0) bufferSize=1;
1200 buffer[0] =
new char[bufferSize];
1201 buffer[1] =
new char[bufferSize];
1205 MPI_Pack(&sendTwo, 1, MPI_CHAR, buffer[0], bufferSize, &position,
1209 MPI_Pack(&sourcePublish, 1, MPI_INT, buffer[0], bufferSize, &position,
1211 MPI_Pack(&destPublish, 1, MPI_INT, buffer[0], bufferSize, &position,
1215 packEntries<ignorePublic>(sourcePairs, *source_, buffer[0], type,
1216 bufferSize, &position, sourcePublish);
1219 packEntries<ignorePublic>(destPairs, *target_, buffer[0], type,
1220 bufferSize, &position, destPublish);
1224 if(sendTwo|| includeSelf_)
1225 unpackCreateRemote(buffer[0], sourcePairs, destPairs, rank, sourcePublish,
1226 destPublish, bufferSize, sendTwo, includeSelf_);
1228 neighbourIds.erase(rank);
1230 if(neighbourIds.size()==0)
1232 Dune::dvverb<<rank<<
": Sending messages in a ring"<<std::endl;
1234 for(
int proc=1; proc<procs; proc++) {
1236 char* p_out = buffer[1-(proc%2)];
1237 char* p_in = buffer[proc%2];
1241 MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
1243 MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
1244 commTag_, comm_, &status);
1246 MPI_Recv(p_in, bufferSize, MPI_PACKED, (rank+procs-1)%procs,
1247 commTag_, comm_, &status);
1248 MPI_Ssend(p_out, bufferSize, MPI_PACKED, (rank+1)%procs,
1254 int remoteProc = (rank+procs-proc)%procs;
1256 unpackCreateRemote(p_in, sourcePairs, destPairs, remoteProc, sourcePublish,
1257 destPublish, bufferSize, sendTwo);
1264 MPI_Request* requests=
new MPI_Request[neighbourIds.size()];
1265 MPI_Request* req=requests;
1267 typedef typename std::set<int>::size_type size_type;
1268 size_type noNeighbours=neighbourIds.size();
1271 for(std::set<int>::iterator neighbour=neighbourIds.begin();
1272 neighbour!= neighbourIds.end(); ++neighbour) {
1274 MPI_Issend(buffer[0], position , MPI_PACKED, *neighbour, commTag_, comm_, req++);
1279 for(size_type received=0; received <noNeighbours; ++received)
1283 MPI_Probe(MPI_ANY_SOURCE, commTag_, comm_, &status);
1284 int remoteProc=status.MPI_SOURCE;
1286 MPI_Get_count(&status, MPI_PACKED, &size);
1288 MPI_Recv(buffer[1], size, MPI_PACKED, remoteProc,
1289 commTag_, comm_, &status);
1291 unpackCreateRemote(buffer[1], sourcePairs, destPairs, remoteProc, sourcePublish,
1292 destPublish, bufferSize, sendTwo);
1295 MPI_Status* statuses =
new MPI_Status[neighbourIds.size()];
1297 if(MPI_ERR_IN_STATUS==MPI_Waitall(neighbourIds.size(), requests, statuses)) {
1298 for(size_type i=0; i < neighbourIds.size(); ++i)
1299 if(statuses[i].MPI_ERROR!=MPI_SUCCESS) {
1300 std::cerr<<rank<<
": MPI_Error occurred while receiving message."<<std::endl;
1301 MPI_Abort(comm_, 999);
1310 if(destPairs!=sourcePairs)
1313 delete[] sourcePairs;
1319 template<
typename T,
typename A>
1320 inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& remote,
1330 if(remoteEntries==0)
1334 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1336 GlobalIndex oldGlobal=index.global();
1337 int n_in=0, localIndex=0;
1340 while(localIndex<localEntries) {
1341 if(local[localIndex]->global()==index.global()) {
1342 int oldLocalIndex=localIndex;
1344 while(localIndex<localEntries &&
1345 local[localIndex]->global()==index.global()) {
1346 if(!fromOurSelf || index.local().attribute() !=
1347 local[localIndex]->local().attribute())
1349 remote.push_back(RemoteIndex(index.local().attribute(),
1350 local[localIndex]));
1355 if((++n_in) < remoteEntries) {
1356 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1358 if(index.global()==oldGlobal)
1360 localIndex=oldLocalIndex;
1362 oldGlobal=index.global();
1370 if (local[localIndex]->global()<index.global()) {
1375 if((++n_in) < remoteEntries) {
1376 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1378 oldGlobal=index.global();
1386 while(++n_in < remoteEntries)
1387 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1392 template<
typename T,
typename A>
1393 inline void RemoteIndices<T,A>::unpackIndices(RemoteIndexList& send,
1394 RemoteIndexList& receive,
1396 PairType** localSource,
1397 int localSourceEntries,
1398 PairType** localDest,
1399 int localDestEntries,
1405 int n_in=0, sourceIndex=0, destIndex=0;
1408 while(n_in<remoteEntries && (sourceIndex<localSourceEntries || destIndex<localDestEntries)) {
1411 MPI_Unpack(p_in, bufferSize, position, &index, 1,
1416 while(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()<index.global())
1419 while(destIndex<localDestEntries && localDest[destIndex]->global()<index.global())
1423 if(sourceIndex<localSourceEntries && localSource[sourceIndex]->global()==index.global())
1424 send.push_back(RemoteIndex(index.local().attribute(),
1425 localSource[sourceIndex]));
1427 if(destIndex < localDestEntries && localDest[destIndex]->global() == index.global())
1428 receive.push_back(RemoteIndex(index.local().attribute(),
1429 localDest[sourceIndex]));
1434 template<
typename T,
typename A>
1435 inline void RemoteIndices<T,A>::free()
1437 typedef typename RemoteIndexMap::iterator Iterator;
1438 Iterator lend = remoteIndices_.end();
1439 for(Iterator lists=remoteIndices_.begin(); lists != lend; ++lists) {
1440 if(lists->second.first==lists->second.second) {
1442 delete lists->second.first;
1444 delete lists->second.first;
1445 delete lists->second.second;
1448 remoteIndices_.clear();
1452 template<
typename T,
typename A>
1453 inline int RemoteIndices<T,A>::neighbours()
const
1455 return remoteIndices_.size();
1458 template<
typename T,
typename A>
1459 template<
bool ignorePublic>
1460 inline void RemoteIndices<T,A>::rebuild()
1464 ignorePublic!=publicIgnored || !
1468 buildRemote<ignorePublic>(includeSelf);
1470 sourceSeqNo_ = source_->seqNo();
1471 destSeqNo_ = target_->seqNo();
1473 publicIgnored=ignorePublic;
1479 template<
typename T,
typename A>
1480 inline bool RemoteIndices<T,A>::isSynced()
const
1482 return sourceSeqNo_==source_->seqNo() && destSeqNo_ ==target_->seqNo();
1485 template<
typename T,
typename A>
1486 template<
bool mode,
bool send>
1487 RemoteIndexListModifier<T,A,mode> RemoteIndices<T,A>::getModifier(
int process)
1493 sourceSeqNo_ = source_->seqNo();
1494 destSeqNo_ = target_->seqNo();
1496 typename RemoteIndexMap::iterator found = remoteIndices_.find(process);
1498 if(found == remoteIndices_.end())
1500 if(source_ != target_)
1501 found = remoteIndices_.insert(found, std::make_pair(process,
1502 std::make_pair(
new RemoteIndexList(),
1503 new RemoteIndexList())));
1505 RemoteIndexList* rlist =
new RemoteIndexList();
1506 found = remoteIndices_.insert(found,
1507 std::make_pair(process,
1508 std::make_pair(rlist, rlist)));
1515 return RemoteIndexListModifier<T,A,mode>(*source_, *(found->second.first));
1517 return RemoteIndexListModifier<T,A,mode>(*target_, *(found->second.second));
1520 template<
typename T,
typename A>
1521 inline typename RemoteIndices<T,A>::const_iterator
1522 RemoteIndices<T,A>::find(
int proc)
const
1524 return remoteIndices_.find(proc);
1527 template<
typename T,
typename A>
1528 inline typename RemoteIndices<T,A>::const_iterator
1529 RemoteIndices<T,A>::begin()
const
1531 return remoteIndices_.begin();
1534 template<
typename T,
typename A>
1535 inline typename RemoteIndices<T,A>::const_iterator
1536 RemoteIndices<T,A>::end()
const
1538 return remoteIndices_.end();
1542 template<
typename T,
typename A>
1545 if(neighbours()!=ri.neighbours())
1548 typedef RemoteIndexList RList;
1549 typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
1551 const const_iterator rend = remoteIndices_.end();
1553 for(const_iterator rindex = remoteIndices_.begin(), rindex1=ri.remoteIndices_.begin(); rindex!=rend; ++rindex, ++rindex1) {
1554 if(rindex->first != rindex1->first)
1556 if(*(rindex->second.first) != *(rindex1->second.first))
1558 if(*(rindex->second.second) != *(rindex1->second.second))
1564 template<
class T,
class A,
bool mode>
1565 RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(
const ParallelIndexSet& indexSet,
1566 RemoteIndexList& rList)
1567 : rList_(&rList), indexSet_(&indexSet), iter_(rList.beginModify()), end_(rList.end()), first_(true)
1569 if(MODIFYINDEXSET) {
1571 for(ConstIterator iter=iter_; iter != end_; ++iter)
1572 glist_.push_back(iter->localIndexPair().global());
1573 giter_ = glist_.beginModify();
1577 template<
typename T,
typename A,
bool mode>
1578 RemoteIndexListModifier<T,A,mode>::RemoteIndexListModifier(
const RemoteIndexListModifier<T,A,mode>& other)
1579 : rList_(other.rList_), indexSet_(other.indexSet_),
1580 glist_(other.glist_), iter_(other.iter_), giter_(other.giter_), end_(other.end_),
1581 first_(other.first_), last_(other.last_)
1584 template<
typename T,
typename A,
bool mode>
1585 inline void RemoteIndexListModifier<T,A,mode>::repairLocalIndexPointers()
1587 if(MODIFYINDEXSET) {
1589 #ifdef DUNE_ISTL_WITH_CHECKING
1590 if(indexSet_->state()!=
GROUND)
1591 DUNE_THROW(InvalidIndexSetState,
"Index has to be in ground mode for repairing pointers to indices");
1593 typedef typename ParallelIndexSet::const_iterator IndexIterator;
1594 typedef typename GlobalList::const_iterator GlobalIterator;
1595 typedef typename RemoteIndexList::iterator Iterator;
1596 GlobalIterator giter = glist_.begin();
1597 IndexIterator index = indexSet_->begin();
1599 for(Iterator iter=rList_->begin(); iter != end_; ++iter) {
1600 while(index->global()<*giter) {
1602 #ifdef DUNE_ISTL_WITH_CHECKING
1603 if(index == indexSet_->end())
1604 DUNE_THROW(InvalidPosition,
"No such global index in set!");
1608 #ifdef DUNE_ISTL_WITH_CHECKING
1609 if(index->global() != *giter)
1610 DUNE_THROW(InvalidPosition,
"No such global index in set!");
1612 iter->localIndex_ = &(*index);
1617 template<
typename T,
typename A,
bool mode>
1618 inline void RemoteIndexListModifier<T,A,mode>::insert(
const RemoteIndex& index)
1620 static_assert(!mode,
"Not allowed if the mode indicates that new indices"
1621 "might be added to the underlying index set. Use "
1622 "insert(const RemoteIndex&, const GlobalIndex&) instead");
1624 #ifdef DUNE_ISTL_WITH_CHECKING
1625 if(!first_ && index.localIndexPair().global()<last_)
1626 DUNE_THROW(InvalidPosition,
"Modifcation of remote indices have to occur with ascending global index.");
1629 while(iter_ != end_ && iter_->localIndexPair().global() < index.localIndexPair().global()) {
1634 assert(iter_==end_ || iter_->localIndexPair().global() != index.localIndexPair().global());
1635 iter_.insert(index);
1636 last_ = index.localIndexPair().global();
1640 template<
typename T,
typename A,
bool mode>
1641 inline void RemoteIndexListModifier<T,A,mode>::insert(
const RemoteIndex& index,
const GlobalIndex& global)
1643 static_assert(mode,
"Not allowed if the mode indicates that no new indices"
1644 "might be added to the underlying index set. Use "
1645 "insert(const RemoteIndex&) instead");
1646 #ifdef DUNE_ISTL_WITH_CHECKING
1647 if(!first_ && global<last_)
1648 DUNE_THROW(InvalidPosition,
"Modification of remote indices have to occur with ascending global index.");
1651 while(iter_ != end_ && *giter_ < global) {
1657 assert(iter_->localIndexPair().global() != global);
1658 iter_.insert(index);
1659 giter_.insert(global);
1665 template<
typename T,
typename A,
bool mode>
1666 bool RemoteIndexListModifier<T,A,mode>::remove(
const GlobalIndex& global)
1668 #ifdef DUNE_ISTL_WITH_CHECKING
1669 if(!first_ && global<last_)
1670 DUNE_THROW(InvalidPosition,
"Modifcation of remote indices have to occur with ascending global index.");
1675 if(MODIFYINDEXSET) {
1677 while(iter_!=end_ && *giter_< global) {
1681 if(*giter_ == global) {
1687 while(iter_!=end_ && iter_->localIndexPair().global() < global)
1690 if(iter_->localIndexPair().global()==global) {
1701 template<
typename T,
typename A>
1703 inline typename RemoteIndices<T,A>::CollectiveIteratorT RemoteIndices<T,A>::iterator()
const
1705 return CollectiveIterator<T,A>(remoteIndices_, send);
1708 template<
typename T,
typename A>
1709 inline MPI_Comm RemoteIndices<T,A>::communicator()
const
1715 template<
typename T,
typename A>
1716 CollectiveIterator<T,A>::CollectiveIterator(
const RemoteIndexMap& pmap,
bool send)
1718 typedef typename RemoteIndexMap::const_iterator const_iterator;
1720 const const_iterator end=pmap.end();
1721 for(const_iterator process=pmap.begin(); process != end; ++process) {
1722 const RemoteIndexList* list = send ? process->second.first : process->second.second;
1723 typedef typename RemoteIndexList::const_iterator iterator;
1724 map_.insert(std::make_pair(process->first,
1725 std::pair<iterator, const iterator>(list->begin(), list->end())));
1729 template<
typename T,
typename A>
1730 inline void CollectiveIterator<T,A>::advance(
const GlobalIndex& index)
1732 typedef typename Map::iterator iterator;
1733 typedef typename Map::const_iterator const_iterator;
1734 const const_iterator end = map_.end();
1736 for(iterator iter = map_.begin(); iter != end;) {
1738 typename RemoteIndexList::const_iterator current = iter->second.first;
1739 typename RemoteIndexList::const_iterator rend = iter->second.second;
1740 RemoteIndex remoteIndex;
1742 remoteIndex = *current;
1744 while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
1745 ++(iter->second.first);
1748 if(iter->second.first == iter->second.second)
1758 template<
typename T,
typename A>
1759 inline void CollectiveIterator<T,A>::advance(
const GlobalIndex& index,
1760 const Attribute& attribute)
1762 typedef typename Map::iterator iterator;
1763 typedef typename Map::const_iterator const_iterator;
1764 const const_iterator end = map_.end();
1766 for(iterator iter = map_.begin(); iter != end;) {
1768 typename RemoteIndexList::const_iterator current = iter->second.first;
1769 typename RemoteIndexList::const_iterator rend = iter->second.second;
1770 RemoteIndex remoteIndex;
1772 remoteIndex = *current;
1775 while(iter->second.first!=iter->second.second && iter->second.first->localIndexPair().global()<index)
1776 ++(iter->second.first);
1779 while(iter->second.first!=iter->second.second
1780 && iter->second.first->localIndexPair().global()==index
1781 && iter->second.first->localIndexPair().local().attribute()<attribute)
1782 ++(iter->second.first);
1785 if(iter->second.first == iter->second.second)
1792 attribute_=attribute;
1796 template<
typename T,
typename A>
1797 inline CollectiveIterator<T,A>& CollectiveIterator<T,A>::operator++()
1799 typedef typename Map::iterator iterator;
1800 typedef typename Map::const_iterator const_iterator;
1801 const const_iterator end = map_.end();
1803 for(iterator iter = map_.begin(); iter != end;) {
1805 typename RemoteIndexList::const_iterator current = iter->second.first;
1806 typename RemoteIndexList::const_iterator rend = iter->second.second;
1809 if(iter->second.first->localIndexPair().global()==index_ &&
1810 (noattribute || iter->second.first->localIndexPair().local().attribute() == attribute_))
1811 ++(iter->second.first);
1814 if(iter->second.first == iter->second.second)
1823 template<
typename T,
typename A>
1824 inline bool CollectiveIterator<T,A>::empty()
1826 return map_.empty();
1829 template<
typename T,
typename A>
1830 inline typename CollectiveIterator<T,A>::iterator
1831 CollectiveIterator<T,A>::begin()
1834 return iterator(map_.begin(), map_.end(), index_);
1836 return iterator(map_.begin(), map_.end(), index_,
1840 template<
typename T,
typename A>
1841 inline typename CollectiveIterator<T,A>::iterator
1842 CollectiveIterator<T,A>::end()
1844 return iterator(map_.end(), map_.end(), index_);
1847 template<
typename TG,
typename TA>
1848 inline std::ostream&
operator<<(std::ostream& os,
const RemoteIndex<TG,TA>& index)
1850 os<<
"[global="<<index.localIndexPair().global()<<
", remote attribute="<<index.attribute()<<
" local attribute="<<index.localIndexPair().local().attribute()<<
"]";
1854 template<
typename T,
typename A>
1855 inline std::ostream&
operator<<(std::ostream& os,
const RemoteIndices<T,A>& indices)
1858 MPI_Comm_rank(indices.comm_, &rank);
1860 typedef typename RemoteIndices<T,A>::RemoteIndexList RList;
1861 typedef typename std::map<int,std::pair<RList*,RList*> >::const_iterator const_iterator;
1863 const const_iterator rend = indices.remoteIndices_.end();
1865 for(const_iterator rindex = indices.remoteIndices_.begin(); rindex!=rend; ++rindex) {
1866 os<<rank<<
": Prozess "<<rindex->first<<
":";
1868 if(!rindex->second.first->empty()) {
1871 const typename RList::const_iterator send= rindex->second.first->end();
1873 for(
typename RList::const_iterator index = rindex->second.first->begin();
1874 index != send; ++index)
1878 if(!rindex->second.second->empty()) {
1879 os<<rank<<
": Prozess "<<rindex->first<<
": "<<
"receive: ";
1881 for(
const auto& index : *(rindex->second.second))
1884 os<<std::endl<<std::flush;