Dune Core Modules (2.9.1)

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 (C) 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_INDICESSYNCER_HH
6#define DUNE_INDICESSYNCER_HH
7
8#include "indexset.hh"
9#include "remoteindices.hh"
11#include <dune/common/sllist.hh>
12#include <cassert>
13#include <cmath>
14#include <limits>
15#include <algorithm>
16#include <functional>
17#include <map>
18#include <tuple>
19
20#if HAVE_MPI
21namespace Dune
22{
39 template<typename T>
41 {
42 public:
43
46
49
52
54 typedef typename ParallelIndexSet::LocalIndex::Attribute Attribute;
55
60
71 RemoteIndices& remoteIndices);
72
82 void sync();
83
94 template<typename T1>
95 void sync(T1& numberer);
96
97 private:
98
100 ParallelIndexSet& indexSet_;
101
103 RemoteIndices& remoteIndices_;
104
106 char** sendBuffers_;
107
109 char* receiveBuffer_;
110
112 std::size_t* sendBufferSizes_;
113
115 int receiveBufferSize_; // int because of MPI
116
120 struct MessageInformation
121 {
122 MessageInformation()
123 : publish(), pairs()
124 {}
131 int pairs;
132 };
133
137 class DefaultNumberer
138 {
139 public:
145 std::size_t operator()([[maybe_unused]] const GlobalIndex& global)
146 {
148 }
149 };
150
152 MPI_Datatype datatype_;
153
155 int rank_;
156
161 typedef SLList<std::pair<GlobalIndex,Attribute>, typename RemoteIndices::Allocator> GlobalIndexList;
162
164 typedef typename GlobalIndexList::ModifyIterator GlobalIndexModifier;
165
170 GlobalIndexIterator;
171
173 typedef std::map<int, GlobalIndexList> GlobalIndicesMap;
174
183 GlobalIndicesMap globalMap_;
184
189
193 typedef typename BoolList::iterator BoolIterator;
194
196 typedef typename BoolList::ModifyIterator BoolListModifier;
197
199 typedef std::map<int,BoolList> BoolMap;
200
205 BoolMap oldMap_;
206
208 std::map<int,MessageInformation> infoSend_;
209
211 typedef typename RemoteIndices::RemoteIndexList RemoteIndexList;
212
214 typedef typename RemoteIndexList::ModifyIterator RemoteIndexModifier;
215
218
220 typedef typename RemoteIndexList::iterator RemoteIndexIterator;
221
223 typedef typename RemoteIndexList::const_iterator ConstRemoteIndexIterator;
224
226 typedef std::tuple<RemoteIndexModifier,GlobalIndexModifier,BoolListModifier,
227 const ConstRemoteIndexIterator> IteratorTuple;
228
236 class Iterators
237 {
238 friend class IndicesSyncer<T>;
239 public:
249 Iterators(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
250 BoolList& booleans);
251
255 Iterators();
256
260 Iterators& operator++();
261
267 void insert(const RemoteIndex& index,
268 const std::pair<GlobalIndex,Attribute>& global);
269
274 RemoteIndex& remoteIndex() const;
275
280 std::pair<GlobalIndex,Attribute>& globalIndexPair() const;
281
282 Attribute& attribute() const;
283
289 bool isOld() const;
290
300 void reset(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
301 BoolList& booleans);
302
308 bool isNotAtEnd() const;
309
315 bool isAtEnd() const;
316
317 private:
327 IteratorTuple iterators_;
328 };
329
331 typedef std::map<int,Iterators> IteratorsMap;
332
344 IteratorsMap iteratorsMap_;
345
347 void calculateMessageSizes();
348
356 void packAndSend(int destination, char* buffer, std::size_t bufferSize, MPI_Request& req);
357
362 template<typename T1>
363 void recvAndUnpack(T1& numberer);
364
368 void registerMessageDatatype();
369
373 void insertIntoRemoteIndexList(int process,
374 const std::pair<GlobalIndex,Attribute>& global,
375 char attribute);
376
380 void resetIteratorsMap();
381
386 bool checkReset();
387
396 bool checkReset(const Iterators& iterators, RemoteIndexList& rlist, GlobalIndexList& gList,
397 BoolList& bList);
398 };
399
400 template<typename TG, typename TA>
401 bool operator<(const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
402 const std::pair<TG,TA>& i2)
403 {
404 return i1.global() < i2.first ||
405 (i1.global() == i2.first && i1.local().attribute()<i2.second);
406 }
407
408 template<typename TG, typename TA>
409 bool operator<(const std::pair<TG,TA>& i1,
410 const IndexPair<TG,ParallelLocalIndex<TA> >& i2)
411 {
412 return i1.first < i2.global() ||
413 (i1.first == i2.global() && i1.second<i2.local().attribute());
414 }
415
416 template<typename TG, typename TA>
417 bool operator==(const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
418 const std::pair<TG,TA>& i2)
419 {
420 return (i1.global() == i2.first && i1.local().attribute()==i2.second);
421 }
422
423 template<typename TG, typename TA>
424 bool operator!=(const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
425 const std::pair<TG,TA>& i2)
426 {
427 return (i1.global() != i2.first || i1.local().attribute()!=i2.second);
428 }
429
430 template<typename TG, typename TA>
431 bool operator==(const std::pair<TG,TA>& i2,
432 const IndexPair<TG,ParallelLocalIndex<TA> >& i1)
433 {
434 return (i1.global() == i2.first && i1.local().attribute()==i2.second);
435 }
436
437 template<typename TG, typename TA>
438 bool operator!=(const std::pair<TG,TA>& i2,
439 const IndexPair<TG,ParallelLocalIndex<TA> >& i1)
440 {
441 return (i1.global() != i2.first || i1.local().attribute()!=i2.second);
442 }
443
460 template<typename T, typename A, typename A1>
461 void storeGlobalIndicesOfRemoteIndices(std::map<int,SLList<std::pair<typename T::GlobalIndex, typename T::LocalIndex::Attribute>,A> >& globalMap,
462 const RemoteIndices<T,A1>& remoteIndices)
463 {
464 for(auto remote = remoteIndices.begin(), end =remoteIndices.end(); remote != end; ++remote) {
465 typedef typename RemoteIndices<T,A1>::RemoteIndexList RemoteIndexList;
467 GlobalIndexList& global = globalMap[remote->first];
468 RemoteIndexList& rList = *(remote->second.first);
469
470 for(auto index = rList.begin(), riEnd = rList.end();
471 index != riEnd; ++index) {
472 global.push_back(std::make_pair(index->localIndexPair().global(),
473 index->localIndexPair().local().attribute()));
474 }
475 }
476 }
477
486 template<typename T, typename A, typename A1>
487 inline void repairLocalIndexPointers(std::map<int,
488 SLList<std::pair<typename T::GlobalIndex,
489 typename T::LocalIndex::Attribute>,A> >& globalMap,
490 RemoteIndices<T,A1>& remoteIndices,
491 const T& indexSet)
492 {
493 assert(globalMap.size()==static_cast<std::size_t>(remoteIndices.neighbours()));
494 // Repair pointers to index set in remote indices.
495 auto global = globalMap.begin();
496 auto end = remoteIndices.remoteIndices_.end();
497
498 for(auto remote = remoteIndices.remoteIndices_.begin(); remote != end; ++remote, ++global) {
499 assert(remote->first==global->first);
500 assert(remote->second.first->size() == global->second.size());
501
502 auto riEnd = remote->second.first->end();
503 auto rIndex = remote->second.first->begin();
504 auto gIndex = global->second.begin();
505 auto index = indexSet.begin();
506
507 assert(rIndex==riEnd || gIndex != global->second.end());
508 while(rIndex != riEnd) {
509 // Search for the index in the set.
510 assert(gIndex != global->second.end());
511
512 while(!(index->global() == gIndex->first
513 && index->local().attribute() == gIndex->second)) {
514 ++index;
515 // this is only needed for ALU, where there may exist
516 // more entries with the same global index in the remote index set
517 // than in the index set
518 if (index->global() > gIndex->first) {
519 index=indexSet.begin();
520 }
521 }
522
523 assert(index != indexSet.end() && *index == *gIndex);
524
525 rIndex->localIndex_ = &(*index);
526 ++index;
527 ++rIndex;
528 ++gIndex;
529 }
530 }
531 remoteIndices.sourceSeqNo_ = remoteIndices.source_->seqNo();
532 remoteIndices.destSeqNo_ = remoteIndices.target_->seqNo();
533 }
534
535 template<typename T>
537 RemoteIndices& remoteIndices)
538 : indexSet_(indexSet), remoteIndices_(remoteIndices)
539 {
540 // index sets must match.
541 assert(remoteIndices.source_ == remoteIndices.target_);
542 assert(remoteIndices.source_ == &indexSet);
543 MPI_Comm_rank(remoteIndices_.communicator(), &rank_);
544 }
545
546 template<typename T>
548 GlobalIndexList& globalIndices,
549 BoolList& booleans)
550 : iterators_(remoteIndices.beginModify(), globalIndices.beginModify(),
551 booleans.beginModify(), remoteIndices.end())
552 { }
553
554 template<typename T>
556 : iterators_()
557 {}
558
559 template<typename T>
561 {
562 ++(std::get<0>(iterators_));
563 ++(std::get<1>(iterators_));
564 ++(std::get<2>(iterators_));
565 return *this;
566 }
567
568 template<typename T>
570 const std::pair<GlobalIndex,Attribute>& global)
571 {
572 std::get<0>(iterators_).insert(index);
573 std::get<1>(iterators_).insert(global);
574 std::get<2>(iterators_).insert(false);
575 }
576
577 template<typename T>
578 inline typename IndicesSyncer<T>::RemoteIndex&
580 {
581 return *(std::get<0>(iterators_));
582 }
583
584 template<typename T>
585 inline std::pair<typename IndicesSyncer<T>::GlobalIndex,typename IndicesSyncer<T>::Attribute>&
587 {
588 return *(std::get<1>(iterators_));
589 }
590
591 template<typename T>
593 {
594 return *(std::get<2>(iterators_));
595 }
596
597 template<typename T>
599 GlobalIndexList& globalIndices,
600 BoolList& booleans)
601 {
602 std::get<0>(iterators_) = remoteIndices.beginModify();
603 std::get<1>(iterators_) = globalIndices.beginModify();
604 std::get<2>(iterators_) = booleans.beginModify();
605 }
606
607 template<typename T>
609 {
610 return std::get<0>(iterators_) != std::get<3>(iterators_);
611 }
612
613 template<typename T>
615 {
616 return std::get<0>(iterators_) == std::get<3>(iterators_);
617 }
618
619 template<typename T>
621 {
622 MPI_Datatype type[2] = {MPI_INT, MPI_INT};
623 int blocklength[2] = {1,1};
624 MPI_Aint displacement[2];
625 MPI_Aint base;
626
627 // Compute displacement
628 MessageInformation message;
629
630 MPI_Get_address( &(message.publish), displacement);
631 MPI_Get_address( &(message.pairs), displacement+1);
632
633 // Make the displacement relative
634 MPI_Get_address(&message, &base);
635 displacement[0] -= base;
636 displacement[1] -= base;
637
638 MPI_Type_create_struct( 2, blocklength, displacement, type, &datatype_);
639 MPI_Type_commit(&datatype_);
640 }
641
642 template<typename T>
643 void IndicesSyncer<T>::calculateMessageSizes()
644 {
645 auto iEnd = indexSet_.end();
646 auto collIter = remoteIndices_.template iterator<true>();
647
648 for(auto index = indexSet_.begin(); index != iEnd; ++index) {
649 collIter.advance(index->global(), index->local().attribute());
650 if(collIter.empty())
651 break;
652 int knownRemote=0;
653 auto end = collIter.end();
654
655 // Count the remote indices we know.
656 for(auto valid = collIter.begin(); valid != end; ++valid) {
657 ++knownRemote;
658 }
659
660 if(knownRemote>0) {
661 Dune::dverb<<rank_<<": publishing "<<knownRemote<<" for index "<<index->global()<< " for processes ";
662
663 // Update MessageInformation
664 for(auto valid = collIter.begin(); valid != end; ++valid) {
665 ++(infoSend_[valid.process()].publish);
666 (infoSend_[valid.process()].pairs) += knownRemote;
667 Dune::dverb<<valid.process()<<" ";
668 Dune::dverb<<"(publish="<<infoSend_[valid.process()].publish<<", pairs="<<infoSend_[valid.process()].pairs
669 <<") ";
670 }
671 Dune::dverb<<std::endl;
672 }
673 }
674
675 const auto end = infoSend_.end();
676
677 // Now determine the buffersizes needed for each neighbour using MPI_Pack_size
678 MessageInformation dummy;
679
680 auto messageIter= infoSend_.begin();
681 const auto rend = remoteIndices_.end();
682 int neighbour=0;
683
684 for(auto remote = remoteIndices_.begin(); remote != rend; ++remote, ++neighbour) {
685 MessageInformation* message;
686 MessageInformation recv;
687
688 if(messageIter != end && messageIter->first==remote->first) {
689 // We want to send message information to that process
690 message = const_cast<MessageInformation*>(&(messageIter->second));
691 ++messageIter;
692 }else
693 // We do not want to send information but the other process might.
694 message = &dummy;
695
696 sendBufferSizes_[neighbour]=0;
697 int tsize;
698 // The number of indices published
699 MPI_Pack_size(1, MPI_INT,remoteIndices_.communicator(), &tsize);
700 sendBufferSizes_[neighbour] += tsize;
701
702 for(int i=0; i < message->publish; ++i) {
703 // The global index
704 MPI_Pack_size(1, MPITraits<GlobalIndex>::getType(), remoteIndices_.communicator(), &tsize);
705 sendBufferSizes_[neighbour] += tsize;
706 // The attribute in the local index
707 MPI_Pack_size(1, MPI_CHAR, remoteIndices_.communicator(), &tsize);
708 sendBufferSizes_[neighbour] += tsize;
709 // The number of corresponding remote indices
710 MPI_Pack_size(1, MPI_INT, remoteIndices_.communicator(), &tsize);
711 sendBufferSizes_[neighbour] += tsize;
712 }
713 for(int i=0; i < message->pairs; ++i) {
714 // The process of the remote index
715 MPI_Pack_size(1, MPI_INT, remoteIndices_.communicator(), &tsize);
716 sendBufferSizes_[neighbour] += tsize;
717 // The attribute of the remote index
718 MPI_Pack_size(1, MPI_CHAR, remoteIndices_.communicator(), &tsize);
719 sendBufferSizes_[neighbour] += tsize;
720 }
721
722 Dune::dverb<<rank_<<": Buffer (neighbour="<<remote->first<<") size is "<< sendBufferSizes_[neighbour]<<" for publish="<<message->publish<<" pairs="<<message->pairs<<std::endl;
723 }
724
725 }
726
727 template<typename T>
729 {
730 DefaultNumberer numberer;
731 sync(numberer);
732 }
733
734 template<typename T>
735 template<typename T1>
736 void IndicesSyncer<T>::sync(T1& numberer)
737 {
738 // The pointers to the local indices in the remote indices
739 // will become invalid due to the resorting of the index set.
740 // Therefore store the corresponding global indices.
741 // Mark all indices as not added
742 const auto end = remoteIndices_.end();
743
744 // Number of neighbours might change during the syncing.
745 // save the old neighbours
746 std::size_t noOldNeighbours = remoteIndices_.neighbours();
747 int* oldNeighbours = new int[noOldNeighbours];
748 sendBufferSizes_ = new std::size_t[noOldNeighbours];
749 std::size_t neighbourI = 0;
750
751 for(auto remote = remoteIndices_.begin(); remote != end; ++remote, ++neighbourI) {
752 oldNeighbours[neighbourI] = remote->first;
753
754 // Make sure we only have one remote index list.
755 assert(remote->second.first==remote->second.second);
756
757 RemoteIndexList& rList = *(remote->second.first);
758
759 // Store the corresponding global indices.
760 GlobalIndexList& global = globalMap_[remote->first];
761 BoolList& added = oldMap_[remote->first];
762 auto riEnd = rList.end();
763
764 for(auto index = rList.begin();
765 index != riEnd; ++index) {
766 global.push_back(std::make_pair(index->localIndexPair().global(),
767 index->localIndexPair().local().attribute()));
768 added.push_back(true);
769 }
770
771 Iterators iterators(rList, global, added);
772 iteratorsMap_.insert(std::make_pair(remote->first, iterators));
773 assert(checkReset(iteratorsMap_[remote->first], rList,global,added));
774 }
775
776 // Exchange indices with each neighbour
777 calculateMessageSizes();
778
779 // Allocate the buffers
780 receiveBufferSize_=1;
781 sendBuffers_ = new char*[noOldNeighbours];
782
783 for(std::size_t i=0; i<noOldNeighbours; ++i) {
784 sendBuffers_[i] = new char[sendBufferSizes_[i]];
785 receiveBufferSize_ = std::max(receiveBufferSize_, static_cast<int>(sendBufferSizes_[i]));
786 }
787
788 receiveBuffer_=new char[receiveBufferSize_];
789
790 indexSet_.beginResize();
791
792 Dune::dverb<<rank_<<": Neighbours: ";
793
794 for(std::size_t i = 0; i<noOldNeighbours; ++i)
795 Dune::dverb<<oldNeighbours[i]<<" ";
796
797 Dune::dverb<<std::endl;
798
799 MPI_Request* requests = new MPI_Request[noOldNeighbours];
800 MPI_Status* statuses = new MPI_Status[noOldNeighbours];
801
802 // Pack Message data and start the sends
803 for(std::size_t i = 0; i<noOldNeighbours; ++i)
804 packAndSend(oldNeighbours[i], sendBuffers_[i], sendBufferSizes_[i], requests[i]);
805
806 // Probe for incoming messages, receive and unpack them
807 for(std::size_t i = 0; i<noOldNeighbours; ++i)
808 recvAndUnpack(numberer);
809 // }else{
810 // recvAndUnpack(oldNeighbours[i], numberer);
811 // packAndSend(oldNeighbours[i]);
812 // }
813 // }
814
815 delete[] receiveBuffer_;
816
817 // Wait for the completion of the sends
818 // Wait for completion of sends
819 if(MPI_SUCCESS!=MPI_Waitall(noOldNeighbours, requests, statuses)) {
820 std::cerr<<": MPI_Error occurred while sending message"<<std::endl;
821 for(std::size_t i=0; i< noOldNeighbours; i++)
822 if(MPI_SUCCESS!=statuses[i].MPI_ERROR)
823 std::cerr<<"Destination "<<statuses[i].MPI_SOURCE<<" error code: "<<statuses[i].MPI_ERROR<<std::endl;
824 }
825
826 delete[] statuses;
827 delete[] requests;
828
829 for(std::size_t i=0; i<noOldNeighbours; ++i)
830 delete[] sendBuffers_[i];
831
832 delete[] sendBuffers_;
833 delete[] sendBufferSizes_;
834
835 // No need for the iterator tuples any more
836 iteratorsMap_.clear();
837
838 indexSet_.endResize();
839
840 delete[] oldNeighbours;
841
842 repairLocalIndexPointers(globalMap_, remoteIndices_, indexSet_);
843
844 oldMap_.clear();
845 globalMap_.clear();
846
847 // update the sequence number
848 remoteIndices_.sourceSeqNo_ = remoteIndices_.destSeqNo_ = indexSet_.seqNo();
849 }
850
851 template<typename T>
852 void IndicesSyncer<T>::packAndSend(int destination, char* buffer, std::size_t bufferSize, MPI_Request& request)
853 {
854 auto iEnd = indexSet_.end();
855 int bpos = 0;
856 int published = 0;
857 int pairs = 0;
858
859 assert(checkReset());
860
861 // Pack the number of indices we publish
862 MPI_Pack(&(infoSend_[destination].publish), 1, MPI_INT, buffer, bufferSize, &bpos,
863 remoteIndices_.communicator());
864
865 for(auto index = indexSet_.begin(); index != iEnd; ++index) {
866 // Search for corresponding remote indices in all iterator tuples
867 auto iteratorsEnd = iteratorsMap_.end();
868
869 // advance all iterators to a position with global index >= index->global()
870 for(auto iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators) {
871 while(iterators->second.isNotAtEnd() &&
872 iterators->second.globalIndexPair().first < index->global())
873 ++(iterators->second);
874 assert(!iterators->second.isNotAtEnd() || iterators->second.globalIndexPair().first >= index->global());
875 }
876
877 // Add all remote indices positioned at global which were already present before calling sync
878 // to the message.
879 // Count how many remote indices we will send
880 int indices = 0;
881 bool knownRemote = false; // Is the remote process supposed to know this index?
882
883 for(auto iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
884 {
885 std::pair<GlobalIndex,Attribute> p;
886 if (iterators->second.isNotAtEnd())
887 {
888 p = iterators->second.globalIndexPair();
889 }
890
891 if(iterators->second.isNotAtEnd() && iterators->second.isOld()
892 && iterators->second.globalIndexPair().first == index->global()) {
893 indices++;
894 if(destination == iterators->first)
895 knownRemote = true;
896 }
897 }
898
899 if(!knownRemote)
900 // We do not need to send any indices
901 continue;
902
903 Dune::dverb<<rank_<<": sending "<<indices<<" for index "<<index->global()<<" to "<<destination<<std::endl;
904
905
906 // Pack the global index, the attribute and the number
907 MPI_Pack(const_cast<GlobalIndex*>(&(index->global())), 1, MPITraits<GlobalIndex>::getType(), buffer, bufferSize, &bpos,
908 remoteIndices_.communicator());
909
910 char attr = index->local().attribute();
911 MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos,
912 remoteIndices_.communicator());
913
914 // Pack the number of remote indices we send.
915 MPI_Pack(&indices, 1, MPI_INT, buffer, bufferSize, &bpos,
916 remoteIndices_.communicator());
917
918 // Pack the information about the remote indices
919 for(auto iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
920 if(iterators->second.isNotAtEnd() && iterators->second.isOld()
921 && iterators->second.globalIndexPair().first == index->global()) {
922 int process = iterators->first;
923
924 ++pairs;
925 assert(pairs <= infoSend_[destination].pairs);
926 MPI_Pack(&process, 1, MPI_INT, buffer, bufferSize, &bpos,
927 remoteIndices_.communicator());
928 char attr2 = iterators->second.remoteIndex().attribute();
929
930 MPI_Pack(&attr2, 1, MPI_CHAR, buffer, bufferSize, &bpos,
931 remoteIndices_.communicator());
932 --indices;
933 }
934 assert(indices==0);
935 ++published;
936 Dune::dvverb<<" (publish="<<published<<", pairs="<<pairs<<")"<<std::endl;
937 assert(published <= infoSend_[destination].publish);
938 }
939
940 // Make sure we send all expected entries
941 assert(published == infoSend_[destination].publish);
942 assert(pairs == infoSend_[destination].pairs);
943 resetIteratorsMap();
944
945 Dune::dverb << rank_<<": Sending message of "<<bpos<<" bytes to "<<destination<<std::endl;
946
947 MPI_Issend(buffer, bpos, MPI_PACKED, destination, 345, remoteIndices_.communicator(),&request);
948 }
949
950 template<typename T>
951 inline void IndicesSyncer<T>::insertIntoRemoteIndexList(int process,
952 const std::pair<GlobalIndex,Attribute>& globalPair,
953 char attribute)
954 {
955 Dune::dverb<<"Inserting from "<<process<<" "<<globalPair.first<<", "<<
956 globalPair.second<<" "<<attribute<<std::endl;
957
958 resetIteratorsMap();
959
960 // There might be cases where there no remote indices for that process yet
961 typename IteratorsMap::iterator found = iteratorsMap_.find(process);
962
963 if( found == iteratorsMap_.end() ) {
964 Dune::dverb<<"Discovered new neighbour "<<process<<std::endl;
965 RemoteIndexList* rlist = new RemoteIndexList();
966 remoteIndices_.remoteIndices_.insert(std::make_pair(process,std::make_pair(rlist,rlist)));
967 Iterators iterators = Iterators(*rlist, globalMap_[process], oldMap_[process]);
968 found = iteratorsMap_.insert(std::make_pair(process, iterators)).first;
969 }
970
971 Iterators& iterators = found->second;
972
973 // Search for the remote index
974 while(iterators.isNotAtEnd() && iterators.globalIndexPair() < globalPair) {
975 // Increment all iterators
976 ++iterators;
977
978 }
979
980 if(iterators.isAtEnd() || iterators.globalIndexPair() != globalPair) {
981 // The entry is not yet known
982 // Insert in the list and do not change the first iterator.
983 iterators.insert(RemoteIndex(Attribute(attribute)),globalPair);
984 return;
985 }
986
987 // Global indices match
988 bool indexIsThere=false;
989 for(Iterators tmpIterators = iterators;
990 !tmpIterators.isAtEnd() && tmpIterators.globalIndexPair() == globalPair;
991 ++tmpIterators)
992 //entry already exists with the same attribute
993 if(tmpIterators.globalIndexPair().second == attribute) {
994 indexIsThere=true;
995 break;
996 }
997
998 if(!indexIsThere)
999 // The entry is not yet known
1000 // Insert in the list and do not change the first iterator.
1001 iterators.insert(RemoteIndex(Attribute(attribute)),globalPair);
1002 }
1003
1004 template<typename T>
1005 template<typename T1>
1006 void IndicesSyncer<T>::recvAndUnpack(T1& numberer)
1007 {
1008 const ParallelIndexSet& constIndexSet = indexSet_;
1009 auto iEnd = constIndexSet.end();
1010 auto index = constIndexSet.begin();
1011 int bpos = 0;
1012 int publish;
1013
1014 assert(checkReset());
1015
1016 MPI_Status status;
1017
1018 // We have to determine the message size and source before the receive
1019 MPI_Probe(MPI_ANY_SOURCE, 345, remoteIndices_.communicator(), &status);
1020
1021 int source=status.MPI_SOURCE;
1022 int count;
1023 MPI_Get_count(&status, MPI_PACKED, &count);
1024
1025 Dune::dvverb<<rank_<<": Receiving message from "<< source<<" with "<<count<<" bytes"<<std::endl;
1026
1027 if(count>receiveBufferSize_) {
1028 receiveBufferSize_=count;
1029 delete[] receiveBuffer_;
1030 receiveBuffer_ = new char[receiveBufferSize_];
1031 }
1032
1033 MPI_Recv(receiveBuffer_, count, MPI_PACKED, source, 345, remoteIndices_.communicator(), &status);
1034
1035 // How many global entries were published?
1036 MPI_Unpack(receiveBuffer_, count, &bpos, &publish, 1, MPI_INT, remoteIndices_.communicator());
1037
1038 // Now unpack the remote indices and add them.
1039 while(publish>0) {
1040
1041 // Unpack information about the local index on the source process
1042 GlobalIndex global; // global index of the current entry
1043 char sourceAttribute; // Attribute on the source process
1044 int pairs;
1045
1046 MPI_Unpack(receiveBuffer_, count, &bpos, &global, 1, MPITraits<GlobalIndex>::getType(),
1047 remoteIndices_.communicator());
1048 MPI_Unpack(receiveBuffer_, count, &bpos, &sourceAttribute, 1, MPI_CHAR,
1049 remoteIndices_.communicator());
1050 MPI_Unpack(receiveBuffer_, count, &bpos, &pairs, 1, MPI_INT,
1051 remoteIndices_.communicator());
1052
1053 // Insert the entry on the remote process to our
1054 // remote index list
1055 SLList<std::pair<int,Attribute> > sourceAttributeList;
1056 sourceAttributeList.push_back(std::make_pair(source,Attribute(sourceAttribute)));
1057#ifndef NDEBUG
1058 bool foundSelf = false;
1059#endif
1060 Attribute myAttribute=Attribute();
1061
1062 // Unpack the remote indices
1063 for(; pairs>0; --pairs) {
1064 // Unpack the process id that knows the index
1065 int process;
1066 char attribute;
1067 MPI_Unpack(receiveBuffer_, count, &bpos, &process, 1, MPI_INT,
1068 remoteIndices_.communicator());
1069 // Unpack the attribute
1070 MPI_Unpack(receiveBuffer_, count, &bpos, &attribute, 1, MPI_CHAR,
1071 remoteIndices_.communicator());
1072
1073 if(process==rank_) {
1074#ifndef NDEBUG
1075 foundSelf=true;
1076#endif
1077 myAttribute=Attribute(attribute);
1078 // Now we know the local attribute of the global index
1079 //Only add the index if it is unknown.
1080 // Do we know that global index already?
1081 auto pos = std::lower_bound(index, iEnd, IndexPair(global));
1082
1083 if(pos == iEnd || pos->global() != global) {
1084 // no entry with this global index
1085 indexSet_.add(global,
1086 ParallelLocalIndex<Attribute>(numberer(global),
1087 myAttribute, true));
1088 Dune::dvverb << "Adding "<<global<<" "<<myAttribute<<std::endl;
1089 continue;
1090 }
1091
1092 // because of above the global indices match. Add only if the attribute is different
1093 bool indexIsThere = false;
1094 index=pos;
1095
1096 for(; pos->global()==global; ++pos)
1097 if(pos->local().attribute() == myAttribute) {
1098 Dune::dvverb<<"found "<<global<<" "<<myAttribute<<std::endl;
1099 indexIsThere = true;
1100 break;
1101 }
1102
1103 if(!indexIsThere) {
1104 indexSet_.add(global,
1105 ParallelLocalIndex<Attribute>(numberer(global),
1106 myAttribute, true));
1107 Dune::dvverb << "Adding "<<global<<" "<<myAttribute<<std::endl;
1108 }
1109
1110 }else{
1111 sourceAttributeList.push_back(std::make_pair(process,Attribute(attribute)));
1112 }
1113 }
1114 assert(foundSelf);
1115 // Insert remote indices
1116 typedef typename SLList<std::pair<int,Attribute> >::const_iterator Iter;
1117 for(Iter i=sourceAttributeList.begin(), end=sourceAttributeList.end();
1118 i!=end; ++i)
1119 insertIntoRemoteIndexList(i->first, std::make_pair(global, myAttribute),
1120 i->second);
1121 --publish;
1122 }
1123
1124 resetIteratorsMap();
1125 }
1126
1127 template<typename T>
1128 void IndicesSyncer<T>::resetIteratorsMap(){
1129
1130 // Reset iterators in all tuples.
1131 const auto remoteEnd = remoteIndices_.remoteIndices_.end();
1132 auto iterators = iteratorsMap_.begin();
1133 auto global = globalMap_.begin();
1134 auto added = oldMap_.begin();
1135
1136 for(auto remote = remoteIndices_.remoteIndices_.begin();
1137 remote != remoteEnd; ++remote, ++global, ++added, ++iterators) {
1138 iterators->second.reset(*(remote->second.first), global->second, added->second);
1139 }
1140 }
1141
1142 template<typename T>
1143 bool IndicesSyncer<T>::checkReset(const Iterators& iterators, RemoteIndexList& rList, GlobalIndexList& gList,
1144 BoolList& bList){
1145
1146 if(std::get<0>(iterators.iterators_) != rList.begin())
1147 return false;
1148 if(std::get<1>(iterators.iterators_) != gList.begin())
1149 return false;
1150 if(std::get<2>(iterators.iterators_) != bList.begin())
1151 return false;
1152 return true;
1153 }
1154
1155
1156 template<typename T>
1157 bool IndicesSyncer<T>::checkReset(){
1158
1159 // Reset iterators in all tuples.
1160 const auto remoteEnd = remoteIndices_.remoteIndices_.end();
1161 auto iterators = iteratorsMap_.begin();
1162 auto global = globalMap_.begin();
1163 auto added = oldMap_.begin();
1164 bool ret = true;
1165
1166 for(auto remote = remoteIndices_.remoteIndices_.begin();
1167 remote != remoteEnd; ++remote, ++global, ++added, ++iterators) {
1168 if(!checkReset(iterators->second, *(remote->second.first), global->second,
1169 added->second))
1170 ret=false;
1171 }
1172 return ret;
1173 }
1174}
1175
1176#endif
1177#endif
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:73
The indices present on remote processes.
Definition: remoteindices.hh:189
MPI_Comm communicator() const
Get the mpi communicator used.
Definition: remoteindices.hh:1696
const_iterator end() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1529
int neighbours() const
Get the number of processors we share indices with.
Definition: remoteindices.hh:1446
typename std::allocator_traits< A >::template rebind_alloc< RemoteIndex > Allocator
The type of the allocator for the remote index list.
Definition: remoteindices.hh:237
const_iterator begin() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1522
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:59
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:487
ParallelIndexSet::GlobalIndex GlobalIndex
Type of the global index used in the index set.
Definition: indicessyncer.hh:51
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:536
std::size_t operator()(const GlobalIndex &global)
Provide the local 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 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:461
void sync()
Sync the index set.
Definition: indicessyncer.hh:728
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: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
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:637
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:237
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:259
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
DVVerbType dvverb(std::cout)
stream for very verbose output.
Definition: stdstreams.hh:95
DVerbType dverb(std::cout)
Singleton of verbose debug stream.
Definition: stdstreams.hh:116
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 (Nov 21, 23:30, 2024)