Dune Core Modules (2.10.0)

indicessyncer.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// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_COMMON_PARALLEL_INDICESSYNCER_HH
6#define DUNE_COMMON_PARALLEL_INDICESSYNCER_HH
7
8#if HAVE_MPI
9
10#include <algorithm>
11#include <cassert>
12#include <cmath>
13#include <functional>
14#include <limits>
15#include <map>
16#include <tuple>
17
18#include <mpi.h>
19
21#include <dune/common/sllist.hh>
24
25namespace Dune
26{
43 template<typename T>
45 {
46 public:
47
50
53
56
58 typedef typename ParallelIndexSet::LocalIndex::Attribute Attribute;
59
64
75 RemoteIndices& remoteIndices);
76
86 void sync();
87
100 template<typename T1>
101 void sync(T1& numberer, bool useFixedOrder = false);
102
103 private:
104
106 ParallelIndexSet& indexSet_;
107
109 RemoteIndices& remoteIndices_;
110
112 char** sendBuffers_;
113
115 char* receiveBuffer_;
116
118 std::size_t* sendBufferSizes_;
119
121 int receiveBufferSize_; // int because of MPI
122
126 struct MessageInformation
127 {
128 MessageInformation()
129 : publish(), pairs()
130 {}
137 int pairs;
138 };
139
143 class DefaultNumberer
144 {
145 public:
151 std::size_t operator()([[maybe_unused]] const GlobalIndex& global)
152 {
154 }
155 };
156
158 MPI_Datatype datatype_;
159
161 int rank_;
162
167 typedef SLList<std::pair<GlobalIndex,Attribute>, typename RemoteIndices::Allocator> GlobalIndexList;
168
170 typedef typename GlobalIndexList::ModifyIterator GlobalIndexModifier;
171
176 GlobalIndexIterator;
177
179 typedef std::map<int, GlobalIndexList> GlobalIndicesMap;
180
189 GlobalIndicesMap globalMap_;
190
195
199 typedef typename BoolList::iterator BoolIterator;
200
202 typedef typename BoolList::ModifyIterator BoolListModifier;
203
205 typedef std::map<int,BoolList> BoolMap;
206
211 BoolMap oldMap_;
212
214 std::map<int,MessageInformation> infoSend_;
215
217 typedef typename RemoteIndices::RemoteIndexList RemoteIndexList;
218
220 typedef typename RemoteIndexList::ModifyIterator RemoteIndexModifier;
221
224
226 typedef typename RemoteIndexList::iterator RemoteIndexIterator;
227
229 typedef typename RemoteIndexList::const_iterator ConstRemoteIndexIterator;
230
232 typedef std::tuple<RemoteIndexModifier,GlobalIndexModifier,BoolListModifier,
233 const ConstRemoteIndexIterator> IteratorTuple;
234
242 class Iterators
243 {
244 friend class IndicesSyncer<T>;
245 public:
255 Iterators(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
256 BoolList& booleans);
257
261 Iterators();
262
266 Iterators& operator++();
267
273 void insert(const RemoteIndex& index,
274 const std::pair<GlobalIndex,Attribute>& global);
275
280 RemoteIndex& remoteIndex() const;
281
286 std::pair<GlobalIndex,Attribute>& globalIndexPair() const;
287
288 Attribute& attribute() const;
289
295 bool isOld() const;
296
306 void reset(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
307 BoolList& booleans);
308
314 bool isNotAtEnd() const;
315
321 bool isAtEnd() const;
322
323 private:
333 IteratorTuple iterators_;
334 };
335
337 typedef std::map<int,Iterators> IteratorsMap;
338
350 IteratorsMap iteratorsMap_;
351
353 void calculateMessageSizes();
354
362 void packAndSend(int destination, char* buffer, std::size_t bufferSize, MPI_Request& req);
363
370 template<typename T1>
371 void recvAndUnpack(T1& numberer, int hardSource, bool useHardSource);
372
376 void registerMessageDatatype();
377
381 void insertIntoRemoteIndexList(int process,
382 const std::pair<GlobalIndex,Attribute>& global,
383 char attribute);
384
388 void resetIteratorsMap();
389
394 bool checkReset();
395
404 bool checkReset(const Iterators& iterators, RemoteIndexList& rlist, GlobalIndexList& gList,
405 BoolList& bList);
406 };
407
408 template<typename TG, typename TA>
409 bool operator<(const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
410 const std::pair<TG,TA>& i2)
411 {
412 return i1.global() < i2.first ||
413 (i1.global() == i2.first && i1.local().attribute()<i2.second);
414 }
415
416 template<typename TG, typename TA>
417 bool operator<(const std::pair<TG,TA>& i1,
418 const IndexPair<TG,ParallelLocalIndex<TA> >& i2)
419 {
420 return i1.first < i2.global() ||
421 (i1.first == i2.global() && i1.second<i2.local().attribute());
422 }
423
424 template<typename TG, typename TA>
425 bool operator==(const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
426 const std::pair<TG,TA>& i2)
427 {
428 return (i1.global() == i2.first && i1.local().attribute()==i2.second);
429 }
430
431 template<typename TG, typename TA>
432 bool operator!=(const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
433 const std::pair<TG,TA>& i2)
434 {
435 return (i1.global() != i2.first || i1.local().attribute()!=i2.second);
436 }
437
438 template<typename TG, typename TA>
439 bool operator==(const std::pair<TG,TA>& i2,
440 const IndexPair<TG,ParallelLocalIndex<TA> >& i1)
441 {
442 return (i1.global() == i2.first && i1.local().attribute()==i2.second);
443 }
444
445 template<typename TG, typename TA>
446 bool operator!=(const std::pair<TG,TA>& i2,
447 const IndexPair<TG,ParallelLocalIndex<TA> >& i1)
448 {
449 return (i1.global() != i2.first || i1.local().attribute()!=i2.second);
450 }
451
468 template<typename T, typename A, typename A1>
469 void storeGlobalIndicesOfRemoteIndices(std::map<int,SLList<std::pair<typename T::GlobalIndex, typename T::LocalIndex::Attribute>,A> >& globalMap,
470 const RemoteIndices<T,A1>& remoteIndices)
471 {
472 for(auto remote = remoteIndices.begin(), end =remoteIndices.end(); remote != end; ++remote) {
473 typedef typename RemoteIndices<T,A1>::RemoteIndexList RemoteIndexList;
475 GlobalIndexList& global = globalMap[remote->first];
476 RemoteIndexList& rList = *(remote->second.first);
477
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()));
482 }
483 }
484 }
485
494 template<typename T, typename A, typename A1>
495 inline void repairLocalIndexPointers(std::map<int,
496 SLList<std::pair<typename T::GlobalIndex,
497 typename T::LocalIndex::Attribute>,A> >& globalMap,
498 RemoteIndices<T,A1>& remoteIndices,
499 const T& indexSet)
500 {
501 assert(globalMap.size()==static_cast<std::size_t>(remoteIndices.neighbours()));
502 // Repair pointers to index set in remote indices.
503 auto global = globalMap.begin();
504 auto end = remoteIndices.remoteIndices_.end();
505
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());
509
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();
514
515 assert(rIndex==riEnd || gIndex != global->second.end());
516 while(rIndex != riEnd) {
517 // Search for the index in the set.
518 assert(gIndex != global->second.end());
519
520 while(!(index->global() == gIndex->first
521 && index->local().attribute() == gIndex->second)) {
522 ++index;
523 // this is only needed for ALU, where there may exist
524 // more entries with the same global index in the remote index set
525 // than in the index set
526 if (index->global() > gIndex->first) {
527 index=indexSet.begin();
528 }
529 }
530
531 assert(index != indexSet.end() && *index == *gIndex);
532
533 rIndex->localIndex_ = &(*index);
534 ++index;
535 ++rIndex;
536 ++gIndex;
537 }
538 }
539 remoteIndices.sourceSeqNo_ = remoteIndices.source_->seqNo();
540 remoteIndices.destSeqNo_ = remoteIndices.target_->seqNo();
541 }
542
543 template<typename T>
545 RemoteIndices& remoteIndices)
546 : indexSet_(indexSet), remoteIndices_(remoteIndices)
547 {
548 // index sets must match.
549 assert(remoteIndices.source_ == remoteIndices.target_);
550 assert(remoteIndices.source_ == &indexSet);
551 MPI_Comm_rank(remoteIndices_.communicator(), &rank_);
552 }
553
554 template<typename T>
556 GlobalIndexList& globalIndices,
557 BoolList& booleans)
558 : iterators_(remoteIndices.beginModify(), globalIndices.beginModify(),
559 booleans.beginModify(), remoteIndices.end())
560 { }
561
562 template<typename T>
564 : iterators_()
565 {}
566
567 template<typename T>
569 {
570 ++(std::get<0>(iterators_));
571 ++(std::get<1>(iterators_));
572 ++(std::get<2>(iterators_));
573 return *this;
574 }
575
576 template<typename T>
578 const std::pair<GlobalIndex,Attribute>& global)
579 {
580 std::get<0>(iterators_).insert(index);
581 std::get<1>(iterators_).insert(global);
582 std::get<2>(iterators_).insert(false);
583 }
584
585 template<typename T>
586 inline typename IndicesSyncer<T>::RemoteIndex&
588 {
589 return *(std::get<0>(iterators_));
590 }
591
592 template<typename T>
593 inline std::pair<typename IndicesSyncer<T>::GlobalIndex,typename IndicesSyncer<T>::Attribute>&
595 {
596 return *(std::get<1>(iterators_));
597 }
598
599 template<typename T>
601 {
602 return *(std::get<2>(iterators_));
603 }
604
605 template<typename T>
607 GlobalIndexList& globalIndices,
608 BoolList& booleans)
609 {
610 std::get<0>(iterators_) = remoteIndices.beginModify();
611 std::get<1>(iterators_) = globalIndices.beginModify();
612 std::get<2>(iterators_) = booleans.beginModify();
613 }
614
615 template<typename T>
617 {
618 return std::get<0>(iterators_) != std::get<3>(iterators_);
619 }
620
621 template<typename T>
623 {
624 return std::get<0>(iterators_) == std::get<3>(iterators_);
625 }
626
627 template<typename T>
629 {
630 MPI_Datatype type[2] = {MPI_INT, MPI_INT};
631 int blocklength[2] = {1,1};
632 MPI_Aint displacement[2];
633 MPI_Aint base;
634
635 // Compute displacement
636 MessageInformation message;
637
638 MPI_Get_address( &(message.publish), displacement);
639 MPI_Get_address( &(message.pairs), displacement+1);
640
641 // Make the displacement relative
642 MPI_Get_address(&message, &base);
643 displacement[0] -= base;
644 displacement[1] -= base;
645
646 MPI_Type_create_struct( 2, blocklength, displacement, type, &datatype_);
647 MPI_Type_commit(&datatype_);
648 }
649
650 template<typename T>
651 void IndicesSyncer<T>::calculateMessageSizes()
652 {
653 auto iEnd = indexSet_.end();
654 auto collIter = remoteIndices_.template iterator<true>();
655
656 for(auto index = indexSet_.begin(); index != iEnd; ++index) {
657 collIter.advance(index->global(), index->local().attribute());
658 if(collIter.empty())
659 break;
660 int knownRemote=0;
661 auto end = collIter.end();
662
663 // Count the remote indices we know.
664 for(auto valid = collIter.begin(); valid != end; ++valid) {
665 ++knownRemote;
666 }
667
668 if(knownRemote>0) {
669 Dune::dverb<<rank_<<": publishing "<<knownRemote<<" for index "<<index->global()<< " for processes ";
670
671 // Update MessageInformation
672 for(auto valid = collIter.begin(); valid != end; ++valid) {
673 ++(infoSend_[valid.process()].publish);
674 (infoSend_[valid.process()].pairs) += knownRemote;
675 Dune::dverb<<valid.process()<<" ";
676 Dune::dverb<<"(publish="<<infoSend_[valid.process()].publish<<", pairs="<<infoSend_[valid.process()].pairs
677 <<") ";
678 }
679 Dune::dverb<<std::endl;
680 }
681 }
682
683 const auto end = infoSend_.end();
684
685 // Now determine the buffersizes needed for each neighbour using MPI_Pack_size
686 MessageInformation dummy;
687
688 auto messageIter= infoSend_.begin();
689 const auto rend = remoteIndices_.end();
690 int neighbour=0;
691
692 for(auto remote = remoteIndices_.begin(); remote != rend; ++remote, ++neighbour) {
693 MessageInformation* message;
694 MessageInformation recv;
695
696 if(messageIter != end && messageIter->first==remote->first) {
697 // We want to send message information to that process
698 message = const_cast<MessageInformation*>(&(messageIter->second));
699 ++messageIter;
700 }else
701 // We do not want to send information but the other process might.
702 message = &dummy;
703
704 sendBufferSizes_[neighbour]=0;
705 int tsize;
706 // The number of indices published
707 MPI_Pack_size(1, MPI_INT,remoteIndices_.communicator(), &tsize);
708 sendBufferSizes_[neighbour] += tsize;
709
710 for(int i=0; i < message->publish; ++i) {
711 // The global index
712 MPI_Pack_size(1, MPITraits<GlobalIndex>::getType(), remoteIndices_.communicator(), &tsize);
713 sendBufferSizes_[neighbour] += tsize;
714 // The attribute in the local index
715 MPI_Pack_size(1, MPI_CHAR, remoteIndices_.communicator(), &tsize);
716 sendBufferSizes_[neighbour] += tsize;
717 // The number of corresponding remote indices
718 MPI_Pack_size(1, MPI_INT, remoteIndices_.communicator(), &tsize);
719 sendBufferSizes_[neighbour] += tsize;
720 }
721 for(int i=0; i < message->pairs; ++i) {
722 // The process of the remote index
723 MPI_Pack_size(1, MPI_INT, remoteIndices_.communicator(), &tsize);
724 sendBufferSizes_[neighbour] += tsize;
725 // The attribute of the remote index
726 MPI_Pack_size(1, MPI_CHAR, remoteIndices_.communicator(), &tsize);
727 sendBufferSizes_[neighbour] += tsize;
728 }
729
730 Dune::dverb<<rank_<<": Buffer (neighbour="<<remote->first<<") size is "<< sendBufferSizes_[neighbour]<<" for publish="<<message->publish<<" pairs="<<message->pairs<<std::endl;
731 }
732
733 }
734
735 template<typename T>
737 {
738 DefaultNumberer numberer;
739 sync(numberer);
740 }
741
742 template<typename T>
743 template<typename T1>
744 void IndicesSyncer<T>::sync(T1& numberer, bool useFixedOrder)
745 {
746 // The pointers to the local indices in the remote indices
747 // will become invalid due to the resorting of the index set.
748 // Therefore store the corresponding global indices.
749 // Mark all indices as not added
750 const auto end = remoteIndices_.end();
751
752 // Number of neighbours might change during the syncing.
753 // save the old neighbours
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;
758
759 for(auto remote = remoteIndices_.begin(); remote != end; ++remote, ++neighbourI) {
760 oldNeighbours[neighbourI] = remote->first;
761
762 // Make sure we only have one remote index list.
763 assert(remote->second.first==remote->second.second);
764
765 RemoteIndexList& rList = *(remote->second.first);
766
767 // Store the corresponding global indices.
768 GlobalIndexList& global = globalMap_[remote->first];
769 BoolList& added = oldMap_[remote->first];
770 auto riEnd = rList.end();
771
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()));
776 added.push_back(true);
777 }
778
779 Iterators iterators(rList, global, added);
780 iteratorsMap_.insert(std::make_pair(remote->first, iterators));
781 assert(checkReset(iteratorsMap_[remote->first], rList,global,added));
782 }
783
784 // Exchange indices with each neighbour
785 calculateMessageSizes();
786
787 // Allocate the buffers
788 receiveBufferSize_=1;
789 sendBuffers_ = new char*[noOldNeighbours];
790
791 for(std::size_t i=0; i<noOldNeighbours; ++i) {
792 sendBuffers_[i] = new char[sendBufferSizes_[i]];
793 receiveBufferSize_ = std::max(receiveBufferSize_, static_cast<int>(sendBufferSizes_[i]));
794 }
795
796 receiveBuffer_=new char[receiveBufferSize_];
797
798 indexSet_.beginResize();
799
800 Dune::dverb<<rank_<<": Neighbours: ";
801
802 for(std::size_t i = 0; i<noOldNeighbours; ++i)
803 Dune::dverb<<oldNeighbours[i]<<" ";
804
805 Dune::dverb<<std::endl;
806
807 MPI_Request* requests = new MPI_Request[noOldNeighbours];
808 MPI_Status* statuses = new MPI_Status[noOldNeighbours];
809
810 // Pack Message data and start the sends
811 for(std::size_t i = 0; i<noOldNeighbours; ++i)
812 packAndSend(oldNeighbours[i], sendBuffers_[i], sendBufferSizes_[i], requests[i]);
813
814 // Probe for incoming messages, receive and unpack them
815 for(std::size_t i = 0; i<noOldNeighbours; ++i)
816 recvAndUnpack(numberer, oldNeighbours[i], useFixedOrder);
817 // }else{
818 // recvAndUnpack(oldNeighbours[i], numberer);
819 // packAndSend(oldNeighbours[i]);
820 // }
821 // }
822
823 delete[] receiveBuffer_;
824
825 // Wait for the completion of the sends
826 // Wait for completion of sends
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;
832 }
833
834 delete[] statuses;
835 delete[] requests;
836
837 for(std::size_t i=0; i<noOldNeighbours; ++i)
838 delete[] sendBuffers_[i];
839
840 delete[] sendBuffers_;
841 delete[] sendBufferSizes_;
842
843 // No need for the iterator tuples any more
844 iteratorsMap_.clear();
845
846 indexSet_.endResize();
847
848 delete[] oldNeighbours;
849
850 repairLocalIndexPointers(globalMap_, remoteIndices_, indexSet_);
851
852 oldMap_.clear();
853 globalMap_.clear();
854
855 // update the sequence number
856 remoteIndices_.sourceSeqNo_ = remoteIndices_.destSeqNo_ = indexSet_.seqNo();
857 }
858
859 template<typename T>
860 void IndicesSyncer<T>::packAndSend(int destination, char* buffer, std::size_t bufferSize, MPI_Request& request)
861 {
862 auto iEnd = indexSet_.end();
863 int bpos = 0;
864 int published = 0;
865 int pairs = 0;
866
867 assert(checkReset());
868
869 // Pack the number of indices we publish
870 MPI_Pack(&(infoSend_[destination].publish), 1, MPI_INT, buffer, bufferSize, &bpos,
871 remoteIndices_.communicator());
872
873 for(auto index = indexSet_.begin(); index != iEnd; ++index) {
874 // Search for corresponding remote indices in all iterator tuples
875 auto iteratorsEnd = iteratorsMap_.end();
876
877 // advance all iterators to a position with global index >= index->global()
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());
883 }
884
885 // Add all remote indices positioned at global which were already present before calling sync
886 // to the message.
887 // Count how many remote indices we will send
888 int indices = 0;
889 bool knownRemote = false; // Is the remote process supposed to know this index?
890
891 for(auto iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
892 {
893 std::pair<GlobalIndex,Attribute> p;
894 if (iterators->second.isNotAtEnd())
895 {
896 p = iterators->second.globalIndexPair();
897 }
898
899 if(iterators->second.isNotAtEnd() && iterators->second.isOld()
900 && iterators->second.globalIndexPair().first == index->global()) {
901 indices++;
902 if(destination == iterators->first)
903 knownRemote = true;
904 }
905 }
906
907 if(!knownRemote)
908 // We do not need to send any indices
909 continue;
910
911 Dune::dverb<<rank_<<": sending "<<indices<<" for index "<<index->global()<<" to "<<destination<<std::endl;
912
913
914 // Pack the global index, the attribute and the number
915 MPI_Pack(const_cast<GlobalIndex*>(&(index->global())), 1, MPITraits<GlobalIndex>::getType(), buffer, bufferSize, &bpos,
916 remoteIndices_.communicator());
917
918 char attr = index->local().attribute();
919 MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos,
920 remoteIndices_.communicator());
921
922 // Pack the number of remote indices we send.
923 MPI_Pack(&indices, 1, MPI_INT, buffer, bufferSize, &bpos,
924 remoteIndices_.communicator());
925
926 // Pack the information about the remote indices
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;
931
932 ++pairs;
933 assert(pairs <= infoSend_[destination].pairs);
934 MPI_Pack(&process, 1, MPI_INT, buffer, bufferSize, &bpos,
935 remoteIndices_.communicator());
936 char attr2 = iterators->second.remoteIndex().attribute();
937
938 MPI_Pack(&attr2, 1, MPI_CHAR, buffer, bufferSize, &bpos,
939 remoteIndices_.communicator());
940 --indices;
941 }
942 assert(indices==0);
943 ++published;
944 Dune::dvverb<<" (publish="<<published<<", pairs="<<pairs<<")"<<std::endl;
945 assert(published <= infoSend_[destination].publish);
946 }
947
948 // Make sure we send all expected entries
949 assert(published == infoSend_[destination].publish);
950 assert(pairs == infoSend_[destination].pairs);
951 resetIteratorsMap();
952
953 Dune::dverb << rank_<<": Sending message of "<<bpos<<" bytes to "<<destination<<std::endl;
954
955 MPI_Issend(buffer, bpos, MPI_PACKED, destination, 345, remoteIndices_.communicator(),&request);
956 }
957
958 template<typename T>
959 inline void IndicesSyncer<T>::insertIntoRemoteIndexList(int process,
960 const std::pair<GlobalIndex,Attribute>& globalPair,
961 char attribute)
962 {
963 Dune::dverb<<"Inserting from "<<process<<" "<<globalPair.first<<", "<<
964 globalPair.second<<" "<<attribute<<std::endl;
965
966 resetIteratorsMap();
967
968 // There might be cases where there no remote indices for that process yet
969 typename IteratorsMap::iterator found = iteratorsMap_.find(process);
970
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;
977 }
978
979 Iterators& iterators = found->second;
980
981 // Search for the remote index
982 while(iterators.isNotAtEnd() && iterators.globalIndexPair() < globalPair) {
983 // Increment all iterators
984 ++iterators;
985
986 }
987
988 if(iterators.isAtEnd() || iterators.globalIndexPair() != globalPair) {
989 // The entry is not yet known
990 // Insert in the list and do not change the first iterator.
991 iterators.insert(RemoteIndex(Attribute(attribute)),globalPair);
992 return;
993 }
994
995 // Global indices match
996 bool indexIsThere=false;
997 for(Iterators tmpIterators = iterators;
998 !tmpIterators.isAtEnd() && tmpIterators.globalIndexPair() == globalPair;
999 ++tmpIterators)
1000 //entry already exists with the same attribute
1001 if(tmpIterators.globalIndexPair().second == attribute) {
1002 indexIsThere=true;
1003 break;
1004 }
1005
1006 if(!indexIsThere)
1007 // The entry is not yet known
1008 // Insert in the list and do not change the first iterator.
1009 iterators.insert(RemoteIndex(Attribute(attribute)),globalPair);
1010 }
1011
1012 template<typename T>
1013 template<typename T1>
1014 void IndicesSyncer<T>::recvAndUnpack(T1& numberer, int hardSource, bool useHardSource)
1015 {
1016 const ParallelIndexSet& constIndexSet = indexSet_;
1017 auto iEnd = constIndexSet.end();
1018 auto index = constIndexSet.begin();
1019 int bpos = 0;
1020 int publish;
1021
1022 assert(checkReset());
1023
1024 MPI_Status status;
1025
1026 // We have to determine the message size and source before the receive
1027
1028 MPI_Probe(useHardSource ? hardSource : MPI_ANY_SOURCE, 345, remoteIndices_.communicator(), &status);
1029
1030 int source=status.MPI_SOURCE;
1031 int count;
1032 MPI_Get_count(&status, MPI_PACKED, &count);
1033
1034 Dune::dvverb<<rank_<<": Receiving message from "<< source<<" with "<<count<<" bytes"<<std::endl;
1035
1036 if(count>receiveBufferSize_) {
1037 receiveBufferSize_=count;
1038 delete[] receiveBuffer_;
1039 receiveBuffer_ = new char[receiveBufferSize_];
1040 }
1041
1042 MPI_Recv(receiveBuffer_, count, MPI_PACKED, source, 345, remoteIndices_.communicator(), &status);
1043
1044 // How many global entries were published?
1045 MPI_Unpack(receiveBuffer_, count, &bpos, &publish, 1, MPI_INT, remoteIndices_.communicator());
1046
1047 // Now unpack the remote indices and add them.
1048 while(publish>0) {
1049
1050 // Unpack information about the local index on the source process
1051 GlobalIndex global; // global index of the current entry
1052 char sourceAttribute; // Attribute on the source process
1053 int pairs;
1054
1055 MPI_Unpack(receiveBuffer_, count, &bpos, &global, 1, MPITraits<GlobalIndex>::getType(),
1056 remoteIndices_.communicator());
1057 MPI_Unpack(receiveBuffer_, count, &bpos, &sourceAttribute, 1, MPI_CHAR,
1058 remoteIndices_.communicator());
1059 MPI_Unpack(receiveBuffer_, count, &bpos, &pairs, 1, MPI_INT,
1060 remoteIndices_.communicator());
1061
1062 // Insert the entry on the remote process to our
1063 // remote index list
1064 SLList<std::pair<int,Attribute> > sourceAttributeList;
1065 sourceAttributeList.push_back(std::make_pair(source,Attribute(sourceAttribute)));
1066#ifndef NDEBUG
1067 bool foundSelf = false;
1068#endif
1069 Attribute myAttribute=Attribute();
1070
1071 // Unpack the remote indices
1072 for(; pairs>0; --pairs) {
1073 // Unpack the process id that knows the index
1074 int process;
1075 char attribute;
1076 MPI_Unpack(receiveBuffer_, count, &bpos, &process, 1, MPI_INT,
1077 remoteIndices_.communicator());
1078 // Unpack the attribute
1079 MPI_Unpack(receiveBuffer_, count, &bpos, &attribute, 1, MPI_CHAR,
1080 remoteIndices_.communicator());
1081
1082 if(process==rank_) {
1083#ifndef NDEBUG
1084 foundSelf=true;
1085#endif
1086 myAttribute=Attribute(attribute);
1087 // Now we know the local attribute of the global index
1088 //Only add the index if it is unknown.
1089 // Do we know that global index already?
1090 auto pos = std::lower_bound(index, iEnd, IndexPair(global));
1091
1092 if(pos == iEnd || pos->global() != global) {
1093 // no entry with this global index
1094 indexSet_.add(global,
1095 ParallelLocalIndex<Attribute>(numberer(global),
1096 myAttribute, true));
1097 Dune::dvverb << "Adding "<<global<<" "<<myAttribute<<std::endl;
1098 continue;
1099 }
1100
1101 // because of above the global indices match. Add only if the attribute is different
1102 bool indexIsThere = false;
1103 index=pos;
1104
1105 for(; pos->global()==global; ++pos)
1106 if(pos->local().attribute() == myAttribute) {
1107 Dune::dvverb<<"found "<<global<<" "<<myAttribute<<std::endl;
1108 indexIsThere = true;
1109 break;
1110 }
1111
1112 if(!indexIsThere) {
1113 indexSet_.add(global,
1114 ParallelLocalIndex<Attribute>(numberer(global),
1115 myAttribute, true));
1116 Dune::dvverb << "Adding "<<global<<" "<<myAttribute<<std::endl;
1117 }
1118
1119 }else{
1120 sourceAttributeList.push_back(std::make_pair(process,Attribute(attribute)));
1121 }
1122 }
1123 assert(foundSelf);
1124 // Insert remote indices
1125 typedef typename SLList<std::pair<int,Attribute> >::const_iterator Iter;
1126 for(Iter i=sourceAttributeList.begin(), end=sourceAttributeList.end();
1127 i!=end; ++i)
1128 insertIntoRemoteIndexList(i->first, std::make_pair(global, myAttribute),
1129 i->second);
1130 --publish;
1131 }
1132
1133 resetIteratorsMap();
1134 }
1135
1136 template<typename T>
1137 void IndicesSyncer<T>::resetIteratorsMap(){
1138
1139 // Reset iterators in all tuples.
1140 const auto remoteEnd = remoteIndices_.remoteIndices_.end();
1141 auto iterators = iteratorsMap_.begin();
1142 auto global = globalMap_.begin();
1143 auto added = oldMap_.begin();
1144
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);
1148 }
1149 }
1150
1151 template<typename T>
1152 bool IndicesSyncer<T>::checkReset(const Iterators& iterators, RemoteIndexList& rList, GlobalIndexList& gList,
1153 BoolList& bList){
1154
1155 if(std::get<0>(iterators.iterators_) != rList.begin())
1156 return false;
1157 if(std::get<1>(iterators.iterators_) != gList.begin())
1158 return false;
1159 if(std::get<2>(iterators.iterators_) != bList.begin())
1160 return false;
1161 return true;
1162 }
1163
1164
1165 template<typename T>
1166 bool IndicesSyncer<T>::checkReset(){
1167
1168 // Reset iterators in all tuples.
1169 const auto remoteEnd = remoteIndices_.remoteIndices_.end();
1170 auto iterators = iteratorsMap_.begin();
1171 auto global = globalMap_.begin();
1172 auto added = oldMap_.begin();
1173 bool ret = true;
1174
1175 for(auto remote = remoteIndices_.remoteIndices_.begin();
1176 remote != remoteEnd; ++remote, ++global, ++added, ++iterators) {
1177 if(!checkReset(iterators->second, *(remote->second.first), global->second,
1178 added->second))
1179 ret=false;
1180 }
1181 return ret;
1182 }
1183}
1184
1185#endif // HAVE_MPI
1186#endif // DUNE_COMMON_PARALLEL_INDICESSYNCER_HH
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
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
Provides a map between global and local indices.
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.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 27, 23:30, 2024)