Dune Core Modules (2.3.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// $Id$
4#ifndef DUNE_INDICESSYNCER_HH
5#define DUNE_INDICESSYNCER_HH
6
7#include "indexset.hh"
8#include "remoteindices.hh"
10#include <dune/common/tuples.hh>
11#include <dune/common/sllist.hh>
12#include <dune/common/unused.hh>
13#include <cassert>
14#include <cmath>
15#include <limits>
16#include <algorithm>
17#include <functional>
18#include <map>
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()(const GlobalIndex& global)
146 {
147 DUNE_UNUSED_PARAMETER(global);
148 return std::numeric_limits<size_t>::max();
149 }
150 };
151
153 MPI_Datatype datatype_;
154
156 int rank_;
157
162 typedef SLList<std::pair<GlobalIndex,Attribute>, typename RemoteIndices::Allocator> GlobalIndexList;
163
165 typedef typename GlobalIndexList::ModifyIterator GlobalIndexModifier;
166
171 GlobalIndexIterator;
172
174 typedef std::map<int, GlobalIndexList> GlobalIndicesMap;
175
184 GlobalIndicesMap globalMap_;
185
190
194 typedef typename BoolList::iterator BoolIterator;
195
197 typedef typename BoolList::ModifyIterator BoolListModifier;
198
200 typedef std::map<int,BoolList> BoolMap;
201
206 BoolMap oldMap_;
207
209 std::map<int,MessageInformation> infoSend_;
210
212 typedef typename RemoteIndices::RemoteIndexList RemoteIndexList;
213
215 typedef typename RemoteIndexList::ModifyIterator RemoteIndexModifier;
216
219
221 typedef typename RemoteIndexList::iterator RemoteIndexIterator;
222
224 typedef typename RemoteIndexList::const_iterator ConstRemoteIndexIterator;
225
227 typedef Dune::tuple<RemoteIndexModifier,GlobalIndexModifier,BoolListModifier,
228 const ConstRemoteIndexIterator> IteratorTuple;
229
237 class Iterators
238 {
239 friend class IndicesSyncer<T>;
240 public:
250 Iterators(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
251 BoolList& booleans);
252
256 Iterators();
257
261 Iterators& operator++();
262
268 void insert(const RemoteIndex& index,
269 const std::pair<GlobalIndex,Attribute>& global);
270
275 RemoteIndex& remoteIndex() const;
276
281 std::pair<GlobalIndex,Attribute>& globalIndexPair() const;
282
283 Attribute& attribute() const;
284
290 bool isOld() const;
291
301 void reset(RemoteIndexList& remoteIndices, GlobalIndexList& globalIndices,
302 BoolList& booleans);
303
309 bool isNotAtEnd() const;
310
316 bool isAtEnd() const;
317
318 private:
328 IteratorTuple iterators_;
329 };
330
332 typedef std::map<int,Iterators> IteratorsMap;
333
345 IteratorsMap iteratorsMap_;
346
348 void calculateMessageSizes();
349
357 void packAndSend(int destination, char* buffer, std::size_t bufferSize, MPI_Request& req);
358
363 template<typename T1>
364 void recvAndUnpack(T1& numberer);
365
369 void registerMessageDatatype();
370
374 void insertIntoRemoteIndexList(int process,
375 const std::pair<GlobalIndex,Attribute>& global,
376 char attribute);
377
381 void resetIteratorsMap();
382
387 bool checkReset();
388
397 bool checkReset(const Iterators& iterators, RemoteIndexList& rlist, GlobalIndexList& gList,
398 BoolList& bList);
399 };
400
401 template<typename TG, typename TA>
402 bool operator<(const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
403 const std::pair<TG,TA>& i2)
404 {
405 return i1.global() < i2.first ||
406 (i1.global() == i2.first && i1.local().attribute()<i2.second);
407 }
408
409 template<typename TG, typename TA>
410 bool operator<(const std::pair<TG,TA>& i1,
411 const IndexPair<TG,ParallelLocalIndex<TA> >& i2)
412 {
413 return i1.first < i2.global() ||
414 (i1.first == i2.global() && i1.second<i2.local().attribute());
415 }
416
417 template<typename TG, typename TA>
418 bool operator==(const IndexPair<TG,ParallelLocalIndex<TA> >& i1,
419 const std::pair<TG,TA>& i2)
420 {
421 return (i1.global() == i2.first && i1.local().attribute()==i2.second);
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 std::pair<TG,TA>& i2,
433 const IndexPair<TG,ParallelLocalIndex<TA> >& i1)
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
461 template<typename T, typename A, typename A1>
462 void storeGlobalIndicesOfRemoteIndices(std::map<int,SLList<std::pair<typename T::GlobalIndex, typename T::LocalIndex::Attribute>,A> >& globalMap,
463 const RemoteIndices<T,A1>& remoteIndices,
464 const T& indexSet)
465 {
466 typedef typename RemoteIndices<T,A1>::const_iterator RemoteIterator;
467
468 for(RemoteIterator remote = remoteIndices.begin(), end =remoteIndices.end(); remote != end; ++remote) {
469 typedef typename RemoteIndices<T,A1>::RemoteIndexList RemoteIndexList;
470 typedef typename RemoteIndexList::const_iterator RemoteIndexIterator;
472 GlobalIndexList& global = globalMap[remote->first];
473 RemoteIndexList& rList = *(remote->second.first);
474
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()));
479 }
480 }
481 }
482
491 template<typename T, typename A, typename A1>
492 inline void repairLocalIndexPointers(std::map<int,
493 SLList<std::pair<typename T::GlobalIndex,
494 typename T::LocalIndex::Attribute>,A> >& globalMap,
495 RemoteIndices<T,A1>& remoteIndices,
496 const T& indexSet)
497 {
498 typedef typename RemoteIndices<T,A1>::RemoteIndexMap::iterator RemoteIterator;
499 typedef typename RemoteIndices<T,A1>::RemoteIndexList::iterator RemoteIndexIterator;
500 typedef typename T::GlobalIndex GlobalIndex;
501 typedef typename T::LocalIndex::Attribute Attribute;
502 typedef std::pair<GlobalIndex,Attribute> GlobalIndexPair;
503 typedef SLList<GlobalIndexPair,A> GlobalIndexPairList;
504 typedef typename GlobalIndexPairList::iterator GlobalIndexIterator;
505
506 assert(globalMap.size()==static_cast<std::size_t>(remoteIndices.neighbours()));
507 // Repair pointers to index set in remote indices.
508 typename std::map<int,GlobalIndexPairList>::iterator global = globalMap.begin();
509 RemoteIterator end = remoteIndices.remoteIndices_.end();
510
511 for(RemoteIterator remote = remoteIndices.remoteIndices_.begin(); remote != end; ++remote, ++global) {
512 typedef typename T::const_iterator IndexIterator;
513
514 assert(remote->first==global->first);
515 assert(remote->second.first->size() == global->second.size());
516
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();
521
522 assert(rIndex==riEnd || gIndex != global->second.end());
523 while(rIndex != riEnd) {
524 // Search for the index in the set.
525 assert(gIndex != global->second.end());
526
527 while(!(index->global() == gIndex->first
528 && index->local().attribute() == gIndex->second)) {
529 ++index;
530 // this is only needed for ALU, where there may exist
531 // more entries with the same global index in the remote index set
532 // than in the index set
533 if (index->global() > gIndex->first) {
534 index=indexSet.begin();
535 }
536 }
537
538 assert(index != indexSet.end() && *index == *gIndex);
539
540 rIndex->localIndex_ = &(*index);
541 ++index;
542 ++rIndex;
543 ++gIndex;
544 }
545 }
546 remoteIndices.sourceSeqNo_ = remoteIndices.source_->seqNo();
547 remoteIndices.destSeqNo_ = remoteIndices.target_->seqNo();
548 }
549
550 template<typename T>
552 RemoteIndices& remoteIndices)
553 : indexSet_(indexSet), remoteIndices_(remoteIndices)
554 {
555 // index sets must match.
556 assert(remoteIndices.source_ == remoteIndices.target_);
557 assert(remoteIndices.source_ == &indexSet);
558 MPI_Comm_rank(remoteIndices_.communicator(), &rank_);
559 }
560
561 template<typename T>
563 GlobalIndexList& globalIndices,
564 BoolList& booleans)
565 : iterators_(remoteIndices.beginModify(), globalIndices.beginModify(),
566 booleans.beginModify(), remoteIndices.end())
567 { }
568
569 template<typename T>
571 : iterators_()
572 {}
573
574 template<typename T>
576 {
577 ++(get<0>(iterators_));
578 ++(get<1>(iterators_));
579 ++(get<2>(iterators_));
580 return *this;
581 }
582
583 template<typename T>
585 const std::pair<GlobalIndex,Attribute>& global)
586 {
587 get<0>(iterators_).insert(index);
588 get<1>(iterators_).insert(global);
589 get<2>(iterators_).insert(false);
590 }
591
592 template<typename T>
593 inline typename IndicesSyncer<T>::RemoteIndex&
595 {
596 return *(get<0>(iterators_));
597 }
598
599 template<typename T>
600 inline std::pair<typename IndicesSyncer<T>::GlobalIndex,typename IndicesSyncer<T>::Attribute>&
602 {
603 return *(get<1>(iterators_));
604 }
605
606 template<typename T>
608 {
609 return *(get<2>(iterators_));
610 }
611
612 template<typename T>
614 GlobalIndexList& globalIndices,
615 BoolList& booleans)
616 {
617 get<0>(iterators_) = remoteIndices.beginModify();
618 get<1>(iterators_) = globalIndices.beginModify();
619 get<2>(iterators_) = booleans.beginModify();
620 }
621
622 template<typename T>
624 {
625 return get<0>(iterators_)!=get<3>(iterators_);
626 }
627
628 template<typename T>
630 {
631 return get<0>(iterators_)==get<3>(iterators_);
632 }
633
634 template<typename T>
636 {
637 MPI_Datatype type[2] = {MPI_INT, MPI_INT};
638 int blocklength[2] = {1,1};
639 MPI_Aint displacement[2];
640 MPI_Aint base;
641
642 // Compute displacement
643 MessageInformation message;
644
645 MPI_Address( &(message.publish), displacement);
646 MPI_Address( &(message.pairs), displacement+1);
647
648 // Make the displacement relative
649 MPI_Address(&message, &base);
650 displacement[0] -= base;
651 displacement[1] -= base;
652
653 MPI_Type_struct( 2, blocklength, displacement, type, &datatype_);
654 MPI_Type_commit(&datatype_);
655 }
656
657 template<typename T>
658 void IndicesSyncer<T>::calculateMessageSizes()
659 {
660 typedef typename ParallelIndexSet::const_iterator IndexIterator;
661 typedef CollectiveIterator<T,typename RemoteIndices::Allocator> CollectiveIterator;
662
663 IndexIterator iEnd = indexSet_.end();
664 CollectiveIterator collIter = remoteIndices_.template iterator<true>();
665
666 for(IndexIterator index = indexSet_.begin(); index != iEnd; ++index) {
667 collIter.advance(index->global(), index->local().attribute());
668 if(collIter.empty())
669 break;
670 int knownRemote=0;
671
672 typedef typename CollectiveIterator::iterator ValidIterator;
673 ValidIterator end = collIter.end();
674
675 // Count the remote indices we know.
676 for(ValidIterator valid = collIter.begin(); valid != end; ++valid) {
677 ++knownRemote;
678 }
679
680 if(knownRemote>0) {
681 Dune::dverb<<rank_<<": publishing "<<knownRemote<<" for index "<<index->global()<< " for processes ";
682
683 // Update MessageInformation
684 for(ValidIterator valid = collIter.begin(); valid != end; ++valid) {
685 ++(infoSend_[valid.process()].publish);
686 (infoSend_[valid.process()].pairs) += knownRemote;
687 Dune::dverb<<valid.process()<<" ";
688 Dune::dverb<<"(publish="<<infoSend_[valid.process()].publish<<", pairs="<<infoSend_[valid.process()].pairs
689 <<") ";
690 }
691 Dune::dverb<<std::endl;
692 }
693 }
694
695 typedef typename std::map<int,MessageInformation>::const_iterator
696 MessageIterator;
697
698 const MessageIterator end = infoSend_.end();
699
700 // registerMessageDatatype();
701
702 // Now determine the buffersizes needed for each neighbour using MPI_Pack_size
703 MessageInformation dummy;
704
705 MessageIterator messageIter= infoSend_.begin();
706
707 typedef typename RemoteIndices::RemoteIndexMap::const_iterator RemoteIterator;
708 const RemoteIterator rend = remoteIndices_.end();
709 int neighbour=0;
710
711 for(RemoteIterator remote = remoteIndices_.begin(); remote != rend; ++remote, ++neighbour) {
712 MessageInformation* message;
713 MessageInformation recv;
714
715 if(messageIter != end && messageIter->first==remote->first) {
716 // We want to send message information to that process
717 message = const_cast<MessageInformation*>(&(messageIter->second));
718 ++messageIter;
719 }else
720 // We do not want to send information but the other process might.
721 message = &dummy;
722
723 sendBufferSizes_[neighbour]=0;
724 int tsize;
725 // The number of indices published
726 MPI_Pack_size(1, MPI_INT,remoteIndices_.communicator(), &tsize);
727 sendBufferSizes_[neighbour] += tsize;
728
729 for(int i=0; i < message->publish; ++i) {
730 // The global index
731 MPI_Pack_size(1, MPITraits<GlobalIndex>::getType(), remoteIndices_.communicator(), &tsize);
732 sendBufferSizes_[neighbour] += tsize;
733 // The attribute in the local index
734 MPI_Pack_size(1, MPI_CHAR, remoteIndices_.communicator(), &tsize);
735 sendBufferSizes_[neighbour] += tsize;
736 // The number of corresponding remote indices
737 MPI_Pack_size(1, MPI_INT, remoteIndices_.communicator(), &tsize);
738 sendBufferSizes_[neighbour] += tsize;
739 }
740 for(int i=0; i < message->pairs; ++i) {
741 // The process of the remote index
742 MPI_Pack_size(1, MPI_INT, remoteIndices_.communicator(), &tsize);
743 sendBufferSizes_[neighbour] += tsize;
744 // The attribute of the remote index
745 MPI_Pack_size(1, MPI_CHAR, remoteIndices_.communicator(), &tsize);
746 sendBufferSizes_[neighbour] += tsize;
747 }
748
749 Dune::dverb<<rank_<<": Buffer (neighbour="<<remote->first<<") size is "<< sendBufferSizes_[neighbour]<<" for publish="<<message->publish<<" pairs="<<message->pairs<<std::endl;
750 }
751
752 }
753
754 template<typename T>
756 {
757 DefaultNumberer numberer;
758 sync(numberer);
759 }
760
761 template<typename T>
762 template<typename T1>
763 void IndicesSyncer<T>::sync(T1& numberer)
764 {
765
766 // The pointers to the local indices in the remote indices
767 // will become invalid due to the resorting of the index set.
768 // Therefore store the corresponding global indices.
769 // Mark all indices as not added
770
771 typedef typename RemoteIndices::RemoteIndexMap::const_iterator
772 RemoteIterator;
773
774 const RemoteIterator end = remoteIndices_.end();
775
776 // Number of neighbours might change during the syncing.
777 // save the old neighbours
778 std::size_t noOldNeighbours = remoteIndices_.neighbours();
779 int* oldNeighbours = new int[noOldNeighbours];
780 sendBufferSizes_ = new std::size_t[noOldNeighbours];
781 std::size_t i=0;
782
783 for(RemoteIterator remote = remoteIndices_.begin(); remote != end; ++remote, ++i) {
784 typedef typename RemoteIndices::RemoteIndexList::const_iterator
786
787 oldNeighbours[i]=remote->first;
788
789 // Make sure we only have one remote index list.
790 assert(remote->second.first==remote->second.second);
791
792 RemoteIndexList& rList = *(remote->second.first);
793
794 // Store the corresponding global indices.
795 GlobalIndexList& global = globalMap_[remote->first];
796 BoolList& added = oldMap_[remote->first];
797 RemoteIndexIterator riEnd = rList.end();
798
799 for(RemoteIndexIterator index = rList.begin();
800 index != riEnd; ++index) {
801 global.push_back(std::make_pair(index->localIndexPair().global(),
802 index->localIndexPair().local().attribute()));
803 added.push_back(true);
804 }
805
806 Iterators iterators(rList, global, added);
807 iteratorsMap_.insert(std::make_pair(remote->first, iterators));
808 assert(checkReset(iteratorsMap_[remote->first], rList,global,added));
809 }
810
811 // Exchange indices with each neighbour
812 calculateMessageSizes();
813
814 // Allocate the buffers
815 receiveBufferSize_=1;
816 sendBuffers_ = new char*[noOldNeighbours];
817
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]));
821 }
822
823 receiveBuffer_=new char[receiveBufferSize_];
824
825 indexSet_.beginResize();
826
827 Dune::dverb<<rank_<<": Neighbours: ";
828
829 for(i = 0; i<noOldNeighbours; ++i)
830 Dune::dverb<<oldNeighbours[i]<<" ";
831
832 Dune::dverb<<std::endl;
833
834 MPI_Request* requests = new MPI_Request[noOldNeighbours];
835 MPI_Status* statuses = new MPI_Status[noOldNeighbours];
836
837 // Pack Message data and start the sends
838 for(i = 0; i<noOldNeighbours; ++i)
839 packAndSend(oldNeighbours[i], sendBuffers_[i], sendBufferSizes_[i], requests[i]);
840
841 // Probe for incoming messages, receive and unpack them
842 for(i = 0; i<noOldNeighbours; ++i)
843 recvAndUnpack(numberer);
844 // }else{
845 // recvAndUnpack(oldNeighbours[i], numberer);
846 // packAndSend(oldNeighbours[i]);
847 // }
848 // }
849
850 delete[] receiveBuffer_;
851
852 // Wait for the completion of the sends
853 // Wait for completion of sends
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;
859 }
860
861 delete[] statuses;
862 delete[] requests;
863
864 for(std::size_t i=0; i<noOldNeighbours; ++i)
865 delete[] sendBuffers_[i];
866
867 delete[] sendBuffers_;
868 delete[] sendBufferSizes_;
869
870 // No need for the iterator tuples any more
871 iteratorsMap_.clear();
872
873 indexSet_.endResize();
874
875 delete[] oldNeighbours;
876
877 repairLocalIndexPointers(globalMap_, remoteIndices_, indexSet_);
878
879 oldMap_.clear();
880 globalMap_.clear();
881
882 // update the sequence number
883 remoteIndices_.sourceSeqNo_ = remoteIndices_.destSeqNo_ = indexSet_.seqNo();
884 }
885
886 template<typename T>
887 void IndicesSyncer<T>::packAndSend(int destination, char* buffer, std::size_t bufferSize, MPI_Request& request)
888 {
889 typedef typename ParallelIndexSet::const_iterator IndexIterator;
890
891 IndexIterator iEnd = indexSet_.end();
892 int bpos = 0;
893 int published = 0;
894 int pairs = 0;
895
896 assert(checkReset());
897
898 // Pack the number of indices we publish
899 MPI_Pack(&(infoSend_[destination].publish), 1, MPI_INT, buffer, bufferSize, &bpos,
900 remoteIndices_.communicator());
901
902 for(IndexIterator index = indexSet_.begin(); index != iEnd; ++index) {
903 // Search for corresponding remote indices in all iterator tuples
904 typedef typename IteratorsMap::iterator Iterator;
905 Iterator iteratorsEnd = iteratorsMap_.end();
906
907 // advance all iterators to a position with global index >= index->global()
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());
913 }
914
915 // Add all remote indices positioned at global which were already present before calling sync
916 // to the message.
917 // Count how many remote indices we will send
918 int indices = 0;
919 bool knownRemote = false; // Is the remote process supposed to know this index?
920
921 for(Iterator iterators = iteratorsMap_.begin(); iteratorsEnd != iterators; ++iterators)
922 {
923 std::pair<GlobalIndex,Attribute> p;
924 if (iterators->second.isNotAtEnd())
925 {
926 p = iterators->second.globalIndexPair();
927 }
928
929 if(iterators->second.isNotAtEnd() && iterators->second.isOld()
930 && iterators->second.globalIndexPair().first == index->global()) {
931 indices++;
932 if(destination == iterators->first)
933 knownRemote = true;
934 }
935 }
936
937 if(!knownRemote)
938 // We do not need to send any indices
939 continue;
940
941 Dune::dverb<<rank_<<": sending "<<indices<<" for index "<<index->global()<<" to "<<destination<<std::endl;
942
943
944 // Pack the global index, the attribute and the number
945 MPI_Pack(const_cast<GlobalIndex*>(&(index->global())), 1, MPITraits<GlobalIndex>::getType(), buffer, bufferSize, &bpos,
946 remoteIndices_.communicator());
947
948 char attr = index->local().attribute();
949 MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos,
950 remoteIndices_.communicator());
951
952 // Pack the number of remote indices we send.
953 MPI_Pack(&indices, 1, MPI_INT, buffer, bufferSize, &bpos,
954 remoteIndices_.communicator());
955
956 // Pack the information about the remote indices
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;
961
962 ++pairs;
963 assert(pairs <= infoSend_[destination].pairs);
964 MPI_Pack(&process, 1, MPI_INT, buffer, bufferSize, &bpos,
965 remoteIndices_.communicator());
966 char attr = iterators->second.remoteIndex().attribute();
967
968 MPI_Pack(&attr, 1, MPI_CHAR, buffer, bufferSize, &bpos,
969 remoteIndices_.communicator());
970 --indices;
971 }
972 assert(indices==0);
973 ++published;
974 Dune::dvverb<<" (publish="<<published<<", pairs="<<pairs<<")"<<std::endl;
975 assert(published <= infoSend_[destination].publish);
976 }
977
978 // Make sure we send all expected entries
979 assert(published == infoSend_[destination].publish);
980 assert(pairs == infoSend_[destination].pairs);
981 resetIteratorsMap();
982
983 Dune::dverb << rank_<<": Sending message of "<<bpos<<" bytes to "<<destination<<std::endl;
984
985 MPI_Issend(buffer, bpos, MPI_PACKED, destination, 345, remoteIndices_.communicator(),&request);
986 }
987
988 template<typename T>
989 inline void IndicesSyncer<T>::insertIntoRemoteIndexList(int process,
990 const std::pair<GlobalIndex,Attribute>& globalPair,
991 char attribute)
992 {
993 Dune::dverb<<"Inserting from "<<process<<" "<<globalPair.first<<", "<<
994 globalPair.second<<" "<<attribute<<std::endl;
995
996 resetIteratorsMap();
997
998 // There might be cases where there no remote indices for that process yet
999 typename IteratorsMap::iterator found = iteratorsMap_.find(process);
1000
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;
1007 }
1008
1009 Iterators& iterators = found->second;
1010
1011 // Search for the remote index
1012 while(iterators.isNotAtEnd() && iterators.globalIndexPair() < globalPair) {
1013 // Increment all iterators
1014 ++iterators;
1015
1016 }
1017
1018 if(iterators.isAtEnd() || iterators.globalIndexPair() != globalPair) {
1019 // The entry is not yet known
1020 // Insert in the the list and do not change the first iterator.
1021 iterators.insert(RemoteIndex(Attribute(attribute)),globalPair);
1022 return;
1023 }
1024
1025 // Global indices match
1026 bool indexIsThere=false;
1027 for(Iterators tmpIterators = iterators;
1028 !tmpIterators.isAtEnd() && tmpIterators.globalIndexPair() == globalPair;
1029 ++tmpIterators)
1030 //entry already exists with the same attribute
1031 if(tmpIterators.globalIndexPair().second == attribute) {
1032 indexIsThere=true;
1033 break;
1034 }
1035
1036 if(!indexIsThere)
1037 // The entry is not yet known
1038 // Insert in the the list and do not change the first iterator.
1039 iterators.insert(RemoteIndex(Attribute(attribute)),globalPair);
1040 }
1041
1042 template<typename T>
1043 template<typename T1>
1044 void IndicesSyncer<T>::recvAndUnpack(T1& numberer)
1045 {
1046 typedef typename ParallelIndexSet::const_iterator IndexIterator;
1047
1048 IndexIterator iEnd = indexSet_.end();
1049 IndexIterator index = indexSet_.begin();
1050 int bpos = 0;
1051 int publish;
1052
1053 assert(checkReset());
1054
1055 MPI_Status status;
1056
1057 // We have to determine the message size and source before the receive
1058 MPI_Probe(MPI_ANY_SOURCE, 345, remoteIndices_.communicator(), &status);
1059
1060 int source=status.MPI_SOURCE;
1061 int count;
1062 MPI_Get_count(&status, MPI_PACKED, &count);
1063
1064 Dune::dvverb<<rank_<<": Receiving message from "<< source<<" with "<<count<<" bytes"<<std::endl;
1065
1066 if(count>receiveBufferSize_) {
1067 receiveBufferSize_=count;
1068 delete[] receiveBuffer_;
1069 receiveBuffer_ = new char[receiveBufferSize_];
1070 }
1071
1072 MPI_Recv(receiveBuffer_, count, MPI_PACKED, source, 345, remoteIndices_.communicator(), &status);
1073
1074 // How many global entries were published?
1075 MPI_Unpack(receiveBuffer_, count, &bpos, &publish, 1, MPI_INT, remoteIndices_.communicator());
1076
1077 // Now unpack the remote indices and add them.
1078 while(publish>0) {
1079
1080 // Unpack information about the local index on the source process
1081 GlobalIndex global; // global index of the current entry
1082 char sourceAttribute; // Attribute on the source process
1083 int pairs;
1084
1085 MPI_Unpack(receiveBuffer_, count, &bpos, &global, 1, MPITraits<GlobalIndex>::getType(),
1086 remoteIndices_.communicator());
1087 MPI_Unpack(receiveBuffer_, count, &bpos, &sourceAttribute, 1, MPI_CHAR,
1088 remoteIndices_.communicator());
1089 MPI_Unpack(receiveBuffer_, count, &bpos, &pairs, 1, MPI_INT,
1090 remoteIndices_.communicator());
1091
1092 // Insert the entry on the remote process to our
1093 // remote index list
1094 SLList<std::pair<int,Attribute> > sourceAttributeList;
1095 sourceAttributeList.push_back(std::make_pair(source,Attribute(sourceAttribute)));
1096#ifndef NDEBUG
1097 bool foundSelf = false;
1098#endif
1099 Attribute myAttribute=Attribute();
1100
1101 // Unpack the remote indices
1102 for(; pairs>0; --pairs) {
1103 // Unpack the process id that knows the index
1104 int process;
1105 char attribute;
1106 MPI_Unpack(receiveBuffer_, count, &bpos, &process, 1, MPI_INT,
1107 remoteIndices_.communicator());
1108 // Unpack the attribute
1109 MPI_Unpack(receiveBuffer_, count, &bpos, &attribute, 1, MPI_CHAR,
1110 remoteIndices_.communicator());
1111
1112 if(process==rank_) {
1113#ifndef NDEBUG
1114 foundSelf=true;
1115#endif
1116 myAttribute=Attribute(attribute);
1117 // Now we know the local attribute of the global index
1118 //Only add the index if it is unknown.
1119 // Do we know that global index already?
1120 IndexIterator pos = std::lower_bound(index, iEnd, IndexPair(global));
1121
1122 if(pos == iEnd || pos->global() != global) {
1123 // no entry with this global index
1124 indexSet_.add(global,
1125 ParallelLocalIndex<Attribute>(numberer(global),
1126 myAttribute, true));
1127 Dune::dvverb << "Adding "<<global<<" "<<myAttribute<<std::endl;
1128 continue;
1129 }
1130
1131 // because of above the global indices match. Add only if the attribute is different
1132 bool indexIsThere = false;
1133 index=pos;
1134
1135 for(; pos->global()==global; ++pos)
1136 if(pos->local().attribute() == myAttribute) {
1137 Dune::dvverb<<"found "<<global<<" "<<myAttribute<<std::endl;
1138 indexIsThere = true;
1139 break;
1140 }
1141
1142 if(!indexIsThere) {
1143 indexSet_.add(global,
1144 ParallelLocalIndex<Attribute>(numberer(global),
1145 myAttribute, true));
1146 Dune::dvverb << "Adding "<<global<<" "<<myAttribute<<std::endl;
1147 }
1148
1149 }else{
1150 sourceAttributeList.push_back(std::make_pair(process,Attribute(attribute)));
1151 }
1152 }
1153 assert(foundSelf);
1154 // Insert remote indices
1155 typedef typename SLList<std::pair<int,Attribute> >::const_iterator Iter;
1156 for(Iter i=sourceAttributeList.begin(), end=sourceAttributeList.end();
1157 i!=end; ++i)
1158 insertIntoRemoteIndexList(i->first, std::make_pair(global, myAttribute),
1159 i->second);
1160 --publish;
1161 }
1162
1163 resetIteratorsMap();
1164 }
1165
1166 template<typename T>
1167 void IndicesSyncer<T>::resetIteratorsMap(){
1168
1169 // Reset iterators in all tuples.
1170 typedef typename IteratorsMap::iterator Iterator;
1171 typedef typename RemoteIndices::RemoteIndexMap::iterator
1172 RemoteIterator;
1173 typedef typename GlobalIndicesMap::iterator GlobalIterator;
1174 typedef typename BoolMap::iterator BoolIterator;
1175
1176 const RemoteIterator remoteEnd = remoteIndices_.remoteIndices_.end();
1177 Iterator iterators = iteratorsMap_.begin();
1178 GlobalIterator global = globalMap_.begin();
1179 BoolIterator added = oldMap_.begin();
1180
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);
1184 }
1185 }
1186
1187 template<typename T>
1188 bool IndicesSyncer<T>::checkReset(const Iterators& iterators, RemoteIndexList& rList, GlobalIndexList& gList,
1189 BoolList& bList){
1190
1191 if(get<0>(iterators.iterators_)!=rList.begin())
1192 return false;
1193 if(get<1>(iterators.iterators_)!=gList.begin())
1194 return false;
1195 if(get<2>(iterators.iterators_)!=bList.begin())
1196 return false;
1197 return true;
1198 }
1199
1200
1201 template<typename T>
1202 bool IndicesSyncer<T>::checkReset(){
1203
1204 // Reset iterators in all tuples.
1205 typedef typename IteratorsMap::iterator Iterator;
1206 typedef typename RemoteIndices::RemoteIndexMap::iterator
1207 RemoteIterator;
1208 typedef typename GlobalIndicesMap::iterator GlobalIterator;
1209 typedef typename BoolMap::iterator BoolIterator;
1210
1211 const RemoteIterator remoteEnd = remoteIndices_.remoteIndices_.end();
1212 Iterator iterators = iteratorsMap_.begin();
1213 GlobalIterator global = globalMap_.begin();
1214 BoolIterator added = oldMap_.begin();
1215 bool ret = true;
1216
1217 for(RemoteIterator remote = remoteIndices_.remoteIndices_.begin();
1218 remote != remoteEnd; ++remote, ++global, ++added, ++iterators) {
1219 if(!checkReset(iterators->second, *(remote->second.first), global->second,
1220 added->second))
1221 ret=false;
1222 }
1223 return ret;
1224 }
1225}
1226
1227#endif
1228#endif
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)