13 #include "Teuchos_StandardParameterEntryValidators.hpp"
24 if (sendType == DISTRIBUTOR_ISEND) {
27 else if (sendType == DISTRIBUTOR_SEND) {
30 else if (sendType == DISTRIBUTOR_ALLTOALL) {
33 #if defined(HAVE_TPETRA_MPI)
34 else if (sendType == DISTRIBUTOR_IALLTOFEWV) {
38 #if defined(HAVE_TPETRACORE_MPI_ADVANCE)
39 else if (sendType == DISTRIBUTOR_MPIADVANCE_ALLTOALL) {
40 return "MpiAdvanceAlltoall";
42 else if (sendType == DISTRIBUTOR_MPIADVANCE_NBRALLTOALLV) {
43 return "MpiAdvanceNbralltoallv";
47 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Invalid "
48 "EDistributorSendType enum value " << sendType <<
".");
55 if (s ==
"Isend")
return DISTRIBUTOR_ISEND;
56 if (s ==
"Send")
return DISTRIBUTOR_SEND;
57 if (s ==
"Alltoall")
return DISTRIBUTOR_ALLTOALL;
58 #if defined(HAVE_TPETRA_MPI)
59 if (s ==
"Ialltofewv")
return DISTRIBUTOR_IALLTOFEWV;
61 #if defined(HAVE_TPETRACORE_MPI_ADVANCE)
62 if (s ==
"MpiAdvanceAlltoall")
return DISTRIBUTOR_MPIADVANCE_ALLTOALL;
63 if (s ==
"MpiAdvanceNbralltoallv")
return DISTRIBUTOR_MPIADVANCE_NBRALLTOALLV;
65 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Invalid string to convert to EDistributorSendType enum value: " << s);
71 if (std::find(valids.begin(), valids.end(), s) == valids.end()) {
72 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Invalid string for EDistributorSendType enum value: " << s);
81 case Details::DISTRIBUTOR_NOT_INITIALIZED:
82 return "Not initialized yet";
83 case Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_SENDS:
84 return "By createFromSends";
85 case Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_RECVS:
86 return "By createFromRecvs";
87 case Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_SENDS_N_RECVS:
88 return "By createFromSendsAndRecvs";
89 case Details::DISTRIBUTOR_INITIALIZED_BY_REVERSE:
90 return "By createReverseDistributor";
91 case Details::DISTRIBUTOR_INITIALIZED_BY_COPY:
92 return "By copy constructor";
98 DistributorPlan::DistributorPlan(Teuchos::RCP<
const Teuchos::Comm<int>> comm)
100 #if defined(HAVE_TPETRACORE_MPI_ADVANCE)
101 mpixComm_(Teuchos::null),
103 howInitialized_(DISTRIBUTOR_NOT_INITIALIZED),
104 reversePlan_(Teuchos::null),
106 sendMessageToSelf_(false),
107 numSendsToOtherProcs_(0),
110 totalReceiveLength_(0)
113 DistributorPlan::DistributorPlan(
const DistributorPlan& otherPlan)
114 : comm_(otherPlan.comm_),
115 #if defined(HAVE_TPETRACORE_MPI_ADVANCE)
116 mpixComm_(otherPlan.mpixComm_),
118 howInitialized_(DISTRIBUTOR_INITIALIZED_BY_COPY),
119 reversePlan_(otherPlan.reversePlan_),
120 sendType_(otherPlan.sendType_),
121 sendMessageToSelf_(otherPlan.sendMessageToSelf_),
122 numSendsToOtherProcs_(otherPlan.numSendsToOtherProcs_),
123 procIdsToSendTo_(otherPlan.procIdsToSendTo_),
124 startsTo_(otherPlan.startsTo_),
125 lengthsTo_(otherPlan.lengthsTo_),
126 maxSendLength_(otherPlan.maxSendLength_),
127 indicesTo_(otherPlan.indicesTo_),
128 numReceives_(otherPlan.numReceives_),
129 totalReceiveLength_(otherPlan.totalReceiveLength_),
130 lengthsFrom_(otherPlan.lengthsFrom_),
131 procsFrom_(otherPlan.procsFrom_),
132 startsFrom_(otherPlan.startsFrom_),
133 indicesFrom_(otherPlan.indicesFrom_)
134 #if defined(HAVE_TPETRACORE_MPI)
136 roots_(otherPlan.roots_)
140 size_t DistributorPlan::createFromSends(
const Teuchos::ArrayView<const int>& exportProcIDs) {
141 using Teuchos::outArg;
142 using Teuchos::REDUCE_MAX;
143 using Teuchos::reduceAll;
145 const char rawPrefix[] =
"Tpetra::DistributorPlan::createFromSends";
147 const size_t numExports = exportProcIDs.size();
148 const int myProcID = comm_->getRank();
149 const int numProcs = comm_->getSize();
193 for (
size_t i = 0; i < numExports; ++i) {
194 const int exportID = exportProcIDs[i];
195 if (exportID >= numProcs || exportID < 0) {
201 reduceAll<int, int> (*comm_, REDUCE_MAX, badID, outArg (gbl_badID));
202 TEUCHOS_TEST_FOR_EXCEPTION
203 (gbl_badID >= 0, std::runtime_error, rawPrefix <<
"Proc "
204 << gbl_badID <<
", perhaps among other processes, got a bad "
220 Teuchos::Array<size_t> starts (numProcs + 1, 0);
223 size_t numActive = 0;
224 int needSendBuff = 0;
226 for (
size_t i = 0; i < numExports; ++i) {
227 const int exportID = exportProcIDs[i];
242 if (needSendBuff == 0 && starts[exportID] > 1 &&
243 exportID != exportProcIDs[i-1]) {
250 #if defined(HAVE_TPETRA_PRINT_EFFICIENCY_WARNINGS)
252 int global_needSendBuff;
253 reduceAll<int, int> (*comm_, REDUCE_MAX, needSendBuff,
254 outArg (global_needSendBuff));
256 global_needSendBuff != 0,
257 "::createFromSends: Grouping export IDs together by process rank often "
258 "improves performance.");
264 if (starts[myProcID] != 0) {
265 sendMessageToSelf_ =
true;
268 sendMessageToSelf_ =
false;
271 if (! needSendBuff) {
273 numSendsToOtherProcs_ = 0;
276 for (
int i = 0; i < numProcs; ++i) {
278 ++numSendsToOtherProcs_;
284 indicesTo_.resize(0);
287 procIdsToSendTo_.assign(numSendsToOtherProcs_,0);
288 startsTo_.assign(numSendsToOtherProcs_,0);
289 lengthsTo_.assign(numSendsToOtherProcs_,0);
296 size_t procIndex = 0;
297 for (
size_t i = 0; i < numSendsToOtherProcs_; ++i) {
298 while (exportProcIDs[procIndex] < 0) {
301 startsTo_[i] = procIndex;
302 int procID = exportProcIDs[procIndex];
303 procIdsToSendTo_[i] = procID;
304 procIndex += starts[procID];
309 if (numSendsToOtherProcs_ > 0) {
310 sort2(procIdsToSendTo_.begin(), procIdsToSendTo_.end(), startsTo_.begin());
314 for (
size_t i = 0; i < numSendsToOtherProcs_; ++i) {
315 int procID = procIdsToSendTo_[i];
316 lengthsTo_[i] = starts[procID];
317 if ((procID != myProcID) && (lengthsTo_[i] > maxSendLength_)) {
318 maxSendLength_ = lengthsTo_[i];
329 if (starts[0] == 0 ) {
330 numSendsToOtherProcs_ = 0;
333 numSendsToOtherProcs_ = 1;
335 for (Teuchos::Array<size_t>::iterator i=starts.begin()+1,
337 i != starts.end(); ++i)
339 if (*i != 0) ++numSendsToOtherProcs_;
345 for (Teuchos::Array<size_t>::reverse_iterator ip1=starts.rbegin(),
347 i != starts.rend(); ++i)
356 indicesTo_.resize(numActive);
358 for (
size_t i = 0; i < numExports; ++i) {
359 if (exportProcIDs[i] >= 0) {
361 indicesTo_[starts[exportProcIDs[i]]] = i;
363 ++starts[exportProcIDs[i]];
375 for (
int proc = numProcs-1; proc != 0; --proc) {
376 starts[proc] = starts[proc-1];
379 starts[numProcs] = numActive;
386 procIdsToSendTo_.resize(numSendsToOtherProcs_);
387 startsTo_.resize(numSendsToOtherProcs_);
388 lengthsTo_.resize(numSendsToOtherProcs_);
395 for (
int proc = 0; proc < numProcs; ++proc ) {
396 if (starts[proc+1] != starts[proc]) {
397 lengthsTo_[snd] = starts[proc+1] - starts[proc];
398 startsTo_[snd] = starts[proc];
400 if ((proc != myProcID) && (lengthsTo_[snd] > maxSendLength_)) {
401 maxSendLength_ = lengthsTo_[snd];
403 procIdsToSendTo_[snd] = proc;
409 if (sendMessageToSelf_) {
410 --numSendsToOtherProcs_;
416 #if defined(HAVE_TPETRA_MPI)
417 maybeInitializeRoots();
422 howInitialized_ = Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_SENDS;
424 return totalReceiveLength_;
427 void DistributorPlan::createFromRecvs(
const Teuchos::ArrayView<const int>& remoteProcIDs)
429 *
this = *getReversePlan();
430 howInitialized_ = Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_RECVS;
433 void DistributorPlan::createFromSendsAndRecvs(
const Teuchos::ArrayView<const int>& exportProcIDs,
434 const Teuchos::ArrayView<const int>& remoteProcIDs)
443 howInitialized_ = Tpetra::Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_SENDS_N_RECVS;
445 int myProcID = comm_->getRank ();
446 int numProcs = comm_->getSize();
448 const size_t numExportIDs = exportProcIDs.size();
449 Teuchos::Array<size_t> starts (numProcs + 1, 0);
451 size_t numActive = 0;
452 int needSendBuff = 0;
454 for(
size_t i = 0; i < numExportIDs; i++ )
456 if( needSendBuff==0 && i && (exportProcIDs[i] < exportProcIDs[i-1]) )
458 if( exportProcIDs[i] >= 0 )
460 ++starts[ exportProcIDs[i] ];
465 sendMessageToSelf_ = ( starts[myProcID] != 0 ) ? 1 : 0;
467 numSendsToOtherProcs_ = 0;
471 if (starts[0] == 0 ) {
472 numSendsToOtherProcs_ = 0;
475 numSendsToOtherProcs_ = 1;
477 for (Teuchos::Array<size_t>::iterator i=starts.begin()+1,
479 i != starts.end(); ++i)
481 if (*i != 0) ++numSendsToOtherProcs_;
487 for (Teuchos::Array<size_t>::reverse_iterator ip1=starts.rbegin(),
489 i != starts.rend(); ++i)
498 indicesTo_.resize(numActive);
500 for (
size_t i = 0; i < numExportIDs; ++i) {
501 if (exportProcIDs[i] >= 0) {
503 indicesTo_[starts[exportProcIDs[i]]] = i;
505 ++starts[exportProcIDs[i]];
508 for (
int proc = numProcs-1; proc != 0; --proc) {
509 starts[proc] = starts[proc-1];
512 starts[numProcs] = numActive;
513 procIdsToSendTo_.resize(numSendsToOtherProcs_);
514 startsTo_.resize(numSendsToOtherProcs_);
515 lengthsTo_.resize(numSendsToOtherProcs_);
518 for (
int proc = 0; proc < numProcs; ++proc ) {
519 if (starts[proc+1] != starts[proc]) {
520 lengthsTo_[snd] = starts[proc+1] - starts[proc];
521 startsTo_[snd] = starts[proc];
523 if ((proc != myProcID) && (lengthsTo_[snd] > maxSendLength_)) {
524 maxSendLength_ = lengthsTo_[snd];
526 procIdsToSendTo_[snd] = proc;
533 numSendsToOtherProcs_ = 0;
536 for (
int i = 0; i < numProcs; ++i) {
538 ++numSendsToOtherProcs_;
544 indicesTo_.resize(0);
547 procIdsToSendTo_.assign(numSendsToOtherProcs_,0);
548 startsTo_.assign(numSendsToOtherProcs_,0);
549 lengthsTo_.assign(numSendsToOtherProcs_,0);
556 size_t procIndex = 0;
557 for (
size_t i = 0; i < numSendsToOtherProcs_; ++i) {
558 while (exportProcIDs[procIndex] < 0) {
561 startsTo_[i] = procIndex;
562 int procID = exportProcIDs[procIndex];
563 procIdsToSendTo_[i] = procID;
564 procIndex += starts[procID];
569 if (numSendsToOtherProcs_ > 0) {
570 sort2(procIdsToSendTo_.begin(), procIdsToSendTo_.end(), startsTo_.begin());
574 for (
size_t i = 0; i < numSendsToOtherProcs_; ++i) {
575 int procID = procIdsToSendTo_[i];
576 lengthsTo_[i] = starts[procID];
577 if ((procID != myProcID) && (lengthsTo_[i] > maxSendLength_)) {
578 maxSendLength_ = lengthsTo_[i];
584 numSendsToOtherProcs_ -= sendMessageToSelf_;
585 std::vector<int> recv_list;
586 recv_list.reserve(numSendsToOtherProcs_);
589 for(
int i=0; i<remoteProcIDs.size(); i++) {
590 if(remoteProcIDs[i]>last_pid) {
591 recv_list.push_back(remoteProcIDs[i]);
592 last_pid = remoteProcIDs[i];
594 else if (remoteProcIDs[i]<last_pid)
595 throw std::runtime_error(
"Tpetra::Distributor:::createFromSendsAndRecvs expected RemotePIDs to be in sorted order");
597 numReceives_ = recv_list.size();
599 procsFrom_.assign(numReceives_,0);
600 lengthsFrom_.assign(numReceives_,0);
601 indicesFrom_.assign(numReceives_,0);
602 startsFrom_.assign(numReceives_,0);
604 for(
size_t i=0,j=0; i<numReceives_; ++i) {
606 procsFrom_[i] = recv_list[i];
608 for( ; j<(size_t)remoteProcIDs.size() &&
609 remoteProcIDs[jlast]==remoteProcIDs[j] ; j++){;}
610 lengthsFrom_[i] = j-jlast;
612 totalReceiveLength_ = remoteProcIDs.size();
613 indicesFrom_.clear ();
614 numReceives_-=sendMessageToSelf_;
616 #if defined(HAVE_TPETRA_MPI)
617 maybeInitializeRoots();
621 Teuchos::RCP<DistributorPlan> DistributorPlan::getReversePlan()
const {
622 if (reversePlan_.is_null()) createReversePlan();
626 void DistributorPlan::createReversePlan()
const
628 reversePlan_ = Teuchos::rcp(
new DistributorPlan(comm_));
629 reversePlan_->howInitialized_ = Details::DISTRIBUTOR_INITIALIZED_BY_REVERSE;
630 reversePlan_->sendType_ = sendType_;
632 #if defined(HAVE_TPETRACORE_MPI)
636 if (DISTRIBUTOR_IALLTOFEWV == sendType_) {
637 if (Behavior::verbose()) {
638 std::stringstream ss;
639 ss << __FILE__ <<
":" << __LINE__ <<
" WARNING (Ialltofewv send type): using default for reversed Ialltofewv\n";
640 std::cerr << ss.str();
651 size_t totalSendLength =
652 std::accumulate(lengthsTo_.begin(), lengthsTo_.end(), 0);
657 size_t maxReceiveLength = 0;
658 const int myProcID = comm_->getRank();
659 for (
size_t i=0; i < numReceives_; ++i) {
660 if (procsFrom_[i] != myProcID) {
662 if (lengthsFrom_[i] > maxReceiveLength) {
663 maxReceiveLength = lengthsFrom_[i];
668 reversePlan_->sendMessageToSelf_ = sendMessageToSelf_;
669 reversePlan_->numSendsToOtherProcs_ = numReceives_;
670 reversePlan_->procIdsToSendTo_ = procsFrom_;
671 reversePlan_->startsTo_ = startsFrom_;
672 reversePlan_->lengthsTo_ = lengthsFrom_;
673 reversePlan_->maxSendLength_ = maxReceiveLength;
674 reversePlan_->indicesTo_ = indicesFrom_;
675 reversePlan_->numReceives_ = numSendsToOtherProcs_;
676 reversePlan_->totalReceiveLength_ = totalSendLength;
677 reversePlan_->lengthsFrom_ = lengthsTo_;
678 reversePlan_->procsFrom_ = procIdsToSendTo_;
679 reversePlan_->startsFrom_ = startsTo_;
680 reversePlan_->indicesFrom_ = indicesTo_;
682 #if defined(HAVE_TPETRACORE_MPI_ADVANCE)
684 reversePlan_->initializeMpiAdvance();
687 #if defined(HAVE_TPETRA_MPI)
688 reversePlan_->maybeInitializeRoots();
692 void DistributorPlan::computeReceives()
694 using Teuchos::Array;
695 using Teuchos::ArrayRCP;
697 using Teuchos::CommStatus;
698 using Teuchos::CommRequest;
699 using Teuchos::ireceive;
702 using Teuchos::REDUCE_SUM;
703 using Teuchos::receive;
704 using Teuchos::reduce;
705 using Teuchos::scatter;
707 using Teuchos::waitAll;
709 const int myRank = comm_->getRank();
710 const int numProcs = comm_->getSize();
712 const int mpiTag = DEFAULT_MPI_TAG;
720 Array<int> toProcsFromMe (numProcs, 0);
721 #ifdef HAVE_TPETRA_DEBUG
722 bool counting_error =
false;
723 #endif // HAVE_TPETRA_DEBUG
724 for (
size_t i = 0; i < (numSendsToOtherProcs_ + (sendMessageToSelf_ ? 1 : 0)); ++i) {
725 #ifdef HAVE_TPETRA_DEBUG
726 if (toProcsFromMe[procIdsToSendTo_[i]] != 0) {
727 counting_error =
true;
729 #endif // HAVE_TPETRA_DEBUG
730 toProcsFromMe[procIdsToSendTo_[i]] = 1;
732 #ifdef HAVE_TPETRA_DEBUG
735 "Tpetra::Distributor::computeReceives: There was an error on at least "
736 "one process in counting the number of messages send by that process to "
737 "the other processs. Please report this bug to the Tpetra developers.",
739 #endif // HAVE_TPETRA_DEBUG
794 Array<int> numRecvsOnEachProc;
795 if (myRank == root) {
796 numRecvsOnEachProc.resize (numProcs);
798 int numReceivesAsInt = 0;
799 reduce<int, int> (toProcsFromMe.getRawPtr (),
800 numRecvsOnEachProc.getRawPtr (),
801 numProcs, REDUCE_SUM, root, *comm_);
802 scatter<int, int> (numRecvsOnEachProc.getRawPtr (), 1,
803 &numReceivesAsInt, 1, root, *comm_);
804 numReceives_ =
static_cast<size_t> (numReceivesAsInt);
810 lengthsFrom_.assign (numReceives_, 0);
811 procsFrom_.assign (numReceives_, 0);
827 const size_t actualNumReceives = numReceives_ - (sendMessageToSelf_ ? 1 : 0);
833 Array<RCP<CommRequest<int> > > requests (actualNumReceives);
834 Array<ArrayRCP<size_t> > lengthsFromBuffers (actualNumReceives);
835 Array<RCP<CommStatus<int> > > statuses (actualNumReceives);
840 const int anySourceProc = MPI_ANY_SOURCE;
842 const int anySourceProc = -1;
846 for (
size_t i = 0; i < actualNumReceives; ++i) {
851 lengthsFromBuffers[i].resize (1);
852 lengthsFromBuffers[i][0] = as<size_t> (0);
853 requests[i] = ireceive<int, size_t> (lengthsFromBuffers[i], anySourceProc,
865 for (
size_t i = 0; i < numSendsToOtherProcs_ + (sendMessageToSelf_ ? 1 : 0); ++i) {
866 if (procIdsToSendTo_[i] != myRank) {
870 const size_t*
const lengthsTo_i = &lengthsTo_[i];
871 send<int, size_t> (lengthsTo_i, 1, as<int> (procIdsToSendTo_[i]), mpiTag, *comm_);
880 lengthsFrom_[numReceives_-1] = lengthsTo_[i];
881 procsFrom_[numReceives_-1] = myRank;
891 waitAll (*comm_, requests (), statuses ());
892 for (
size_t i = 0; i < actualNumReceives; ++i) {
893 lengthsFrom_[i] = *lengthsFromBuffers[i];
894 procsFrom_[i] = statuses[i]->getSourceRank ();
900 sort2 (procsFrom_.begin(), procsFrom_.end(), lengthsFrom_.begin());
903 totalReceiveLength_ =
904 std::accumulate (lengthsFrom_.begin (), lengthsFrom_.end (), 0);
905 indicesFrom_.clear ();
907 startsFrom_.clear ();
908 startsFrom_.reserve (numReceives_);
909 for (
size_t i = 0, j = 0; i < numReceives_; ++i) {
910 startsFrom_.push_back(j);
911 j += lengthsFrom_[i];
914 if (sendMessageToSelf_) {
919 void DistributorPlan::setParameterList(
const Teuchos::RCP<Teuchos::ParameterList>& plist)
921 using Teuchos::FancyOStream;
922 using Teuchos::getIntegralValue;
923 using Teuchos::ParameterList;
924 using Teuchos::parameterList;
928 if (! plist.is_null()) {
929 RCP<const ParameterList> validParams = getValidParameters ();
930 plist->validateParametersAndSetDefaults (*validParams);
933 getIntegralValue<Details::EDistributorSendType> (*plist,
"Send type");
936 sendType_ = sendType;
938 #if defined(HAVE_TPETRACORE_MPI_ADVANCE)
939 initializeMpiAdvance();
944 this->setMyParamList (plist);
946 #if defined(HAVE_TPETRA_MPI)
947 maybeInitializeRoots();
954 Teuchos::Array<std::string> sendTypes;
955 sendTypes.push_back (
"Isend");
956 sendTypes.push_back (
"Send");
957 sendTypes.push_back (
"Alltoall");
958 #if defined(HAVE_TPETRA_MPI)
959 sendTypes.push_back (
"Ialltofewv");
961 #if defined(HAVE_TPETRACORE_MPI_ADVANCE)
962 sendTypes.push_back (
"MpiAdvanceAlltoall");
963 sendTypes.push_back (
"MpiAdvanceNbralltoallv");
969 Teuchos::Array<EDistributorSendType> res;
970 res.push_back (DISTRIBUTOR_ISEND);
971 res.push_back (DISTRIBUTOR_SEND);
972 res.push_back (DISTRIBUTOR_ALLTOALL);
973 #if defined(HAVE_TPETRA_MPI)
974 res.push_back (DISTRIBUTOR_IALLTOFEWV);
976 #if defined(HAVE_TPETRACORE_MPI_ADVANCE)
977 res.push_back (DISTRIBUTOR_MPIADVANCE_ALLTOALL);
978 res.push_back (DISTRIBUTOR_MPIADVANCE_NBRALLTOALLV);
983 Teuchos::RCP<const Teuchos::ParameterList>
984 DistributorPlan::getValidParameters()
const
986 using Teuchos::Array;
987 using Teuchos::ParameterList;
988 using Teuchos::parameterList;
990 using Teuchos::setStringToIntegralParameter;
992 Array<std::string> sendTypes = distributorSendTypes ();
993 const Array<Details::EDistributorSendType> sendTypeEnums = distributorSendTypeEnums ();
995 const std::string validatedSendType = validSendTypeOrThrow(Behavior::defaultSendType());
997 RCP<ParameterList> plist = parameterList (
"Tpetra::Distributor");
999 setStringToIntegralParameter<Details::EDistributorSendType> (
"Send type",
1000 validatedSendType,
"When using MPI, the variant of send to use in "
1001 "do[Reverse]Posts()", sendTypes(), sendTypeEnums(), plist.getRawPtr());
1002 plist->set (
"Timer Label",
"",
"Label for Time Monitor output");
1004 return Teuchos::rcp_const_cast<
const ParameterList> (plist);
1007 #if defined(HAVE_TPETRACORE_MPI_ADVANCE)
1010 struct MpixCommDeallocator {
1011 void free(MPIX_Comm **comm)
const {
1012 MPIX_Comm_free(comm);
1016 void DistributorPlan::initializeMpiAdvance() {
1019 TEUCHOS_ASSERT(mpixComm_.is_null());
1022 Teuchos::RCP<const Teuchos::MpiComm<int> > mpiComm = Teuchos::rcp_dynamic_cast<
const Teuchos::MpiComm<int> >(comm_);
1023 Teuchos::RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > rawComm = mpiComm->getRawMpiComm();
1025 if (sendType_ == DISTRIBUTOR_MPIADVANCE_ALLTOALL ||
1026 sendType_ == DISTRIBUTOR_MPIADVANCE_NBRALLTOALLV ) {
1027 MPIX_Comm **mpixComm =
new(MPIX_Comm*);
1028 err = MPIX_Comm_init(mpixComm, (*rawComm)());
1029 mpixComm_ = Teuchos::RCP(mpixComm,
1030 MpixCommDeallocator(),
1035 TEUCHOS_ASSERT(err == 0);
1039 #if defined(HAVE_TPETRA_MPI)
1041 void DistributorPlan::maybeInitializeRoots() {
1044 if (DISTRIBUTOR_IALLTOFEWV != sendType_) {
1049 ProfilingRegion region_maybeInitializeRoots (
"Tpetra::DistributorPlan::maybeInitializeRoots");
1052 const int numRecvs = (int)(getNumReceives() + (hasSelfMessage() ? 1 : 0));
1053 std::vector<int> sendbuf(comm_->getSize(), numRecvs);
1054 std::vector<int> recvbuf(comm_->getSize());
1059 Teuchos::RCP<const Teuchos::MpiComm<int> > mpiComm = Teuchos::rcp_dynamic_cast<
const Teuchos::MpiComm<int> >(comm_);
1060 Teuchos::RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > rawComm = mpiComm->getRawMpiComm();
1061 MPI_Comm comm = (*rawComm)();
1062 MPI_Alltoall(sendbuf.data(), 1, MPI_INT, recvbuf.data(), 1, MPI_INT, comm);
1065 for (
size_t root = 0; root < recvbuf.size(); ++root) {
1066 if (recvbuf[root] > 0) {
1067 roots_.push_back(root);
1073 int slow = !getIndicesTo().is_null() ? 1 : 0;
1074 MPI_Allreduce(MPI_IN_PLACE, &slow, 1, MPI_INT, MPI_LOR, comm);
1079 std::stringstream ss;
1080 ss << __FILE__ <<
":" << __LINE__ <<
" " << comm_->getRank() <<
": WARNING: Ialltoallv send mode set, at least one rank's data is not grouped by rank. Setting to \"Send\"" << std::endl;
1081 std::cerr << ss.str();
1086 sendType_ = DISTRIBUTOR_SEND;
1092 if (roots_.size() * roots_.size() >= size_t(comm_->getSize())) {
1094 std::stringstream ss;
1095 ss << __FILE__ <<
":" << __LINE__ <<
" " << comm_->getRank() <<
": WARNING (Ialltoallv send type): too many roots (" << roots_.size() <<
") for " << comm_->getSize() <<
" ranks. Setting send-type to \"Send\"" << std::endl;
1096 std::cerr << ss.str();
1099 sendType_ = DISTRIBUTOR_SEND;
1102 #endif // HAVE_TPETRA_MPI
1104 DistributorPlan::SubViewLimits DistributorPlan::getImportViewLimits(
size_t numPackets)
const {
1105 const size_t actualNumReceives = getNumReceives() + (hasSelfMessage() ? 1 : 0);
1107 IndexView importStarts(actualNumReceives);
1108 IndexView importLengths(actualNumReceives);
1111 for (
size_t i = 0; i < actualNumReceives; ++i) {
1112 importStarts[i] = offset;
1113 offset += getLengthsFrom()[i] * numPackets;
1114 importLengths[i] = getLengthsFrom()[i] * numPackets;
1116 return std::make_pair(importStarts, importLengths);
1119 DistributorPlan::SubViewLimits DistributorPlan::getImportViewLimits(
const Teuchos::ArrayView<const size_t> &numImportPacketsPerLID)
const {
1121 const size_t actualNumReceives = getNumReceives() + (hasSelfMessage() ? 1 : 0);
1123 IndexView importStarts(actualNumReceives);
1124 IndexView importLengths(actualNumReceives);
1127 size_t curLIDoffset = 0;
1128 for (
size_t i = 0; i < actualNumReceives; ++i) {
1129 size_t totalPacketsFrom_i = 0;
1130 for (
size_t j = 0; j < getLengthsFrom()[i]; ++j) {
1131 totalPacketsFrom_i += numImportPacketsPerLID[curLIDoffset + j];
1133 curLIDoffset += getLengthsFrom()[i];
1134 importStarts[i] = offset;
1135 offset += totalPacketsFrom_i;
1136 importLengths[i] = totalPacketsFrom_i;
1138 return std::make_pair(importStarts, importLengths);
1142 DistributorPlan::SubViewLimits DistributorPlan::getExportViewLimits(
size_t numPackets)
const {
1143 if (getIndicesTo().is_null()) {
1145 const size_t actualNumSends = getNumSends() + (hasSelfMessage() ? 1 : 0);
1146 IndexView exportStarts(actualNumSends);
1147 IndexView exportLengths(actualNumSends);
1148 for (
size_t pp = 0; pp < actualNumSends; ++pp) {
1149 exportStarts[pp] = getStartsTo()[pp] * numPackets;
1150 exportLengths[pp] = getLengthsTo()[pp] * numPackets;
1152 return std::make_pair(exportStarts, exportLengths);
1154 const size_t numIndices = getIndicesTo().size();
1155 IndexView exportStarts(numIndices);
1156 IndexView exportLengths(numIndices);
1157 for (
size_t j = 0; j < numIndices; ++j) {
1158 exportStarts[j] = getIndicesTo()[j]*numPackets;
1159 exportLengths[j] = numPackets;
1161 return std::make_pair(exportStarts, exportLengths);
1165 DistributorPlan::SubViewLimits DistributorPlan::getExportViewLimits(
const Teuchos::ArrayView<const size_t> &numExportPacketsPerLID)
const {
1166 if (getIndicesTo().is_null()) {
1167 const size_t actualNumSends = getNumSends() + (hasSelfMessage() ? 1 : 0);
1168 IndexView exportStarts(actualNumSends);
1169 IndexView exportLengths(actualNumSends);
1171 for (
size_t pp = 0; pp < actualNumSends; ++pp) {
1172 size_t numPackets = 0;
1173 for (
size_t j = getStartsTo()[pp];
1174 j < getStartsTo()[pp] + getLengthsTo()[pp]; ++j) {
1175 numPackets += numExportPacketsPerLID[j];
1177 exportStarts[pp] = offset;
1178 offset += numPackets;
1179 exportLengths[pp] = numPackets;
1181 return std::make_pair(exportStarts, exportLengths);
1183 const size_t numIndices = getIndicesTo().size();
1184 IndexView exportStarts(numIndices);
1185 IndexView exportLengths(numIndices);
1187 for (
size_t j = 0; j < numIndices; ++j) {
1188 exportStarts[j] = offset;
1189 offset += numExportPacketsPerLID[j];
1190 exportLengths[j] = numExportPacketsPerLID[j];
1192 return std::make_pair(exportStarts, exportLengths);
Declaration of Tpetra::Details::Profiling, a scope guard for Kokkos Profiling.
EDistributorHowInitialized
Enum indicating how and whether a Distributor was initialized.
static bool debug()
Whether Tpetra is in debug mode.
const std::string & validSendTypeOrThrow(const std::string &s)
Valid enum values of distributor send types.
std::string DistributorSendTypeEnumToString(EDistributorSendType sendType)
Convert an EDistributorSendType enum value to a string.
#define TPETRA_EFFICIENCY_WARNING(throw_exception_test, msg)
Print or throw an efficency warning.
EDistributorSendType DistributorSendTypeStringToEnum(const std::string_view s)
Convert a string to an EDistributorSendType. Throw on error.
Teuchos::Array< EDistributorSendType > distributorSendTypeEnums()
Valid enum values of distributor send types.
void sort2(const IT1 &first1, const IT1 &last1, const IT2 &first2, const bool stableSort=false)
Sort the first array, and apply the resulting permutation to the second array.
static bool verbose()
Whether Tpetra is in verbose mode.
std::string DistributorHowInitializedEnumToString(EDistributorHowInitialized how)
Convert an EDistributorHowInitialized enum value to a string.
Stand-alone utility functions and macros.
Teuchos::Array< std::string > distributorSendTypes()
Valid string values for Distributor's "Send type" parameter.
#define SHARED_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg, comm)
Test for exception, with reduction over the given communicator.
EDistributorSendType
The type of MPI send that Distributor should use.
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.