4#ifndef DUNE_INDICESSYNCER_HH
5#define DUNE_INDICESSYNCER_HH
54 typedef typename ParallelIndexSet::LocalIndex::Attribute
Attribute;
95 void sync(T1& numberer);
109 char* receiveBuffer_;
112 std::size_t* sendBufferSizes_;
115 int receiveBufferSize_;
120 struct MessageInformation
137 class DefaultNumberer
148 return std::numeric_limits<size_t>::max();
153 MPI_Datatype datatype_;
174 typedef std::map<int, GlobalIndexList> GlobalIndicesMap;
184 GlobalIndicesMap globalMap_;
200 typedef std::map<int,BoolList> BoolMap;
209 std::map<int,MessageInformation> infoSend_;
227 typedef Dune::tuple<RemoteIndexModifier,GlobalIndexModifier,BoolListModifier,
228 const ConstRemoteIndexIterator> IteratorTuple;
250 Iterators(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
261 Iterators& operator++();
269 const std::pair<GlobalIndex,Attribute>& global);
281 std::pair<GlobalIndex,Attribute>& globalIndexPair()
const;
301 void reset(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
309 bool isNotAtEnd()
const;
316 bool isAtEnd()
const;
328 IteratorTuple iterators_;
332 typedef std::map<int,Iterators> IteratorsMap;
345 IteratorsMap iteratorsMap_;
348 void calculateMessageSizes();
357 void packAndSend(
int destination,
char* buffer, std::size_t bufferSize, MPI_Request& req);
363 template<
typename T1>
364 void recvAndUnpack(T1& numberer);
369 void registerMessageDatatype();
374 void insertIntoRemoteIndexList(
int process,
375 const std::pair<GlobalIndex,Attribute>& global,
381 void resetIteratorsMap();
397 bool checkReset(
const Iterators& iterators, RemoteIndexList& rlist, GlobalIndexList& gList,
401 template<
typename TG,
typename TA>
402 bool operator<(
const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
403 const std::pair<TG,TA>& i2)
405 return i1.global() < i2.first ||
406 (i1.global() == i2.first && i1.local().attribute()<i2.second);
409 template<
typename TG,
typename TA>
410 bool operator<(
const std::pair<TG,TA>& i1,
411 const IndexPair<TG,ParallelLocalIndex<TA> >& i2)
413 return i1.first < i2.global() ||
414 (i1.first == i2.global() && i1.second<i2.local().attribute());
417 template<
typename TG,
typename TA>
418 bool operator==(
const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
419 const std::pair<TG,TA>& i2)
421 return (i1.global() == i2.first && i1.local().attribute()==i2.second);
424 template<
typename TG,
typename TA>
425 bool operator!=(
const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
426 const std::pair<TG,TA>& i2)
428 return (i1.global() != i2.first || i1.local().attribute()!=i2.second);
431 template<
typename TG,
typename TA>
433 const IndexPair<TG,ParallelLocalIndex<TA> >& i1)
435 return (i1.global() == i2.first && i1.local().attribute()==i2.second);
438 template<
typename TG,
typename TA>
440 const IndexPair<TG,ParallelLocalIndex<TA> >& i1)
442 return (i1.global() != i2.first || i1.local().attribute()!=i2.second);
461 template<
typename T,
typename A,
typename A1>
466 typedef typename RemoteIndices<T,A1>::const_iterator RemoteIterator;
468 for(RemoteIterator remote = remoteIndices.
begin(), end =remoteIndices.
end(); remote != end; ++remote) {
470 typedef typename RemoteIndexList::const_iterator RemoteIndexIterator;
472 GlobalIndexList& global = globalMap[remote->first];
473 RemoteIndexList& rList = *(remote->second.first);
475 for(RemoteIndexIterator index = rList.begin(), riEnd = rList.end();
476 index != riEnd; ++index) {
477 global.push_back(std::make_pair(index->localIndexPair().global(),
478 index->localIndexPair().local().attribute()));
491 template<
typename T,
typename A,
typename A1>
493 SLList<std::pair<
typename T::GlobalIndex,
494 typename T::LocalIndex::Attribute>,A> >& globalMap,
500 typedef typename T::GlobalIndex GlobalIndex;
501 typedef typename T::LocalIndex::Attribute Attribute;
502 typedef std::pair<GlobalIndex,Attribute> GlobalIndexPair;
504 typedef typename GlobalIndexPairList::iterator GlobalIndexIterator;
506 assert(globalMap.size()==
static_cast<std::size_t
>(remoteIndices.
neighbours()));
508 typename std::map<int,GlobalIndexPairList>::iterator global = globalMap.begin();
509 RemoteIterator end = remoteIndices.remoteIndices_.end();
511 for(RemoteIterator remote = remoteIndices.remoteIndices_.begin(); remote != end; ++remote, ++global) {
512 typedef typename T::const_iterator IndexIterator;
514 assert(remote->first==global->first);
515 assert(remote->second.first->size() == global->second.size());
517 RemoteIndexIterator riEnd = remote->second.first->end();
518 RemoteIndexIterator rIndex = remote->second.first->begin();
519 GlobalIndexIterator gIndex = global->second.begin();
520 IndexIterator index = indexSet.begin();
522 assert(rIndex==riEnd || gIndex != global->second.end());
523 while(rIndex != riEnd) {
525 assert(gIndex != global->second.end());
527 while(!(index->global() == gIndex->first
528 && index->local().attribute() == gIndex->second)) {
533 if (index->global() > gIndex->first) {
534 index=indexSet.begin();
538 assert(index != indexSet.end() && *index == *gIndex);
540 rIndex->localIndex_ = &(*index);
546 remoteIndices.sourceSeqNo_ = remoteIndices.source_->
seqNo();
547 remoteIndices.destSeqNo_ = remoteIndices.target_->
seqNo();
553 : indexSet_(indexSet), remoteIndices_(remoteIndices)
556 assert(remoteIndices.source_ == remoteIndices.target_);
557 assert(remoteIndices.source_ == &indexSet);
565 : iterators_(remoteIndices.beginModify(), globalIndices.beginModify(),
566 booleans.beginModify(), remoteIndices.end())
577 ++(get<0>(iterators_));
578 ++(get<1>(iterators_));
579 ++(get<2>(iterators_));
585 const std::pair<GlobalIndex,Attribute>& global)
587 get<0>(iterators_).insert(index);
588 get<1>(iterators_).insert(global);
589 get<2>(iterators_).insert(
false);
596 return *(get<0>(iterators_));
603 return *(get<1>(iterators_));
609 return *(get<2>(iterators_));
625 return get<0>(iterators_)!=get<3>(iterators_);
631 return get<0>(iterators_)==get<3>(iterators_);
637 MPI_Datatype type[2] = {MPI_INT, MPI_INT};
638 int blocklength[2] = {1,1};
639 MPI_Aint displacement[2];
643 MessageInformation message;
645 MPI_Address( &(message.publish), displacement);
646 MPI_Address( &(message.pairs), displacement+1);
649 MPI_Address(&message, &base);
650 displacement[0] -= base;
651 displacement[1] -= base;
653 MPI_Type_struct( 2, blocklength, displacement, type, &datatype_);
654 MPI_Type_commit(&datatype_);
658 void IndicesSyncer<T>::calculateMessageSizes()
661 typedef CollectiveIterator<T,typename RemoteIndices::Allocator> CollectiveIterator;
663 IndexIterator iEnd = indexSet_.end();
664 CollectiveIterator collIter = remoteIndices_.template iterator<true>();
666 for(IndexIterator index = indexSet_.begin(); index != iEnd; ++index) {
667 collIter.advance(index->global(), index->local().attribute());
672 typedef typename CollectiveIterator::iterator ValidIterator;
673 ValidIterator end = collIter.end();
676 for(ValidIterator valid = collIter.begin(); valid != end; ++valid) {
681 Dune::dverb<<rank_<<
": publishing "<<knownRemote<<
" for index "<<index->global()<<
" for processes ";
684 for(ValidIterator valid = collIter.begin(); valid != end; ++valid) {
685 ++(infoSend_[valid.process()].publish);
686 (infoSend_[valid.process()].pairs) += knownRemote;
688 Dune::dverb<<
"(publish="<<infoSend_[valid.process()].publish<<
", pairs="<<infoSend_[valid.process()].pairs
695 typedef typename std::map<int,MessageInformation>::const_iterator
698 const MessageIterator end = infoSend_.end();
703 MessageInformation dummy;
705 MessageIterator messageIter= infoSend_.begin();
707 typedef typename RemoteIndices::RemoteIndexMap::const_iterator RemoteIterator;
708 const RemoteIterator rend = remoteIndices_.
end();
711 for(RemoteIterator remote = remoteIndices_.
begin(); remote != rend; ++remote, ++neighbour) {
712 MessageInformation* message;
713 MessageInformation recv;
715 if(messageIter != end && messageIter->first==remote->first) {
717 message =
const_cast<MessageInformation*
>(&(messageIter->second));
723 sendBufferSizes_[neighbour]=0;
726 MPI_Pack_size(1, MPI_INT,remoteIndices_.
communicator(), &tsize);
727 sendBufferSizes_[neighbour] += tsize;
729 for(
int i=0; i < message->publish; ++i) {
731 MPI_Pack_size(1, MPITraits<GlobalIndex>::getType(), remoteIndices_.
communicator(), &tsize);
732 sendBufferSizes_[neighbour] += tsize;
734 MPI_Pack_size(1, MPI_CHAR, remoteIndices_.
communicator(), &tsize);
735 sendBufferSizes_[neighbour] += tsize;
737 MPI_Pack_size(1, MPI_INT, remoteIndices_.
communicator(), &tsize);
738 sendBufferSizes_[neighbour] += tsize;
740 for(
int i=0; i < message->pairs; ++i) {
742 MPI_Pack_size(1, MPI_INT, remoteIndices_.
communicator(), &tsize);
743 sendBufferSizes_[neighbour] += tsize;
745 MPI_Pack_size(1, MPI_CHAR, remoteIndices_.
communicator(), &tsize);
746 sendBufferSizes_[neighbour] += tsize;
749 Dune::dverb<<rank_<<
": Buffer (neighbour="<<remote->first<<
") size is "<< sendBufferSizes_[neighbour]<<
" for publish="<<message->publish<<
" pairs="<<message->pairs<<std::endl;
757 DefaultNumberer numberer;
762 template<
typename T1>
771 typedef typename RemoteIndices::RemoteIndexMap::const_iterator
774 const RemoteIterator end = remoteIndices_.
end();
778 std::size_t noOldNeighbours = remoteIndices_.
neighbours();
779 int* oldNeighbours =
new int[noOldNeighbours];
780 sendBufferSizes_ =
new std::size_t[noOldNeighbours];
783 for(RemoteIterator remote = remoteIndices_.
begin(); remote != end; ++remote, ++i) {
784 typedef typename RemoteIndices::RemoteIndexList::const_iterator
787 oldNeighbours[i]=remote->first;
790 assert(remote->second.first==remote->second.second);
796 BoolList& added = oldMap_[remote->first];
800 index != riEnd; ++index) {
801 global.
push_back(std::make_pair(index->localIndexPair().global(),
802 index->localIndexPair().local().attribute()));
806 Iterators iterators(rList, global, added);
807 iteratorsMap_.insert(std::make_pair(remote->first, iterators));
808 assert(checkReset(iteratorsMap_[remote->first], rList,global,added));
812 calculateMessageSizes();
815 receiveBufferSize_=1;
816 sendBuffers_ =
new char*[noOldNeighbours];
818 for(std::size_t i=0; i<noOldNeighbours; ++i) {
819 sendBuffers_[i] =
new char[sendBufferSizes_[i]];
820 receiveBufferSize_ = std::max(receiveBufferSize_,
static_cast<int>(sendBufferSizes_[i]));
823 receiveBuffer_=
new char[receiveBufferSize_];
825 indexSet_.beginResize();
829 for(i = 0; i<noOldNeighbours; ++i)
834 MPI_Request* requests =
new MPI_Request[noOldNeighbours];
835 MPI_Status* statuses =
new MPI_Status[noOldNeighbours];
838 for(i = 0; i<noOldNeighbours; ++i)
839 packAndSend(oldNeighbours[i], sendBuffers_[i], sendBufferSizes_[i], requests[i]);
842 for(i = 0; i<noOldNeighbours; ++i)
843 recvAndUnpack(numberer);
850 delete[] receiveBuffer_;
854 if(MPI_SUCCESS!=MPI_Waitall(noOldNeighbours, requests, statuses)) {
855 std::cerr<<
": MPI_Error occurred while sending message"<<std::endl;
856 for(i=0; i< noOldNeighbours; i++)
857 if(MPI_SUCCESS!=statuses[i].MPI_ERROR)
858 std::cerr<<
"Destination "<<statuses[i].MPI_SOURCE<<
" error code: "<<statuses[i].MPI_ERROR<<std::endl;
864 for(std::size_t i=0; i<noOldNeighbours; ++i)
865 delete[] sendBuffers_[i];
867 delete[] sendBuffers_;
868 delete[] sendBufferSizes_;
871 iteratorsMap_.clear();
873 indexSet_.endResize();
875 delete[] oldNeighbours;
883 remoteIndices_.sourceSeqNo_ = remoteIndices_.destSeqNo_ = indexSet_.seqNo();
891 IndexIterator iEnd = indexSet_.end();
896 assert(checkReset());
899 MPI_Pack(&(infoSend_[destination].publish), 1, MPI_INT, buffer, bufferSize, &bpos,
902 for(IndexIterator index = indexSet_.begin(); index != iEnd; ++index) {
904 typedef typename IteratorsMap::iterator Iterator;
905 Iterator iteratorsEnd = iteratorsMap_.end();
908 for(Iterator iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators) {
909 while(iterators->second.isNotAtEnd() &&
910 iterators->second.globalIndexPair().first < index->global())
911 ++(iterators->second);
912 assert(!iterators->second.isNotAtEnd() || iterators->second.globalIndexPair().first >= index->global());
919 bool knownRemote =
false;
921 for(Iterator iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
923 std::pair<GlobalIndex,Attribute> p;
924 if (iterators->second.isNotAtEnd())
926 p = iterators->second.globalIndexPair();
929 if(iterators->second.isNotAtEnd() && iterators->second.isOld()
930 && iterators->second.globalIndexPair().first == index->global()) {
932 if(destination == iterators->first)
941 Dune::dverb<<rank_<<
": sending "<<indices<<
" for index "<<index->global()<<
" to "<<destination<<std::endl;
945 MPI_Pack(
const_cast<GlobalIndex*
>(&(index->global())), 1, MPITraits<GlobalIndex>::getType(), buffer, bufferSize, &bpos,
948 char attr = index->local().attribute();
949 MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos,
953 MPI_Pack(&indices, 1, MPI_INT, buffer, bufferSize, &bpos,
957 for(Iterator iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
958 if(iterators->second.isNotAtEnd() && iterators->second.isOld()
959 && iterators->second.globalIndexPair().first == index->global()) {
960 int process = iterators->first;
963 assert(pairs <= infoSend_[destination].pairs);
964 MPI_Pack(&process, 1, MPI_INT, buffer, bufferSize, &bpos,
966 char attr = iterators->second.remoteIndex().attribute();
968 MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos,
974 Dune::dvverb<<
" (publish="<<published<<
", pairs="<<pairs<<
")"<<std::endl;
975 assert(published <= infoSend_[destination].publish);
979 assert(published == infoSend_[destination].publish);
980 assert(pairs == infoSend_[destination].pairs);
983 Dune::dverb << rank_<<
": Sending message of "<<bpos<<
" bytes to "<<destination<<std::endl;
985 MPI_Issend(buffer, bpos, MPI_PACKED, destination, 345, remoteIndices_.
communicator(),&request);
989 inline void IndicesSyncer<T>::insertIntoRemoteIndexList(
int process,
990 const std::pair<GlobalIndex,Attribute>& globalPair,
993 Dune::dverb<<
"Inserting from "<<process<<
" "<<globalPair.first<<
", "<<
994 globalPair.second<<
" "<<attribute<<std::endl;
999 typename IteratorsMap::iterator found = iteratorsMap_.find(process);
1001 if( found == iteratorsMap_.end() ) {
1002 Dune::dverb<<
"Discovered new neighbour "<<process<<std::endl;
1003 RemoteIndexList* rlist =
new RemoteIndexList();
1004 remoteIndices_.remoteIndices_.insert(std::make_pair(process,std::make_pair(rlist,rlist)));
1005 Iterators iterators = Iterators(*rlist, globalMap_[process], oldMap_[process]);
1006 found = iteratorsMap_.insert(std::make_pair(process, iterators)).first;
1009 Iterators& iterators = found->second;
1012 while(iterators.isNotAtEnd() && iterators.globalIndexPair() < globalPair) {
1018 if(iterators.isAtEnd() || iterators.globalIndexPair() != globalPair) {
1021 iterators.insert(RemoteIndex(
Attribute(attribute)),globalPair);
1026 bool indexIsThere=
false;
1027 for(Iterators tmpIterators = iterators;
1028 !tmpIterators.isAtEnd() && tmpIterators.globalIndexPair() == globalPair;
1031 if(tmpIterators.globalIndexPair().second == attribute) {
1039 iterators.insert(RemoteIndex(
Attribute(attribute)),globalPair);
1042 template<
typename T>
1043 template<
typename T1>
1044 void IndicesSyncer<T>::recvAndUnpack(T1& numberer)
1048 IndexIterator iEnd = indexSet_.end();
1049 IndexIterator index = indexSet_.begin();
1053 assert(checkReset());
1058 MPI_Probe(MPI_ANY_SOURCE, 345, remoteIndices_.
communicator(), &status);
1060 int source=status.MPI_SOURCE;
1062 MPI_Get_count(&status, MPI_PACKED, &count);
1064 Dune::dvverb<<rank_<<
": Receiving message from "<< source<<
" with "<<count<<
" bytes"<<std::endl;
1066 if(count>receiveBufferSize_) {
1067 receiveBufferSize_=count;
1068 delete[] receiveBuffer_;
1069 receiveBuffer_ =
new char[receiveBufferSize_];
1072 MPI_Recv(receiveBuffer_, count, MPI_PACKED, source, 345, remoteIndices_.
communicator(), &status);
1075 MPI_Unpack(receiveBuffer_, count, &bpos, &publish, 1, MPI_INT, remoteIndices_.
communicator());
1082 char sourceAttribute;
1085 MPI_Unpack(receiveBuffer_, count, &bpos, &global, 1, MPITraits<GlobalIndex>::getType(),
1087 MPI_Unpack(receiveBuffer_, count, &bpos, &sourceAttribute, 1, MPI_CHAR,
1089 MPI_Unpack(receiveBuffer_, count, &bpos, &pairs, 1, MPI_INT,
1094 SLList<std::pair<int,Attribute> > sourceAttributeList;
1095 sourceAttributeList.push_back(std::make_pair(source,
Attribute(sourceAttribute)));
1097 bool foundSelf =
false;
1102 for(; pairs>0; --pairs) {
1106 MPI_Unpack(receiveBuffer_, count, &bpos, &process, 1, MPI_INT,
1109 MPI_Unpack(receiveBuffer_, count, &bpos, &attribute, 1, MPI_CHAR,
1112 if(process==rank_) {
1120 IndexIterator pos = std::lower_bound(index, iEnd,
IndexPair(global));
1122 if(pos == iEnd || pos->global() != global) {
1124 indexSet_.add(global,
1125 ParallelLocalIndex<Attribute>(numberer(global),
1126 myAttribute,
true));
1127 Dune::dvverb <<
"Adding "<<global<<
" "<<myAttribute<<std::endl;
1132 bool indexIsThere =
false;
1135 for(; pos->global()==global; ++pos)
1136 if(pos->local().attribute() == myAttribute) {
1137 Dune::dvverb<<
"found "<<global<<
" "<<myAttribute<<std::endl;
1138 indexIsThere =
true;
1143 indexSet_.add(global,
1144 ParallelLocalIndex<Attribute>(numberer(global),
1145 myAttribute,
true));
1146 Dune::dvverb <<
"Adding "<<global<<
" "<<myAttribute<<std::endl;
1150 sourceAttributeList.push_back(std::make_pair(process,
Attribute(attribute)));
1155 typedef typename SLList<std::pair<int,Attribute> >::const_iterator Iter;
1156 for(Iter i=sourceAttributeList.begin(), end=sourceAttributeList.end();
1158 insertIntoRemoteIndexList(i->first, std::make_pair(global, myAttribute),
1163 resetIteratorsMap();
1166 template<
typename T>
1167 void IndicesSyncer<T>::resetIteratorsMap(){
1170 typedef typename IteratorsMap::iterator Iterator;
1171 typedef typename RemoteIndices::RemoteIndexMap::iterator
1173 typedef typename GlobalIndicesMap::iterator GlobalIterator;
1174 typedef typename BoolMap::iterator BoolIterator;
1176 const RemoteIterator remoteEnd = remoteIndices_.remoteIndices_.end();
1177 Iterator iterators = iteratorsMap_.begin();
1178 GlobalIterator global = globalMap_.begin();
1179 BoolIterator added = oldMap_.begin();
1181 for(RemoteIterator remote = remoteIndices_.remoteIndices_.begin();
1182 remote != remoteEnd; ++remote, ++global, ++added, ++iterators) {
1183 iterators->second.reset(*(remote->second.first), global->second, added->second);
1187 template<
typename T>
1188 bool IndicesSyncer<T>::checkReset(
const Iterators& iterators, RemoteIndexList& rList, GlobalIndexList& gList,
1191 if(get<0>(iterators.iterators_)!=rList.begin())
1193 if(get<1>(iterators.iterators_)!=gList.begin())
1195 if(get<2>(iterators.iterators_)!=bList.begin())
1201 template<
typename T>
1202 bool IndicesSyncer<T>::checkReset(){
1205 typedef typename IteratorsMap::iterator Iterator;
1206 typedef typename RemoteIndices::RemoteIndexMap::iterator
1208 typedef typename GlobalIndicesMap::iterator GlobalIterator;
1209 typedef typename BoolMap::iterator BoolIterator;
1211 const RemoteIterator remoteEnd = remoteIndices_.remoteIndices_.end();
1212 Iterator iterators = iteratorsMap_.begin();
1213 GlobalIterator global = globalMap_.begin();
1214 BoolIterator added = oldMap_.begin();
1217 for(RemoteIterator remote = remoteIndices_.remoteIndices_.begin();
1218 remote != remoteEnd; ++remote, ++global, ++added, ++iterators) {
1219 if(!checkReset(iterators->second, *(remote->second.first), global->second,
A constant random access iterator for the Dune::ArrayList class.
Definition: arraylist.hh:380
A pair consisting of a global and local index.
Definition: indexset.hh:85
Class for recomputing missing indices of a distributed index set.
Definition: indicessyncer.hh:41
Information about an index residing on another processor.
Definition: remoteindices.hh:68
The indices present on remote processes.
Definition: remoteindices.hh:183
MPI_Comm communicator() const
Get the mpi communicator used.
Definition: remoteindices.hh:1698
const_iterator end() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1525
int neighbours() const
Get the number of processors we share indices with.
Definition: remoteindices.hh:1442
const_iterator begin() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1518
A::template rebind< RemoteIndex >::other Allocator
The type of the allocator for the remote index list.
Definition: remoteindices.hh:231
A mutable iterator for the SLList.
Definition: sllist.hh:270
A single linked list.
Definition: sllist.hh:43
A Tuple of objects.
Definition: tuples.hh:294
Dune::RemoteIndices< ParallelIndexSet > RemoteIndices
Type of the remote indices.
Definition: indicessyncer.hh:59
ArrayList< IndexPair, N >::const_iterator const_iterator
The constant iterator over the pairs.
Definition: indexset.hh:306
int publish
The number of indices we publish for the other process.
Definition: indicessyncer.hh:126
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:492
ParallelIndexSet::GlobalIndex GlobalIndex
Type of the global index used in the index set.
Definition: indicessyncer.hh:51
void storeGlobalIndicesOfRemoteIndices(std::map< int, SLList< std::pair< typename T::GlobalIndex, typename T::LocalIndex::Attribute >, A > > &globalMap, const RemoteIndices< T, A1 > &remoteIndices, const T &indexSet)
Stores the corresponding global indices of the remote index information.
Definition: indicessyncer.hh:462
int pairs
The number of pairs (attribute and process number) we publish to the neighbour process.
Definition: indicessyncer.hh:131
ParallelIndexSet::LocalIndex::Attribute Attribute
Type of the attribute used in the index set.
Definition: indicessyncer.hh:54
IndicesSyncer(ParallelIndexSet &indexSet, RemoteIndices &remoteIndices)
Constructor.
Definition: indicessyncer.hh:551
std::size_t operator()(const GlobalIndex &global)
Provide the lcoal index, always std::numeric_limits<size_t>::max()
Definition: indicessyncer.hh:145
T ParallelIndexSet
The type of the index set.
Definition: indicessyncer.hh:45
TG GlobalIndex
the type of the global index. This type has to provide at least a operator< for sorting.
Definition: indexset.hh:226
void sync()
Sync the index set.
Definition: indicessyncer.hh:755
int seqNo() const
Get the internal sequence number.
ParallelIndexSet::IndexPair IndexPair
The type of the index pair.
Definition: indicessyncer.hh:48
void push_back(const MemberType &item)
Add a new entry to the end of the list.
Definition: sllist.hh:658
SLListIterator< T, A > iterator
The mutable iterator of the list.
Definition: sllist.hh:68
iterator end()
Get an iterator pointing to the end of the list.
Definition: sllist.hh:789
ModifyIterator beginModify()
Get an iterator capable of deleting and inserting elements.
Definition: sllist.hh:802
SLListConstIterator< RemoteIndex, Allocator > const_iterator
The constant iterator of the list.
Definition: sllist.hh:73
SLListModifyIterator< T, A > ModifyIterator
The type of the iterator capable of deletion and insertion.
Definition: sllist.hh:102
iterator begin()
Get an iterator pointing to the first element in the list.
Definition: sllist.hh:777
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:231
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:253
DVVerbType dvverb(std::cout)
stream for very verbose output.
Definition: stdstreams.hh:94
DVerbType dverb(std::cout)
Singleton of verbose debug stream.
Definition: stdstreams.hh:115
Provides a map between global and local indices.
Dune namespace.
Definition: alignment.hh:14
Classes describing a distributed indexset.
Implements a singly linked list together with the necessary iterators.
Standard Dune debug streams.
Fallback implementation of the std::tuple class.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18