5#ifndef DUNE_COMMON_PARALLEL_INDICESSYNCER_HH
6#define DUNE_COMMON_PARALLEL_INDICESSYNCER_HH
58 typedef typename ParallelIndexSet::LocalIndex::Attribute
Attribute;
100 template<
typename T1>
101 void sync(T1& numberer,
bool useFixedOrder =
false);
115 char* receiveBuffer_;
118 std::size_t* sendBufferSizes_;
121 int receiveBufferSize_;
126 struct MessageInformation
143 class DefaultNumberer
158 MPI_Datatype datatype_;
179 typedef std::map<int, GlobalIndexList> GlobalIndicesMap;
189 GlobalIndicesMap globalMap_;
205 typedef std::map<int,BoolList> BoolMap;
214 std::map<int,MessageInformation> infoSend_;
232 typedef std::tuple<RemoteIndexModifier,GlobalIndexModifier,BoolListModifier,
233 const ConstRemoteIndexIterator> IteratorTuple;
255 Iterators(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
266 Iterators& operator++();
274 const std::pair<GlobalIndex,Attribute>& global);
286 std::pair<GlobalIndex,Attribute>& globalIndexPair()
const;
306 void reset(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
314 bool isNotAtEnd()
const;
321 bool isAtEnd()
const;
333 IteratorTuple iterators_;
337 typedef std::map<int,Iterators> IteratorsMap;
350 IteratorsMap iteratorsMap_;
353 void calculateMessageSizes();
362 void packAndSend(
int destination,
char* buffer, std::size_t bufferSize, MPI_Request& req);
370 template<
typename T1>
371 void recvAndUnpack(T1& numberer,
int hardSource,
bool useHardSource);
376 void registerMessageDatatype();
381 void insertIntoRemoteIndexList(
int process,
382 const std::pair<GlobalIndex,Attribute>& global,
388 void resetIteratorsMap();
404 bool checkReset(
const Iterators& iterators, RemoteIndexList& rlist, GlobalIndexList& gList,
408 template<
typename TG,
typename TA>
409 bool operator<(
const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
410 const std::pair<TG,TA>& i2)
412 return i1.global() < i2.first ||
413 (i1.global() == i2.first && i1.local().attribute()<i2.second);
416 template<
typename TG,
typename TA>
417 bool operator<(
const std::pair<TG,TA>& i1,
418 const IndexPair<TG,ParallelLocalIndex<TA> >& i2)
420 return i1.first < i2.global() ||
421 (i1.first == i2.global() && i1.second<i2.local().attribute());
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>
432 bool operator!=(
const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
433 const std::pair<TG,TA>& i2)
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);
445 template<
typename TG,
typename TA>
447 const IndexPair<TG,ParallelLocalIndex<TA> >& i1)
449 return (i1.global() != i2.first || i1.local().attribute()!=i2.second);
468 template<
typename T,
typename A,
typename A1>
472 for(
auto remote = remoteIndices.
begin(), end =remoteIndices.
end(); remote != end; ++remote) {
475 GlobalIndexList& global = globalMap[remote->first];
476 RemoteIndexList& rList = *(remote->second.first);
478 for(
auto index = rList.begin(), riEnd = rList.end();
479 index != riEnd; ++index) {
480 global.push_back(std::make_pair(index->localIndexPair().global(),
481 index->localIndexPair().local().attribute()));
494 template<
typename T,
typename A,
typename A1>
496 SLList<std::pair<
typename T::GlobalIndex,
497 typename T::LocalIndex::Attribute>,A> >& globalMap,
501 assert(globalMap.size()==
static_cast<std::size_t
>(remoteIndices.
neighbours()));
503 auto global = globalMap.begin();
504 auto end = remoteIndices.remoteIndices_.end();
506 for(
auto remote = remoteIndices.remoteIndices_.begin(); remote != end; ++remote, ++global) {
507 assert(remote->first==global->first);
508 assert(remote->second.first->size() == global->second.size());
510 auto riEnd = remote->second.first->end();
511 auto rIndex = remote->second.first->begin();
512 auto gIndex = global->second.begin();
513 auto index = indexSet.begin();
515 assert(rIndex==riEnd || gIndex != global->second.end());
516 while(rIndex != riEnd) {
518 assert(gIndex != global->second.end());
520 while(!(index->global() == gIndex->first
521 && index->local().attribute() == gIndex->second)) {
526 if (index->global() > gIndex->first) {
527 index=indexSet.begin();
531 assert(index != indexSet.end() && *index == *gIndex);
533 rIndex->localIndex_ = &(*index);
539 remoteIndices.sourceSeqNo_ = remoteIndices.source_->
seqNo();
540 remoteIndices.destSeqNo_ = remoteIndices.target_->
seqNo();
546 : indexSet_(indexSet), remoteIndices_(remoteIndices)
549 assert(remoteIndices.source_ == remoteIndices.target_);
550 assert(remoteIndices.source_ == &indexSet);
558 : iterators_(remoteIndices.beginModify(), globalIndices.beginModify(),
559 booleans.beginModify(), remoteIndices.end())
570 ++(std::get<0>(iterators_));
571 ++(std::get<1>(iterators_));
572 ++(std::get<2>(iterators_));
578 const std::pair<GlobalIndex,Attribute>& global)
580 std::get<0>(iterators_).insert(index);
581 std::get<1>(iterators_).insert(global);
582 std::get<2>(iterators_).insert(
false);
589 return *(std::get<0>(iterators_));
596 return *(std::get<1>(iterators_));
602 return *(std::get<2>(iterators_));
610 std::get<0>(iterators_) = remoteIndices.
beginModify();
611 std::get<1>(iterators_) = globalIndices.
beginModify();
618 return std::get<0>(iterators_) != std::get<3>(iterators_);
624 return std::get<0>(iterators_) == std::get<3>(iterators_);
630 MPI_Datatype type[2] = {MPI_INT, MPI_INT};
631 int blocklength[2] = {1,1};
632 MPI_Aint displacement[2];
636 MessageInformation message;
638 MPI_Get_address( &(message.publish), displacement);
639 MPI_Get_address( &(message.pairs), displacement+1);
642 MPI_Get_address(&message, &base);
643 displacement[0] -= base;
644 displacement[1] -= base;
646 MPI_Type_create_struct( 2, blocklength, displacement, type, &datatype_);
647 MPI_Type_commit(&datatype_);
651 void IndicesSyncer<T>::calculateMessageSizes()
653 auto iEnd = indexSet_.end();
654 auto collIter = remoteIndices_.template iterator<true>();
656 for(
auto index = indexSet_.begin(); index != iEnd; ++index) {
657 collIter.advance(index->global(), index->local().attribute());
661 auto end = collIter.end();
664 for(
auto valid = collIter.begin(); valid != end; ++valid) {
669 Dune::dverb<<rank_<<
": publishing "<<knownRemote<<
" for index "<<index->global()<<
" for processes ";
672 for(
auto valid = collIter.begin(); valid != end; ++valid) {
673 ++(infoSend_[valid.process()].publish);
674 (infoSend_[valid.process()].pairs) += knownRemote;
676 Dune::dverb<<
"(publish="<<infoSend_[valid.process()].publish<<
", pairs="<<infoSend_[valid.process()].pairs
683 const auto end = infoSend_.end();
686 MessageInformation dummy;
688 auto messageIter= infoSend_.begin();
689 const auto rend = remoteIndices_.
end();
692 for(
auto remote = remoteIndices_.
begin(); remote != rend; ++remote, ++neighbour) {
693 MessageInformation* message;
694 MessageInformation recv;
696 if(messageIter != end && messageIter->first==remote->first) {
698 message =
const_cast<MessageInformation*
>(&(messageIter->second));
704 sendBufferSizes_[neighbour]=0;
707 MPI_Pack_size(1, MPI_INT,remoteIndices_.
communicator(), &tsize);
708 sendBufferSizes_[neighbour] += tsize;
710 for(
int i=0; i < message->publish; ++i) {
712 MPI_Pack_size(1, MPITraits<GlobalIndex>::getType(), remoteIndices_.
communicator(), &tsize);
713 sendBufferSizes_[neighbour] += tsize;
715 MPI_Pack_size(1, MPI_CHAR, remoteIndices_.
communicator(), &tsize);
716 sendBufferSizes_[neighbour] += tsize;
718 MPI_Pack_size(1, MPI_INT, remoteIndices_.
communicator(), &tsize);
719 sendBufferSizes_[neighbour] += tsize;
721 for(
int i=0; i < message->pairs; ++i) {
723 MPI_Pack_size(1, MPI_INT, remoteIndices_.
communicator(), &tsize);
724 sendBufferSizes_[neighbour] += tsize;
726 MPI_Pack_size(1, MPI_CHAR, remoteIndices_.
communicator(), &tsize);
727 sendBufferSizes_[neighbour] += tsize;
730 Dune::dverb<<rank_<<
": Buffer (neighbour="<<remote->first<<
") size is "<< sendBufferSizes_[neighbour]<<
" for publish="<<message->publish<<
" pairs="<<message->pairs<<std::endl;
738 DefaultNumberer numberer;
743 template<
typename T1>
750 const auto end = remoteIndices_.
end();
754 std::size_t noOldNeighbours = remoteIndices_.
neighbours();
755 int* oldNeighbours =
new int[noOldNeighbours];
756 sendBufferSizes_ =
new std::size_t[noOldNeighbours];
757 std::size_t neighbourI = 0;
759 for(
auto remote = remoteIndices_.
begin(); remote != end; ++remote, ++neighbourI) {
760 oldNeighbours[neighbourI] = remote->first;
763 assert(remote->second.first==remote->second.second);
769 BoolList& added = oldMap_[remote->first];
770 auto riEnd = rList.
end();
772 for(
auto index = rList.
begin();
773 index != riEnd; ++index) {
774 global.
push_back(std::make_pair(index->localIndexPair().global(),
775 index->localIndexPair().local().attribute()));
779 Iterators iterators(rList, global, added);
780 iteratorsMap_.insert(std::make_pair(remote->first, iterators));
781 assert(checkReset(iteratorsMap_[remote->first], rList,global,added));
785 calculateMessageSizes();
788 receiveBufferSize_=1;
789 sendBuffers_ =
new char*[noOldNeighbours];
791 for(std::size_t i=0; i<noOldNeighbours; ++i) {
792 sendBuffers_[i] =
new char[sendBufferSizes_[i]];
793 receiveBufferSize_ = std::max<int>(receiveBufferSize_, sendBufferSizes_[i]);
796 receiveBuffer_=
new char[receiveBufferSize_];
798 indexSet_.beginResize();
802 for(std::size_t i = 0; i<noOldNeighbours; ++i)
807 MPI_Request* requests =
new MPI_Request[noOldNeighbours];
808 MPI_Status* statuses =
new MPI_Status[noOldNeighbours];
811 for(std::size_t i = 0; i<noOldNeighbours; ++i)
812 packAndSend(oldNeighbours[i], sendBuffers_[i], sendBufferSizes_[i], requests[i]);
815 for(std::size_t i = 0; i<noOldNeighbours; ++i)
816 recvAndUnpack(numberer, oldNeighbours[i], useFixedOrder);
823 delete[] receiveBuffer_;
827 if(MPI_SUCCESS!=MPI_Waitall(noOldNeighbours, requests, statuses)) {
828 std::cerr<<
": MPI_Error occurred while sending message"<<std::endl;
829 for(std::size_t i=0; i< noOldNeighbours; i++)
830 if(MPI_SUCCESS!=statuses[i].MPI_ERROR)
831 std::cerr<<
"Destination "<<statuses[i].MPI_SOURCE<<
" error code: "<<statuses[i].MPI_ERROR<<std::endl;
837 for(std::size_t i=0; i<noOldNeighbours; ++i)
838 delete[] sendBuffers_[i];
840 delete[] sendBuffers_;
841 delete[] sendBufferSizes_;
844 iteratorsMap_.clear();
846 indexSet_.endResize();
848 delete[] oldNeighbours;
856 remoteIndices_.sourceSeqNo_ = remoteIndices_.destSeqNo_ = indexSet_.seqNo();
862 auto iEnd = indexSet_.end();
867 assert(checkReset());
870 MPI_Pack(&(infoSend_[destination].publish), 1, MPI_INT, buffer, bufferSize, &bpos,
873 for(
auto index = indexSet_.begin(); index != iEnd; ++index) {
875 auto iteratorsEnd = iteratorsMap_.end();
878 for(
auto iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators) {
879 while(iterators->second.isNotAtEnd() &&
880 iterators->second.globalIndexPair().first < index->global())
881 ++(iterators->second);
882 assert(!iterators->second.isNotAtEnd() || iterators->second.globalIndexPair().first >= index->global());
889 bool knownRemote =
false;
891 for(
auto iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
893 std::pair<GlobalIndex,Attribute> p;
894 if (iterators->second.isNotAtEnd())
896 p = iterators->second.globalIndexPair();
899 if(iterators->second.isNotAtEnd() && iterators->second.isOld()
900 && iterators->second.globalIndexPair().first == index->global()) {
902 if(destination == iterators->first)
911 Dune::dverb<<rank_<<
": sending "<<indices<<
" for index "<<index->global()<<
" to "<<destination<<std::endl;
915 MPI_Pack(
const_cast<GlobalIndex*
>(&(index->global())), 1, MPITraits<GlobalIndex>::getType(), buffer, bufferSize, &bpos,
918 char attr = index->local().attribute();
919 MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos,
923 MPI_Pack(&indices, 1, MPI_INT, buffer, bufferSize, &bpos,
927 for(
auto iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
928 if(iterators->second.isNotAtEnd() && iterators->second.isOld()
929 && iterators->second.globalIndexPair().first == index->global()) {
930 int process = iterators->first;
933 assert(pairs <= infoSend_[destination].pairs);
934 MPI_Pack(&process, 1, MPI_INT, buffer, bufferSize, &bpos,
936 char attr2 = iterators->second.remoteIndex().attribute();
938 MPI_Pack(&attr2, 1, MPI_CHAR, buffer, bufferSize, &bpos,
944 Dune::dvverb<<
" (publish="<<published<<
", pairs="<<pairs<<
")"<<std::endl;
945 assert(published <= infoSend_[destination].publish);
949 assert(published == infoSend_[destination].publish);
950 assert(pairs == infoSend_[destination].pairs);
953 Dune::dverb << rank_<<
": Sending message of "<<bpos<<
" bytes to "<<destination<<std::endl;
955 MPI_Issend(buffer, bpos, MPI_PACKED, destination, 345, remoteIndices_.
communicator(),&request);
959 inline void IndicesSyncer<T>::insertIntoRemoteIndexList(
int process,
960 const std::pair<GlobalIndex,Attribute>& globalPair,
963 Dune::dverb<<
"Inserting from "<<process<<
" "<<globalPair.first<<
", "<<
964 globalPair.second<<
" "<<attribute<<std::endl;
969 typename IteratorsMap::iterator found = iteratorsMap_.find(process);
971 if( found == iteratorsMap_.end() ) {
972 Dune::dverb<<
"Discovered new neighbour "<<process<<std::endl;
973 RemoteIndexList* rlist =
new RemoteIndexList();
974 remoteIndices_.remoteIndices_.insert(std::make_pair(process,std::make_pair(rlist,rlist)));
975 Iterators iterators = Iterators(*rlist, globalMap_[process], oldMap_[process]);
976 found = iteratorsMap_.insert(std::make_pair(process, iterators)).first;
979 Iterators& iterators = found->second;
982 while(iterators.isNotAtEnd() && iterators.globalIndexPair() < globalPair) {
988 if(iterators.isAtEnd() || iterators.globalIndexPair() != globalPair) {
991 iterators.insert(RemoteIndex(
Attribute(attribute)),globalPair);
996 bool indexIsThere=
false;
997 for(Iterators tmpIterators = iterators;
998 !tmpIterators.isAtEnd() && tmpIterators.globalIndexPair() == globalPair;
1001 if(tmpIterators.globalIndexPair().second == attribute) {
1009 iterators.insert(RemoteIndex(
Attribute(attribute)),globalPair);
1012 template<
typename T>
1013 template<
typename T1>
1014 void IndicesSyncer<T>::recvAndUnpack(T1& numberer,
int hardSource,
bool useHardSource)
1017 auto iEnd = constIndexSet.end();
1018 auto index = constIndexSet.begin();
1022 assert(checkReset());
1028 MPI_Probe(useHardSource ? hardSource : MPI_ANY_SOURCE, 345, remoteIndices_.
communicator(), &status);
1030 int source=status.MPI_SOURCE;
1032 MPI_Get_count(&status, MPI_PACKED, &count);
1034 Dune::dvverb<<rank_<<
": Receiving message from "<< source<<
" with "<<count<<
" bytes"<<std::endl;
1036 if(count>receiveBufferSize_) {
1037 receiveBufferSize_=count;
1038 delete[] receiveBuffer_;
1039 receiveBuffer_ =
new char[receiveBufferSize_];
1042 MPI_Recv(receiveBuffer_, count, MPI_PACKED, source, 345, remoteIndices_.
communicator(), &status);
1045 MPI_Unpack(receiveBuffer_, count, &bpos, &publish, 1, MPI_INT, remoteIndices_.
communicator());
1052 char sourceAttribute;
1055 MPI_Unpack(receiveBuffer_, count, &bpos, &global, 1, MPITraits<GlobalIndex>::getType(),
1057 MPI_Unpack(receiveBuffer_, count, &bpos, &sourceAttribute, 1, MPI_CHAR,
1059 MPI_Unpack(receiveBuffer_, count, &bpos, &pairs, 1, MPI_INT,
1064 SLList<std::pair<int,Attribute> > sourceAttributeList;
1065 sourceAttributeList.push_back(std::make_pair(source,
Attribute(sourceAttribute)));
1067 bool foundSelf =
false;
1072 for(; pairs>0; --pairs) {
1076 MPI_Unpack(receiveBuffer_, count, &bpos, &process, 1, MPI_INT,
1079 MPI_Unpack(receiveBuffer_, count, &bpos, &attribute, 1, MPI_CHAR,
1082 if(process==rank_) {
1090 auto pos = std::lower_bound(index, iEnd,
IndexPair(global));
1092 if(pos == iEnd || pos->global() != global) {
1094 indexSet_.add(global,
1095 ParallelLocalIndex<Attribute>(numberer(global),
1096 myAttribute,
true));
1097 Dune::dvverb <<
"Adding "<<global<<
" "<<myAttribute<<std::endl;
1102 bool indexIsThere =
false;
1105 for(; pos->global()==global; ++pos)
1106 if(pos->local().attribute() == myAttribute) {
1107 Dune::dvverb<<
"found "<<global<<
" "<<myAttribute<<std::endl;
1108 indexIsThere =
true;
1113 indexSet_.add(global,
1114 ParallelLocalIndex<Attribute>(numberer(global),
1115 myAttribute,
true));
1116 Dune::dvverb <<
"Adding "<<global<<
" "<<myAttribute<<std::endl;
1120 sourceAttributeList.push_back(std::make_pair(process,
Attribute(attribute)));
1125 typedef typename SLList<std::pair<int,Attribute> >::const_iterator Iter;
1126 for(Iter i=sourceAttributeList.begin(), end=sourceAttributeList.end();
1128 insertIntoRemoteIndexList(i->first, std::make_pair(global, myAttribute),
1133 resetIteratorsMap();
1136 template<
typename T>
1137 void IndicesSyncer<T>::resetIteratorsMap(){
1140 const auto remoteEnd = remoteIndices_.remoteIndices_.end();
1141 auto iterators = iteratorsMap_.begin();
1142 auto global = globalMap_.begin();
1143 auto added = oldMap_.begin();
1145 for(
auto remote = remoteIndices_.remoteIndices_.begin();
1146 remote != remoteEnd; ++remote, ++global, ++added, ++iterators) {
1147 iterators->second.reset(*(remote->second.first), global->second, added->second);
1151 template<
typename T>
1152 bool IndicesSyncer<T>::checkReset(
const Iterators& iterators, RemoteIndexList& rList, GlobalIndexList& gList,
1155 if(std::get<0>(iterators.iterators_) != rList.begin())
1157 if(std::get<1>(iterators.iterators_) != gList.begin())
1159 if(std::get<2>(iterators.iterators_) != bList.begin())
1165 template<
typename T>
1166 bool IndicesSyncer<T>::checkReset(){
1169 const auto remoteEnd = remoteIndices_.remoteIndices_.end();
1170 auto iterators = iteratorsMap_.begin();
1171 auto global = globalMap_.begin();
1172 auto added = oldMap_.begin();
1175 for(
auto remote = remoteIndices_.remoteIndices_.begin();
1176 remote != remoteEnd; ++remote, ++global, ++added, ++iterators) {
1177 if(!checkReset(iterators->second, *(remote->second.first), global->second,
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:45
Information about an index residing on another processor.
Definition: remoteindices.hh:74
The indices present on remote processes.
Definition: remoteindices.hh:190
typename std::allocator_traits< A >::template rebind_alloc< RemoteIndex > Allocator
The type of the allocator for the remote index list.
Definition: remoteindices.hh:238
A mutable iterator for the SLList.
Definition: sllist.hh:271
A single linked list.
Definition: sllist.hh:44
Provides a map between global and local indices.
Dune::RemoteIndices< ParallelIndexSet > RemoteIndices
Type of the remote indices.
Definition: indicessyncer.hh:63
int publish
The number of indices we publish for the other process.
Definition: indicessyncer.hh:132
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:495
MPI_Comm communicator() const
Get the mpi communicator used.
Definition: remoteindices.hh:1695
ParallelIndexSet::GlobalIndex GlobalIndex
Type of the global index used in the index set.
Definition: indicessyncer.hh:55
int pairs
The number of pairs (attribute and process number) we publish to the neighbour process.
Definition: indicessyncer.hh:137
const_iterator end() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1528
ParallelIndexSet::LocalIndex::Attribute Attribute
Type of the attribute used in the index set.
Definition: indicessyncer.hh:58
IndicesSyncer(ParallelIndexSet &indexSet, RemoteIndices &remoteIndices)
Constructor.
Definition: indicessyncer.hh:544
std::size_t operator()(const GlobalIndex &global)
Provide the local index, always std::numeric_limits<size_t>::max()
Definition: indicessyncer.hh:151
T ParallelIndexSet
The type of the index set.
Definition: indicessyncer.hh:49
int neighbours() const
Get the number of processors we share indices with.
Definition: remoteindices.hh:1445
TG GlobalIndex
the type of the global index. This type has to provide at least a operator< for sorting.
Definition: indexset.hh:226
void storeGlobalIndicesOfRemoteIndices(std::map< int, SLList< std::pair< typename T::GlobalIndex, typename T::LocalIndex::Attribute >, A > > &globalMap, const RemoteIndices< T, A1 > &remoteIndices)
Stores the corresponding global indices of the remote index information.
Definition: indicessyncer.hh:469
void sync()
Sync the index set.
Definition: indicessyncer.hh:736
int seqNo() const
Get the internal sequence number.
ParallelIndexSet::IndexPair IndexPair
The type of the index pair.
Definition: indicessyncer.hh:52
const_iterator begin() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1521
void push_back(const MemberType &item)
Add a new entry to the end of the list.
Definition: sllist.hh:643
SLListIterator< T, A > iterator
The mutable iterator of the list.
Definition: sllist.hh:69
iterator end()
Get an iterator pointing to the end of the list.
Definition: sllist.hh:774
ModifyIterator beginModify()
Get an iterator capable of deleting and inserting elements.
Definition: sllist.hh:787
SLListConstIterator< RemoteIndex, Allocator > const_iterator
The constant iterator of the list.
Definition: sllist.hh:74
SLListModifyIterator< T, A > ModifyIterator
The type of the iterator capable of deletion and insertion.
Definition: sllist.hh:103
iterator begin()
Get an iterator pointing to the first element in the list.
Definition: sllist.hh:762
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484
EnableIfInterOperable< T1, T2, bool >::type operator<(const RandomAccessIteratorFacade< T1, V1, R1, D > &lhs, const RandomAccessIteratorFacade< T2, V2, R2, D > &rhs)
Comparison operator.
Definition: iteratorfacades.hh:638
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:238
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:260
DVVerbType dvverb(std::cout)
stream for very verbose output.
Definition: stdstreams.hh:96
DVerbType dverb(std::cout)
Singleton of verbose debug stream.
Definition: stdstreams.hh:117
Dune namespace.
Definition: alignedallocator.hh:13
Classes describing a distributed indexset.
Implements a singly linked list together with the necessary iterators.
Standard Dune debug streams.