13#include "Teuchos_StandardParameterEntryValidators.hpp" 
   25  } 
else if (
sendType == DISTRIBUTOR_SEND) {
 
   27  } 
else if (
sendType == DISTRIBUTOR_ALLTOALL) {
 
   30#if defined(HAVE_TPETRA_MPI) 
   35#if defined(HAVE_TPETRACORE_MPI_ADVANCE) 
   37    return "MpiAdvanceAlltoall";
 
   39    return "MpiAdvanceNbralltoallv";
 
   45                               "EDistributorSendType enum value " 
 
   52  if (
s == 
"Isend") 
return DISTRIBUTOR_ISEND;
 
   53  if (
s == 
"Send") 
return DISTRIBUTOR_SEND;
 
   54  if (
s == 
"Alltoall") 
return DISTRIBUTOR_ALLTOALL;
 
   55#if defined(HAVE_TPETRA_MPI) 
   58#if defined(HAVE_TPETRACORE_MPI_ADVANCE) 
 
   77    case Details::DISTRIBUTOR_NOT_INITIALIZED:
 
   78      return "Not initialized yet";
 
   79    case Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_SENDS:
 
   80      return "By createFromSends";
 
   81    case Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_RECVS:
 
   82      return "By createFromRecvs";
 
   83    case Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_SENDS_N_RECVS:
 
   84      return "By createFromSendsAndRecvs";
 
   85    case Details::DISTRIBUTOR_INITIALIZED_BY_REVERSE:
 
   86      return "By createReverseDistributor";
 
   87    case Details::DISTRIBUTOR_INITIALIZED_BY_COPY:
 
   88      return "By copy constructor";
 
 
   94DistributorPlan::DistributorPlan(Teuchos::RCP<
const Teuchos::Comm<int>> comm)
 
  101  howInitialized_(DISTRIBUTOR_NOT_INITIALIZED)
 
  102  , reversePlan_(Teuchos::
null)
 
  104  , sendMessageToSelf_(
false)
 
  105  , numSendsToOtherProcs_(0)
 
  108  , totalReceiveLength_(0) {
 
  111DistributorPlan::DistributorPlan(
const DistributorPlan& otherPlan)
 
  112  : comm_(otherPlan.comm_)
 
  114#if defined(HAVE_TPETRACORE_MPI_ADVANCE)
 
  115  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)
 
  135  , roots_(otherPlan.roots_)
 
  140size_t DistributorPlan::createFromSends(
const Teuchos::ArrayView<const int>& exportProcIDs) {
 
  142  using Teuchos::outArg;
 
  143  using Teuchos::REDUCE_MAX;
 
  144  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(gbl_badID >= 0, std::runtime_error, rawPrefix << 
"Proc " << gbl_badID << 
", perhaps among other processes, got a bad " 
  218  Teuchos::Array<size_t> starts(numProcs + 1, 0);
 
  221  size_t numActive = 0;
 
  222  int needSendBuff = 0;  
 
  224  for (
size_t i = 0; i < numExports; ++i) {
 
  225    const int exportID = exportProcIDs[i];
 
  240      if (needSendBuff == 0 && starts[exportID] > 1 &&
 
  241          exportID != exportProcIDs[i - 1]) {
 
  248#if defined(HAVE_TPETRA_PRINT_EFFICIENCY_WARNINGS) 
  250    int global_needSendBuff;
 
  251    reduceAll<int, int>(*comm_, REDUCE_MAX, needSendBuff,
 
  252                        outArg(global_needSendBuff));
 
  254        global_needSendBuff != 0,
 
  255        "::createFromSends: Grouping export IDs together by process rank often " 
  256        "improves performance.");
 
  262  if (starts[myProcID] != 0) {
 
  263    sendMessageToSelf_ = 
true;
 
  265    sendMessageToSelf_ = 
false;
 
  270    numSendsToOtherProcs_ = 0;
 
  273    for (
int i = 0; i < numProcs; ++i) {
 
  275        ++numSendsToOtherProcs_;
 
  281    indicesTo_.resize(0);
 
  284    procIdsToSendTo_.assign(numSendsToOtherProcs_, 0);
 
  285    startsTo_.assign(numSendsToOtherProcs_, 0);
 
  286    lengthsTo_.assign(numSendsToOtherProcs_, 0);
 
  293      size_t procIndex = 0;
 
  294      for (
size_t i = 0; i < numSendsToOtherProcs_; ++i) {
 
  295        while (exportProcIDs[procIndex] < 0) {
 
  298        startsTo_[i]        = procIndex;
 
  299        int procID          = exportProcIDs[procIndex];
 
  300        procIdsToSendTo_[i] = procID;
 
  301        procIndex += starts[procID];
 
  306    if (numSendsToOtherProcs_ > 0) {
 
  307      sort2(procIdsToSendTo_.begin(), procIdsToSendTo_.end(), startsTo_.begin());
 
  311    for (
size_t i = 0; i < numSendsToOtherProcs_; ++i) {
 
  312      int procID    = procIdsToSendTo_[i];
 
  313      lengthsTo_[i] = starts[procID];
 
  314      if ((procID != myProcID) && (lengthsTo_[i] > maxSendLength_)) {
 
  315        maxSendLength_ = lengthsTo_[i];
 
  325    if (starts[0] == 0) {
 
  326      numSendsToOtherProcs_ = 0;
 
  328      numSendsToOtherProcs_ = 1;
 
  330    for (Teuchos::Array<size_t>::iterator i   = starts.begin() + 1,
 
  331                                          im1 = starts.begin();
 
  332         i != starts.end(); ++i) {
 
  333      if (*i != 0) ++numSendsToOtherProcs_;
 
  339    for (Teuchos::Array<size_t>::reverse_iterator ip1 = starts.rbegin(),
 
  340                                                  i   = starts.rbegin() + 1;
 
  341         i != starts.rend(); ++i) {
 
  349    indicesTo_.resize(numActive);
 
  351    for (
size_t i = 0; i < numExports; ++i) {
 
  352      if (exportProcIDs[i] >= 0) {
 
  354        indicesTo_[starts[exportProcIDs[i]]] = i;
 
  356        ++starts[exportProcIDs[i]];
 
  368    for (
int proc = numProcs - 1; proc != 0; --proc) {
 
  369      starts[proc] = starts[proc - 1];
 
  372    starts[numProcs] = numActive;
 
  379    procIdsToSendTo_.resize(numSendsToOtherProcs_);
 
  380    startsTo_.resize(numSendsToOtherProcs_);
 
  381    lengthsTo_.resize(numSendsToOtherProcs_);
 
  388    for (
int proc = 0; proc < numProcs; ++proc) {
 
  389      if (starts[proc + 1] != starts[proc]) {
 
  390        lengthsTo_[snd] = starts[proc + 1] - starts[proc];
 
  391        startsTo_[snd]  = starts[proc];
 
  393        if ((proc != myProcID) && (lengthsTo_[snd] > maxSendLength_)) {
 
  394          maxSendLength_ = lengthsTo_[snd];
 
  396        procIdsToSendTo_[snd] = proc;
 
  402  if (sendMessageToSelf_) {
 
  403    --numSendsToOtherProcs_;
 
  409#if defined(HAVE_TPETRA_MPI) 
  410  maybeInitializeRoots();
 
  415  howInitialized_ = Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_SENDS;
 
  417  return totalReceiveLength_;
 
  420void DistributorPlan::createFromRecvs(
const Teuchos::ArrayView<const int>& remoteProcIDs) {
 
  421  *
this           = *getReversePlan();
 
  422  howInitialized_ = Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_RECVS;
 
  425void DistributorPlan::createFromSendsAndRecvs(
const Teuchos::ArrayView<const int>& exportProcIDs,
 
  426                                              const Teuchos::ArrayView<const int>& remoteProcIDs) {
 
  434  howInitialized_ = Tpetra::Details::DISTRIBUTOR_INITIALIZED_BY_CREATE_FROM_SENDS_N_RECVS;
 
  436  int myProcID = comm_->getRank();
 
  437  int numProcs = comm_->getSize();
 
  439  const size_t numExportIDs = exportProcIDs.size();
 
  440  Teuchos::Array<size_t> starts(numProcs + 1, 0);
 
  442  size_t numActive = 0;
 
  443  int needSendBuff = 0;  
 
  445  for (
size_t i = 0; i < numExportIDs; i++) {
 
  446    if (needSendBuff == 0 && i && (exportProcIDs[i] < exportProcIDs[i - 1]))
 
  448    if (exportProcIDs[i] >= 0) {
 
  449      ++starts[exportProcIDs[i]];
 
  454  sendMessageToSelf_ = (starts[myProcID] != 0) ? 1 : 0;
 
  456  numSendsToOtherProcs_ = 0;
 
  460    if (starts[0] == 0) {
 
  461      numSendsToOtherProcs_ = 0;
 
  463      numSendsToOtherProcs_ = 1;
 
  465    for (Teuchos::Array<size_t>::iterator i   = starts.begin() + 1,
 
  466                                          im1 = starts.begin();
 
  467         i != starts.end(); ++i) {
 
  468      if (*i != 0) ++numSendsToOtherProcs_;
 
  474    for (Teuchos::Array<size_t>::reverse_iterator ip1 = starts.rbegin(),
 
  475                                                  i   = starts.rbegin() + 1;
 
  476         i != starts.rend(); ++i) {
 
  484    indicesTo_.resize(numActive);
 
  486    for (
size_t i = 0; i < numExportIDs; ++i) {
 
  487      if (exportProcIDs[i] >= 0) {
 
  489        indicesTo_[starts[exportProcIDs[i]]] = i;
 
  491        ++starts[exportProcIDs[i]];
 
  494    for (
int proc = numProcs - 1; proc != 0; --proc) {
 
  495      starts[proc] = starts[proc - 1];
 
  498    starts[numProcs] = numActive;
 
  499    procIdsToSendTo_.resize(numSendsToOtherProcs_);
 
  500    startsTo_.resize(numSendsToOtherProcs_);
 
  501    lengthsTo_.resize(numSendsToOtherProcs_);
 
  504    for (
int proc = 0; proc < numProcs; ++proc) {
 
  505      if (starts[proc + 1] != starts[proc]) {
 
  506        lengthsTo_[snd] = starts[proc + 1] - starts[proc];
 
  507        startsTo_[snd]  = starts[proc];
 
  509        if ((proc != myProcID) && (lengthsTo_[snd] > maxSendLength_)) {
 
  510          maxSendLength_ = lengthsTo_[snd];
 
  512        procIdsToSendTo_[snd] = proc;
 
  518    numSendsToOtherProcs_ = 0;
 
  521    for (
int i = 0; i < numProcs; ++i) {
 
  523        ++numSendsToOtherProcs_;
 
  529    indicesTo_.resize(0);
 
  532    procIdsToSendTo_.assign(numSendsToOtherProcs_, 0);
 
  533    startsTo_.assign(numSendsToOtherProcs_, 0);
 
  534    lengthsTo_.assign(numSendsToOtherProcs_, 0);
 
  541      size_t procIndex = 0;
 
  542      for (
size_t i = 0; i < numSendsToOtherProcs_; ++i) {
 
  543        while (exportProcIDs[procIndex] < 0) {
 
  546        startsTo_[i]        = procIndex;
 
  547        int procID          = exportProcIDs[procIndex];
 
  548        procIdsToSendTo_[i] = procID;
 
  549        procIndex += starts[procID];
 
  554    if (numSendsToOtherProcs_ > 0) {
 
  555      sort2(procIdsToSendTo_.begin(), procIdsToSendTo_.end(), startsTo_.begin());
 
  559    for (
size_t i = 0; i < numSendsToOtherProcs_; ++i) {
 
  560      int procID    = procIdsToSendTo_[i];
 
  561      lengthsTo_[i] = starts[procID];
 
  562      if ((procID != myProcID) && (lengthsTo_[i] > maxSendLength_)) {
 
  563        maxSendLength_ = lengthsTo_[i];
 
  568  numSendsToOtherProcs_ -= sendMessageToSelf_;
 
  569  std::vector<int> recv_list;
 
  570  recv_list.reserve(numSendsToOtherProcs_);  
 
  573  for (
int i = 0; i < remoteProcIDs.size(); i++) {
 
  574    if (remoteProcIDs[i] > last_pid) {
 
  575      recv_list.push_back(remoteProcIDs[i]);
 
  576      last_pid = remoteProcIDs[i];
 
  577    } 
else if (remoteProcIDs[i] < last_pid)
 
  578      throw std::runtime_error(
"Tpetra::Distributor:::createFromSendsAndRecvs expected RemotePIDs to be in sorted order");
 
  580  numReceives_ = recv_list.size();
 
  582    procsFrom_.assign(numReceives_, 0);
 
  583    lengthsFrom_.assign(numReceives_, 0);
 
  584    indicesFrom_.assign(numReceives_, 0);
 
  585    startsFrom_.assign(numReceives_, 0);
 
  587  for (
size_t i = 0, j = 0; i < numReceives_; ++i) {
 
  589    procsFrom_[i]  = recv_list[i];
 
  591    for (; j < (size_t)remoteProcIDs.size() &&
 
  592           remoteProcIDs[jlast] == remoteProcIDs[j];
 
  596    lengthsFrom_[i] = j - jlast;
 
  598  totalReceiveLength_ = remoteProcIDs.size();
 
  599  indicesFrom_.clear();
 
  600  numReceives_ -= sendMessageToSelf_;
 
  602#if defined(HAVE_TPETRA_MPI) 
  603  maybeInitializeRoots();
 
  607Teuchos::RCP<DistributorPlan> DistributorPlan::getReversePlan()
 const {
 
  608  if (reversePlan_.is_null()) createReversePlan();
 
  612void DistributorPlan::createReversePlan()
 const {
 
  613  reversePlan_                  = Teuchos::rcp(
new DistributorPlan(comm_));
 
  614  reversePlan_->howInitialized_ = Details::DISTRIBUTOR_INITIALIZED_BY_REVERSE;
 
  615  reversePlan_->sendType_       = sendType_;
 
  617#if defined(HAVE_TPETRACORE_MPI) 
  621  if (DISTRIBUTOR_IALLTOFEWV == sendType_) {
 
  622    if (Behavior::verbose()) {
 
  623      std::stringstream ss;
 
  624      ss << __FILE__ << 
":" << __LINE__ << 
" WARNING (Ialltofewv send type): using default for reversed Ialltofewv\n";
 
  625      std::cerr << ss.str();
 
  635  size_t totalSendLength =
 
  636      std::accumulate(lengthsTo_.begin(), lengthsTo_.end(), 0);
 
  641  size_t maxReceiveLength = 0;
 
  642  const int myProcID      = comm_->getRank();
 
  643  for (
size_t i = 0; i < numReceives_; ++i) {
 
  644    if (procsFrom_[i] != myProcID) {
 
  646      if (lengthsFrom_[i] > maxReceiveLength) {
 
  647        maxReceiveLength = lengthsFrom_[i];
 
  652  reversePlan_->sendMessageToSelf_    = sendMessageToSelf_;
 
  653  reversePlan_->numSendsToOtherProcs_ = numReceives_;
 
  654  reversePlan_->procIdsToSendTo_      = procsFrom_;
 
  655  reversePlan_->startsTo_             = startsFrom_;
 
  656  reversePlan_->lengthsTo_            = lengthsFrom_;
 
  657  reversePlan_->maxSendLength_        = maxReceiveLength;
 
  658  reversePlan_->indicesTo_            = indicesFrom_;
 
  659  reversePlan_->numReceives_          = numSendsToOtherProcs_;
 
  660  reversePlan_->totalReceiveLength_   = totalSendLength;
 
  661  reversePlan_->lengthsFrom_          = lengthsTo_;
 
  662  reversePlan_->procsFrom_            = procIdsToSendTo_;
 
  663  reversePlan_->startsFrom_           = startsTo_;
 
  664  reversePlan_->indicesFrom_          = indicesTo_;
 
  666#if defined(HAVE_TPETRACORE_MPI_ADVANCE) 
  668  reversePlan_->initializeMpiAdvance();
 
  671#if defined(HAVE_TPETRA_MPI) 
  672  reversePlan_->maybeInitializeRoots();
 
  676void DistributorPlan::computeReceives() {
 
  677  using Teuchos::Array;
 
  678  using Teuchos::ArrayRCP;
 
  680  using Teuchos::CommRequest;
 
  681  using Teuchos::CommStatus;
 
  682  using Teuchos::ireceive;
 
  685  using Teuchos::receive;
 
  686  using Teuchos::reduce;
 
  687  using Teuchos::REDUCE_SUM;
 
  688  using Teuchos::scatter;
 
  690  using Teuchos::waitAll;
 
  692  const int myRank   = comm_->getRank();
 
  693  const int numProcs = comm_->getSize();
 
  695  const int mpiTag = DEFAULT_MPI_TAG;
 
  703    Array<int> toProcsFromMe(numProcs, 0);
 
  704#ifdef HAVE_TPETRA_DEBUG 
  705    bool counting_error = 
false;
 
  707    for (
size_t i = 0; i < (numSendsToOtherProcs_ + (sendMessageToSelf_ ? 1 : 0)); ++i) {
 
  708#ifdef HAVE_TPETRA_DEBUG 
  709      if (toProcsFromMe[procIdsToSendTo_[i]] != 0) {
 
  710        counting_error = 
true;
 
  713      toProcsFromMe[procIdsToSendTo_[i]] = 1;
 
  715#ifdef HAVE_TPETRA_DEBUG 
  718                              "Tpetra::Distributor::computeReceives: There was an error on at least " 
  719                              "one process in counting the number of messages send by that process to " 
  720                              "the other processs.  Please report this bug to the Tpetra developers.",
 
  777    Array<int> numRecvsOnEachProc;  
 
  778    if (myRank == root) {
 
  779      numRecvsOnEachProc.resize(numProcs);
 
  781    int numReceivesAsInt = 0;  
 
  782    reduce<int, int>(toProcsFromMe.getRawPtr(),
 
  783                     numRecvsOnEachProc.getRawPtr(),
 
  784                     numProcs, REDUCE_SUM, root, *comm_);
 
  785    scatter<int, int>(numRecvsOnEachProc.getRawPtr(), 1,
 
  786                      &numReceivesAsInt, 1, root, *comm_);
 
  787    numReceives_ = 
static_cast<size_t>(numReceivesAsInt);
 
  793  lengthsFrom_.assign(numReceives_, 0);
 
  794  procsFrom_.assign(numReceives_, 0);
 
  810  const size_t actualNumReceives = numReceives_ - (sendMessageToSelf_ ? 1 : 0);
 
  816  Array<RCP<CommRequest<int>>> requests(actualNumReceives);
 
  817  Array<ArrayRCP<size_t>> lengthsFromBuffers(actualNumReceives);
 
  818  Array<RCP<CommStatus<int>>> statuses(actualNumReceives);
 
  823  const int anySourceProc = MPI_ANY_SOURCE;
 
  825  const int anySourceProc = -1;
 
  829  for (
size_t i = 0; i < actualNumReceives; ++i) {
 
  834    lengthsFromBuffers[i].resize(1);
 
  835    lengthsFromBuffers[i][0] = as<size_t>(0);
 
  836    requests[i]              = ireceive<int, size_t>(lengthsFromBuffers[i], anySourceProc,
 
  848  for (
size_t i = 0; i < numSendsToOtherProcs_ + (sendMessageToSelf_ ? 1 : 0); ++i) {
 
  849    if (procIdsToSendTo_[i] != myRank) {
 
  853      const size_t* 
const lengthsTo_i = &lengthsTo_[i];
 
  854      send<int, size_t>(lengthsTo_i, 1, as<int>(procIdsToSendTo_[i]), mpiTag, *comm_);
 
  862      lengthsFrom_[numReceives_ - 1] = lengthsTo_[i];
 
  863      procsFrom_[numReceives_ - 1]   = myRank;
 
  873  waitAll(*comm_, requests(), statuses());
 
  874  for (
size_t i = 0; i < actualNumReceives; ++i) {
 
  875    lengthsFrom_[i] = *lengthsFromBuffers[i];
 
  876    procsFrom_[i]   = statuses[i]->getSourceRank();
 
  882  sort2(procsFrom_.begin(), procsFrom_.end(), lengthsFrom_.begin());
 
  885  totalReceiveLength_ =
 
  886      std::accumulate(lengthsFrom_.begin(), lengthsFrom_.end(), 0);
 
  887  indicesFrom_.clear();
 
  890  startsFrom_.reserve(numReceives_);
 
  891  for (
size_t i = 0, j = 0; i < numReceives_; ++i) {
 
  892    startsFrom_.push_back(j);
 
  893    j += lengthsFrom_[i];
 
  896  if (sendMessageToSelf_) {
 
  901void DistributorPlan::setParameterList(
const Teuchos::RCP<Teuchos::ParameterList>& plist) {
 
  903  using Teuchos::FancyOStream;
 
  904  using Teuchos::getIntegralValue;
 
  905  using Teuchos::ParameterList;
 
  906  using Teuchos::parameterList;
 
  909  if (!plist.is_null()) {
 
  910    RCP<const ParameterList> validParams = getValidParameters();
 
  911    plist->validateParametersAndSetDefaults(*validParams);
 
  914        getIntegralValue<Details::EDistributorSendType>(*plist, 
"Send type");
 
  917    sendType_ = sendType;
 
  919#if defined(HAVE_TPETRACORE_MPI_ADVANCE) 
  920    initializeMpiAdvance();
 
  925    this->setMyParamList(plist);
 
  927#if defined(HAVE_TPETRA_MPI) 
  928    maybeInitializeRoots();
 
  938#if defined(HAVE_TPETRA_MPI) 
  941#if defined(HAVE_TPETRACORE_MPI_ADVANCE) 
  942  sendTypes.push_back(
"MpiAdvanceAlltoall");
 
  943  sendTypes.push_back(
"MpiAdvanceNbralltoallv");
 
 
  949  Teuchos::Array<EDistributorSendType> 
res;
 
  950  res.push_back(DISTRIBUTOR_ISEND);
 
  951  res.push_back(DISTRIBUTOR_SEND);
 
  952  res.push_back(DISTRIBUTOR_ALLTOALL);
 
  953#if defined(HAVE_TPETRA_MPI) 
  956#if defined(HAVE_TPETRACORE_MPI_ADVANCE) 
 
  963Teuchos::RCP<const Teuchos::ParameterList>
 
  964DistributorPlan::getValidParameters()
 const {
 
  965  using Teuchos::Array;
 
  966  using Teuchos::ParameterList;
 
  967  using Teuchos::parameterList;
 
  969  using Teuchos::setStringToIntegralParameter;
 
  974  const std::string 
validatedSendType = validSendTypeOrThrow(Behavior::defaultSendType());
 
  980                                                              "When using MPI, the variant of send to use in " 
  981                                                              "do[Reverse]Posts()",
 
  983  plist->set(
"Timer Label", 
"", 
"Label for Time Monitor output");
 
  985  return Teuchos::rcp_const_cast<const ParameterList>(
plist);
 
  988#if defined(HAVE_TPETRACORE_MPI_ADVANCE) 
  991struct MpixCommDeallocator {
 
  992  void free(MPIX_Comm** comm)
 const {
 
  993    MPIX_Comm_free(comm);
 
  997void DistributorPlan::initializeMpiAdvance() {
 
  999  TEUCHOS_ASSERT(mpixComm_.is_null());
 
 1002  Teuchos::RCP<const Teuchos::MpiComm<int>> mpiComm            = Teuchos::rcp_dynamic_cast<const Teuchos::MpiComm<int>>(comm_);
 
 1003  Teuchos::RCP<const Teuchos::OpaqueWrapper<MPI_Comm>> rawComm = mpiComm->getRawMpiComm();
 
 1005  if (sendType_ == DISTRIBUTOR_MPIADVANCE_ALLTOALL ||
 
 1006      sendType_ == DISTRIBUTOR_MPIADVANCE_NBRALLTOALLV) {
 
 1007    MPIX_Comm** mpixComm = 
new (MPIX_Comm*);
 
 1008    err                  = MPIX_Comm_init(mpixComm, (*rawComm)());
 
 1009    mpixComm_            = Teuchos::RCP(mpixComm,
 
 1010                                        MpixCommDeallocator(),
 
 1015  TEUCHOS_ASSERT(err == 0);
 
 1019#if defined(HAVE_TPETRA_MPI) 
 1021void DistributorPlan::maybeInitializeRoots() {
 
 1023  if (DISTRIBUTOR_IALLTOFEWV != sendType_) {
 
 1028  ProfilingRegion region_maybeInitializeRoots(
"Tpetra::DistributorPlan::maybeInitializeRoots");
 
 1031  const int numRecvs = (int)(getNumReceives() + (hasSelfMessage() ? 1 : 0));
 
 1032  std::vector<int> sendbuf(comm_->getSize(), numRecvs);
 
 1033  std::vector<int> recvbuf(comm_->getSize());
 
 1038  Teuchos::RCP<const Teuchos::MpiComm<int>> mpiComm            = Teuchos::rcp_dynamic_cast<const Teuchos::MpiComm<int>>(comm_);
 
 1039  Teuchos::RCP<const Teuchos::OpaqueWrapper<MPI_Comm>> rawComm = mpiComm->getRawMpiComm();
 
 1040  MPI_Comm comm                                                = (*rawComm)();
 
 1041  MPI_Alltoall(sendbuf.data(), 1, MPI_INT, recvbuf.data(), 1, MPI_INT, comm);
 
 1044  for (
size_t root = 0; root < recvbuf.size(); ++root) {
 
 1045    if (recvbuf[root] > 0) {
 
 1046      roots_.push_back(root);
 
 1052  int slow = !getIndicesTo().is_null() ? 1 : 0;
 
 1053  MPI_Allreduce(MPI_IN_PLACE, &slow, 1, MPI_INT, MPI_LOR, comm);
 
 1057        std::stringstream ss;
 
 1058        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;
 
 1059        std::cerr << ss.str();
 
 1064    sendType_ = DISTRIBUTOR_SEND;
 
 1070  if (roots_.size() * roots_.size() >= 
size_t(comm_->getSize())) {
 
 1072      std::stringstream ss;
 
 1073      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;
 
 1074      std::cerr << ss.str();
 
 1077    sendType_ = DISTRIBUTOR_SEND;
 
 1082DistributorPlan::SubViewLimits DistributorPlan::getImportViewLimits(
size_t numPackets)
 const {
 
 1083  const size_t actualNumReceives = getNumReceives() + (hasSelfMessage() ? 1 : 0);
 
 1085  IndexView importStarts(actualNumReceives);
 
 1086  IndexView importLengths(actualNumReceives);
 
 1089  for (
size_t i = 0; i < actualNumReceives; ++i) {
 
 1090    importStarts[i] = offset;
 
 1091    offset += getLengthsFrom()[i] * numPackets;
 
 1092    importLengths[i] = getLengthsFrom()[i] * numPackets;
 
 1094  return std::make_pair(importStarts, importLengths);
 
 1097DistributorPlan::SubViewLimits DistributorPlan::getImportViewLimits(
const Teuchos::ArrayView<const size_t>& numImportPacketsPerLID)
 const {
 
 1098  const size_t actualNumReceives = getNumReceives() + (hasSelfMessage() ? 1 : 0);
 
 1100  IndexView importStarts(actualNumReceives);
 
 1101  IndexView importLengths(actualNumReceives);
 
 1104  size_t curLIDoffset = 0;
 
 1105  for (
size_t i = 0; i < actualNumReceives; ++i) {
 
 1106    size_t totalPacketsFrom_i = 0;
 
 1107    for (
size_t j = 0; j < getLengthsFrom()[i]; ++j) {
 
 1108      totalPacketsFrom_i += numImportPacketsPerLID[curLIDoffset + j];
 
 1110    curLIDoffset += getLengthsFrom()[i];
 
 1111    importStarts[i] = offset;
 
 1112    offset += totalPacketsFrom_i;
 
 1113    importLengths[i] = totalPacketsFrom_i;
 
 1115  return std::make_pair(importStarts, importLengths);
 
 1118DistributorPlan::SubViewLimits DistributorPlan::getExportViewLimits(
size_t numPackets)
 const {
 
 1119  if (getIndicesTo().is_null()) {
 
 1120    const size_t actualNumSends = getNumSends() + (hasSelfMessage() ? 1 : 0);
 
 1121    IndexView exportStarts(actualNumSends);
 
 1122    IndexView exportLengths(actualNumSends);
 
 1123    for (
size_t pp = 0; pp < actualNumSends; ++pp) {
 
 1124      exportStarts[pp]  = getStartsTo()[pp] * numPackets;
 
 1125      exportLengths[pp] = getLengthsTo()[pp] * numPackets;
 
 1127    return std::make_pair(exportStarts, exportLengths);
 
 1129    const size_t numIndices = getIndicesTo().size();
 
 1130    IndexView exportStarts(numIndices);
 
 1131    IndexView exportLengths(numIndices);
 
 1132    for (
size_t j = 0; j < numIndices; ++j) {
 
 1133      exportStarts[j]  = getIndicesTo()[j] * numPackets;
 
 1134      exportLengths[j] = numPackets;
 
 1136    return std::make_pair(exportStarts, exportLengths);
 
 1140DistributorPlan::SubViewLimits DistributorPlan::getExportViewLimits(
const Teuchos::ArrayView<const size_t>& numExportPacketsPerLID)
 const {
 
 1141  if (getIndicesTo().is_null()) {
 
 1142    const size_t actualNumSends = getNumSends() + (hasSelfMessage() ? 1 : 0);
 
 1143    IndexView exportStarts(actualNumSends);
 
 1144    IndexView exportLengths(actualNumSends);
 
 1146    for (
size_t pp = 0; pp < actualNumSends; ++pp) {
 
 1147      size_t numPackets = 0;
 
 1148      for (
size_t j = getStartsTo()[pp];
 
 1149           j < getStartsTo()[pp] + getLengthsTo()[pp]; ++j) {
 
 1150        numPackets += numExportPacketsPerLID[j];
 
 1152      exportStarts[pp] = offset;
 
 1153      offset += numPackets;
 
 1154      exportLengths[pp] = numPackets;
 
 1156    return std::make_pair(exportStarts, exportLengths);
 
 1158    const size_t numIndices = getIndicesTo().size();
 
 1159    IndexView exportStarts(numIndices);
 
 1160    IndexView exportLengths(numIndices);
 
 1162    for (
size_t j = 0; j < numIndices; ++j) {
 
 1163      exportStarts[j] = offset;
 
 1164      offset += numExportPacketsPerLID[j];
 
 1165      exportLengths[j] = numExportPacketsPerLID[j];
 
 1167    return std::make_pair(exportStarts, exportLengths);
 
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.
 
Declaration of Tpetra::Details::Profiling, a scope guard for Kokkos Profiling.
 
Stand-alone utility functions and macros.
 
#define TPETRA_EFFICIENCY_WARNING(throw_exception_test, msg)
Print or throw an efficency warning.
 
#define SHARED_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg, comm)
Test for exception, with reduction over the given communicator.
 
Struct that holds views of the contents of a CrsMatrix.
 
static bool debug()
Whether Tpetra is in debug mode.
 
static bool verbose()
Whether Tpetra is in verbose mode.
 
Implementation details of Tpetra.
 
std::string DistributorSendTypeEnumToString(EDistributorSendType sendType)
Convert an EDistributorSendType enum value to a string.
 
EDistributorSendType
The type of MPI send that Distributor should use.
 
Teuchos::Array< EDistributorSendType > distributorSendTypeEnums()
Valid enum values of distributor send types.
 
Teuchos::Array< std::string > distributorSendTypes()
Valid string values for Distributor's "Send type" parameter.
 
EDistributorSendType DistributorSendTypeStringToEnum(const std::string_view s)
Convert a string to an EDistributorSendType. Throw on error.
 
const std::string & validSendTypeOrThrow(const std::string &s)
Valid enum values of distributor send types.
 
EDistributorHowInitialized
Enum indicating how and whether a Distributor was initialized.
 
std::string DistributorHowInitializedEnumToString(EDistributorHowInitialized how)
Convert an EDistributorHowInitialized enum value to a string.
 
Namespace Tpetra contains the class and methods constituting the Tpetra library.
 
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.