40 #ifndef TPETRA_DETAILS_DISTRIBUTOR_ACTOR_HPP
41 #define TPETRA_DETAILS_DISTRIBUTOR_ACTOR_HPP
43 #include "Tpetra_Details_DistributorPlan.hpp"
46 #include "Teuchos_Array.hpp"
47 #include "Teuchos_Comm.hpp"
48 #include "Teuchos_RCP.hpp"
49 #include "Teuchos_Time.hpp"
54 template <
class View1,
class View2>
55 constexpr
bool areKokkosViews = Kokkos::Impl::is_view<View1>::value && Kokkos::Impl::is_view<View2>::value;
57 class DistributorActor {
61 DistributorActor(
const DistributorActor& otherActor);
63 template <
class ExpView,
class ImpView>
64 void doPostsAndWaits(
const DistributorPlan& plan,
65 const ExpView &exports,
67 const ImpView &imports);
69 template <
class ExpView,
class ImpView>
70 void doPostsAndWaits(
const DistributorPlan& plan,
71 const ExpView &exports,
72 const Teuchos::ArrayView<const size_t>& numExportPacketsPerLID,
73 const ImpView &imports,
74 const Teuchos::ArrayView<const size_t>& numImportPacketsPerLID);
76 template <
class ExpView,
class ImpView>
77 void doPosts(
const DistributorPlan& plan,
78 const ExpView& exports,
80 const ImpView& imports);
82 template <
class ExpView,
class ImpView>
83 void doPosts(
const DistributorPlan& plan,
84 const ExpView &exports,
85 const Teuchos::ArrayView<const size_t>& numExportPacketsPerLID,
86 const ImpView &imports,
87 const Teuchos::ArrayView<const size_t>& numImportPacketsPerLID);
89 void doWaits(
const DistributorPlan& plan);
92 Teuchos::Array<Teuchos::RCP<Teuchos::CommRequest<int>>> requests_;
94 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
95 Teuchos::RCP<Teuchos::Time> timer_doPosts3KV_;
96 Teuchos::RCP<Teuchos::Time> timer_doPosts4KV_;
97 Teuchos::RCP<Teuchos::Time> timer_doWaits_;
98 Teuchos::RCP<Teuchos::Time> timer_doPosts3KV_recvs_;
99 Teuchos::RCP<Teuchos::Time> timer_doPosts4KV_recvs_;
100 Teuchos::RCP<Teuchos::Time> timer_doPosts3KV_barrier_;
101 Teuchos::RCP<Teuchos::Time> timer_doPosts4KV_barrier_;
102 Teuchos::RCP<Teuchos::Time> timer_doPosts3KV_sends_;
103 Teuchos::RCP<Teuchos::Time> timer_doPosts4KV_sends_;
104 Teuchos::RCP<Teuchos::Time> timer_doPosts3KV_sends_slow_;
105 Teuchos::RCP<Teuchos::Time> timer_doPosts4KV_sends_slow_;
106 Teuchos::RCP<Teuchos::Time> timer_doPosts3KV_sends_fast_;
107 Teuchos::RCP<Teuchos::Time> timer_doPosts4KV_sends_fast_;
114 template <
class ExpView,
class ImpView>
115 void DistributorActor::doPostsAndWaits(
const DistributorPlan& plan,
116 const ExpView& exports,
118 const ImpView& imports)
120 static_assert(areKokkosViews<ExpView, ImpView>,
121 "Data arrays for DistributorActor::doPostsAndWaits must be Kokkos::Views");
122 doPosts(plan, exports, numPackets, imports);
126 template <
class ExpView,
class ImpView>
127 void DistributorActor::doPostsAndWaits(
const DistributorPlan& plan,
128 const ExpView& exports,
129 const Teuchos::ArrayView<const size_t>& numExportPacketsPerLID,
130 const ImpView& imports,
131 const Teuchos::ArrayView<const size_t>& numImportPacketsPerLID)
133 static_assert(areKokkosViews<ExpView, ImpView>,
134 "Data arrays for DistributorActor::doPostsAndWaits must be Kokkos::Views");
135 doPosts(plan, exports, numExportPacketsPerLID, imports, numImportPacketsPerLID);
139 template <
class ExpView,
class ImpView>
140 void DistributorActor::doPosts(
const DistributorPlan& plan,
141 const ExpView& exports,
143 const ImpView& imports)
145 static_assert(areKokkosViews<ExpView, ImpView>,
146 "Data arrays for DistributorActor::doPosts must be Kokkos::Views");
147 using Teuchos::Array;
149 using Teuchos::FancyOStream;
150 using Teuchos::includesVerbLevel;
151 using Teuchos::ireceive;
152 using Teuchos::isend;
153 using Teuchos::readySend;
155 using Teuchos::ssend;
156 using Teuchos::TypeNameTraits;
157 using Teuchos::typeName;
159 using Kokkos::Compat::create_const_view;
160 using Kokkos::Compat::create_view;
161 using Kokkos::Compat::subview_offset;
162 using Kokkos::Compat::deep_copy_offset;
163 typedef Array<size_t>::size_type size_type;
164 typedef ExpView exports_view_type;
165 typedef ImpView imports_view_type;
167 #ifdef KOKKOS_ENABLE_CUDA
169 (! std::is_same<typename ExpView::memory_space, Kokkos::CudaUVMSpace>::value &&
170 ! std::is_same<typename ImpView::memory_space, Kokkos::CudaUVMSpace>::value,
171 "Please do not use Tpetra::Distributor with UVM allocations. "
172 "See Trilinos GitHub issue #1088.");
175 #ifdef KOKKOS_ENABLE_SYCL
177 (! std::is_same<typename ExpView::memory_space, Kokkos::Experimental::SYCLSharedUSMSpace>::value &&
178 ! std::is_same<typename ImpView::memory_space, Kokkos::Experimental::SYCLSharedUSMSpace>::value,
179 "Please do not use Tpetra::Distributor with SharedUSM allocations. "
180 "See Trilinos GitHub issue #1088 (corresponding to CUDA).");
183 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
184 Teuchos::TimeMonitor timeMon (*timer_doPosts3KV_);
187 const int myRank = plan.getComm()->getRank ();
191 const bool doBarrier = plan.barrierBetweenRecvSend();
193 TEUCHOS_TEST_FOR_EXCEPTION(
194 sendType == Details::DISTRIBUTOR_RSEND && ! doBarrier,
196 "Tpetra::Distributor::doPosts(3 args, Kokkos): Ready-send version "
197 "requires a barrier between posting receives and posting ready sends. "
198 "This should have been checked before. "
199 "Please report this bug to the Tpetra developers.");
201 size_t selfReceiveOffset = 0;
208 const int pathTag = 0;
209 const int tag = plan.getTag(pathTag);
211 #ifdef HAVE_TPETRA_DEBUG
212 TEUCHOS_TEST_FOR_EXCEPTION
213 (requests_.size () != 0,
215 "Tpetra::Distributor::doPosts(3 args, Kokkos): Process "
216 << myRank <<
": requests_.size() = " << requests_.size () <<
" != 0.");
232 const size_type actualNumReceives = as<size_type> (plan.getNumReceives()) +
233 as<size_type> (plan.hasSelfMessage() ? 1 : 0);
234 requests_.resize (0);
242 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
243 Teuchos::TimeMonitor timeMonRecvs (*timer_doPosts3KV_recvs_);
246 size_t curBufferOffset = 0;
247 for (size_type i = 0; i < actualNumReceives; ++i) {
248 const size_t curBufLen = plan.getLengthsFrom()[i] * numPackets;
249 if (plan.getProcsFrom()[i] != myRank) {
257 TEUCHOS_TEST_FOR_EXCEPTION(
258 curBufferOffset + curBufLen >
static_cast<size_t> (imports.size ()),
259 std::logic_error,
"Tpetra::Distributor::doPosts(3 args, Kokkos): "
260 "Exceeded size of 'imports' array in packing loop on Process " <<
261 myRank <<
". imports.size() = " << imports.size () <<
" < "
262 "curBufferOffset(" << curBufferOffset <<
") + curBufLen(" <<
264 imports_view_type recvBuf =
265 subview_offset (imports, curBufferOffset, curBufLen);
266 requests_.push_back (ireceive<int> (recvBuf, plan.getProcsFrom()[i],
267 tag, *plan.getComm()));
270 selfReceiveOffset = curBufferOffset;
272 curBufferOffset += curBufLen;
277 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
278 Teuchos::TimeMonitor timeMonBarrier (*timer_doPosts3KV_barrier_);
286 plan.getComm()->barrier ();
289 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
290 Teuchos::TimeMonitor timeMonSends (*timer_doPosts3KV_sends_);
298 size_t numBlocks = plan.getNumSends() + plan.hasSelfMessage();
299 size_t procIndex = 0;
300 while ((procIndex < numBlocks) && (plan.getProcsTo()[procIndex] < myRank)) {
303 if (procIndex == numBlocks) {
308 size_t selfIndex = 0;
310 if (plan.getIndicesTo().is_null()) {
312 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
313 Teuchos::TimeMonitor timeMonSends2 (*timer_doPosts3KV_sends_fast_);
318 for (
size_t i = 0; i < numBlocks; ++i) {
319 size_t p = i + procIndex;
320 if (p > (numBlocks - 1)) {
324 if (plan.getProcsTo()[p] != myRank) {
325 exports_view_type tmpSend = subview_offset(
326 exports, plan.getStartsTo()[p]*numPackets, plan.getLengthsTo()[p]*numPackets);
328 if (sendType == Details::DISTRIBUTOR_SEND) {
330 as<int> (tmpSend.size ()),
331 plan.getProcsTo()[p], tag, *plan.getComm());
333 else if (sendType == Details::DISTRIBUTOR_ISEND) {
334 exports_view_type tmpSendBuf =
335 subview_offset (exports, plan.getStartsTo()[p] * numPackets,
336 plan.getLengthsTo()[p] * numPackets);
337 requests_.push_back (isend<int> (tmpSendBuf, plan.getProcsTo()[p],
338 tag, *plan.getComm()));
340 else if (sendType == Details::DISTRIBUTOR_RSEND) {
341 readySend<int> (tmpSend,
342 as<int> (tmpSend.size ()),
343 plan.getProcsTo()[p], tag, *plan.getComm());
345 else if (sendType == Details::DISTRIBUTOR_SSEND) {
347 as<int> (tmpSend.size ()),
348 plan.getProcsTo()[p], tag, *plan.getComm());
350 TEUCHOS_TEST_FOR_EXCEPTION(
353 "Tpetra::Distributor::doPosts(3 args, Kokkos): "
354 "Invalid send type. We should never get here. "
355 "Please report this bug to the Tpetra developers.");
363 if (plan.hasSelfMessage()) {
371 deep_copy_offset(imports, exports, selfReceiveOffset,
372 plan.getStartsTo()[selfNum]*numPackets,
373 plan.getLengthsTo()[selfNum]*numPackets);
378 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
379 Teuchos::TimeMonitor timeMonSends2 (*timer_doPosts3KV_sends_slow_);
382 typedef typename ExpView::non_const_value_type Packet;
383 typedef typename ExpView::array_layout Layout;
384 typedef typename ExpView::device_type Device;
385 typedef typename ExpView::memory_traits Mem;
386 Kokkos::View<Packet*,Layout,Device,Mem> sendArray (
"sendArray",
387 plan.getMaxSendLength() * numPackets);
391 TEUCHOS_TEST_FOR_EXCEPTION(
392 sendType == Details::DISTRIBUTOR_ISEND,
394 "Tpetra::Distributor::doPosts(3 args, Kokkos): The \"send buffer\" code path "
395 "doesn't currently work with nonblocking sends.");
397 for (
size_t i = 0; i < numBlocks; ++i) {
398 size_t p = i + procIndex;
399 if (p > (numBlocks - 1)) {
403 if (plan.getProcsTo()[p] != myRank) {
404 size_t sendArrayOffset = 0;
405 size_t j = plan.getStartsTo()[p];
406 for (
size_t k = 0; k < plan.getLengthsTo()[p]; ++k, ++j) {
407 deep_copy_offset(sendArray, exports, sendArrayOffset,
408 plan.getIndicesTo()[j]*numPackets, numPackets);
409 sendArrayOffset += numPackets;
412 subview_offset(sendArray,
size_t(0), plan.getLengthsTo()[p]*numPackets);
414 if (sendType == Details::DISTRIBUTOR_SEND) {
416 as<int> (tmpSend.size ()),
417 plan.getProcsTo()[p], tag, *plan.getComm());
419 else if (sendType == Details::DISTRIBUTOR_ISEND) {
420 exports_view_type tmpSendBuf =
421 subview_offset (sendArray,
size_t(0), plan.getLengthsTo()[p] * numPackets);
422 requests_.push_back (isend<int> (tmpSendBuf, plan.getProcsTo()[p],
423 tag, *plan.getComm()));
425 else if (sendType == Details::DISTRIBUTOR_RSEND) {
426 readySend<int> (tmpSend,
427 as<int> (tmpSend.size ()),
428 plan.getProcsTo()[p], tag, *plan.getComm());
430 else if (sendType == Details::DISTRIBUTOR_SSEND) {
432 as<int> (tmpSend.size ()),
433 plan.getProcsTo()[p], tag, *plan.getComm());
436 TEUCHOS_TEST_FOR_EXCEPTION(
439 "Tpetra::Distributor::doPosts(3 args, Kokkos): "
440 "Invalid send type. We should never get here. "
441 "Please report this bug to the Tpetra developers.");
446 selfIndex = plan.getStartsTo()[p];
450 if (plan.hasSelfMessage()) {
451 for (
size_t k = 0; k < plan.getLengthsTo()[selfNum]; ++k) {
452 deep_copy_offset(imports, exports, selfReceiveOffset,
453 plan.getIndicesTo()[selfIndex]*numPackets, numPackets);
455 selfReceiveOffset += numPackets;
461 template <
class ExpView,
class ImpView>
462 void DistributorActor::doPosts(
const DistributorPlan& plan,
463 const ExpView &exports,
464 const Teuchos::ArrayView<const size_t>& numExportPacketsPerLID,
465 const ImpView &imports,
466 const Teuchos::ArrayView<const size_t>& numImportPacketsPerLID)
468 static_assert(areKokkosViews<ExpView, ImpView>,
469 "Data arrays for DistributorActor::doPosts must be Kokkos::Views");
470 using Teuchos::Array;
472 using Teuchos::ireceive;
473 using Teuchos::isend;
474 using Teuchos::readySend;
476 using Teuchos::ssend;
477 using Teuchos::TypeNameTraits;
479 using Kokkos::Compat::create_const_view;
480 using Kokkos::Compat::create_view;
481 using Kokkos::Compat::subview_offset;
482 using Kokkos::Compat::deep_copy_offset;
483 typedef Array<size_t>::size_type size_type;
484 typedef ExpView exports_view_type;
485 typedef ImpView imports_view_type;
487 #ifdef KOKKOS_ENABLE_CUDA
488 static_assert (! std::is_same<typename ExpView::memory_space, Kokkos::CudaUVMSpace>::value &&
489 ! std::is_same<typename ImpView::memory_space, Kokkos::CudaUVMSpace>::value,
490 "Please do not use Tpetra::Distributor with UVM "
491 "allocations. See GitHub issue #1088.");
494 #ifdef KOKKOS_ENABLE_SYCL
495 static_assert (! std::is_same<typename ExpView::memory_space, Kokkos::Experimental::SYCLSharedUSMSpace>::value &&
496 ! std::is_same<typename ImpView::memory_space, Kokkos::Experimental::SYCLSharedUSMSpace>::value,
497 "Please do not use Tpetra::Distributor with SharedUSM "
498 "allocations. See GitHub issue #1088 (corresponding to CUDA).");
501 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
502 Teuchos::TimeMonitor timeMon (*timer_doPosts4KV_);
508 const bool doBarrier = plan.barrierBetweenRecvSend();
510 TEUCHOS_TEST_FOR_EXCEPTION(
511 sendType == Details::DISTRIBUTOR_RSEND && ! doBarrier,
512 std::logic_error,
"Tpetra::Distributor::doPosts(4 args, Kokkos): Ready-send "
513 "version requires a barrier between posting receives and posting ready "
514 "sends. This should have been checked before. "
515 "Please report this bug to the Tpetra developers.");
517 const int myProcID = plan.getComm()->getRank ();
518 size_t selfReceiveOffset = 0;
520 #ifdef HAVE_TEUCHOS_DEBUG
522 size_t totalNumImportPackets = 0;
523 for (size_type ii = 0; ii < numImportPacketsPerLID.size (); ++ii) {
524 totalNumImportPackets += numImportPacketsPerLID[ii];
526 TEUCHOS_TEST_FOR_EXCEPTION(
527 imports.extent (0) < totalNumImportPackets, std::runtime_error,
528 "Tpetra::Distributor::doPosts(4 args, Kokkos): The 'imports' array must have "
529 "enough entries to hold the expected number of import packets. "
530 "imports.extent(0) = " << imports.extent (0) <<
" < "
531 "totalNumImportPackets = " << totalNumImportPackets <<
".");
539 const int pathTag = 1;
540 const int tag = plan.getTag(pathTag);
542 #ifdef HAVE_TEUCHOS_DEBUG
543 TEUCHOS_TEST_FOR_EXCEPTION
544 (requests_.size () != 0, std::logic_error,
"Tpetra::Distributor::"
545 "doPosts(4 args, Kokkos): Process " << myProcID <<
": requests_.size () = "
546 << requests_.size () <<
" != 0.");
561 const size_type actualNumReceives = as<size_type> (plan.getNumReceives()) +
562 as<size_type> (plan.hasSelfMessage() ? 1 : 0);
563 requests_.resize (0);
571 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
572 Teuchos::TimeMonitor timeMonRecvs (*timer_doPosts4KV_recvs_);
575 size_t curBufferOffset = 0;
576 size_t curLIDoffset = 0;
577 for (size_type i = 0; i < actualNumReceives; ++i) {
578 size_t totalPacketsFrom_i = 0;
579 for (
size_t j = 0; j < plan.getLengthsFrom()[i]; ++j) {
580 totalPacketsFrom_i += numImportPacketsPerLID[curLIDoffset+j];
582 curLIDoffset += plan.getLengthsFrom()[i];
583 if (plan.getProcsFrom()[i] != myProcID && totalPacketsFrom_i) {
592 imports_view_type recvBuf =
593 subview_offset (imports, curBufferOffset, totalPacketsFrom_i);
594 requests_.push_back (ireceive<int> (recvBuf, plan.getProcsFrom()[i],
595 tag, *plan.getComm()));
598 selfReceiveOffset = curBufferOffset;
600 curBufferOffset += totalPacketsFrom_i;
605 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
606 Teuchos::TimeMonitor timeMonBarrier (*timer_doPosts4KV_barrier_);
613 plan.getComm()->barrier ();
616 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
617 Teuchos::TimeMonitor timeMonSends (*timer_doPosts4KV_sends_);
622 Array<size_t> sendPacketOffsets(plan.getNumSends(),0), packetsPerSend(plan.getNumSends(),0);
623 size_t maxNumPackets = 0;
624 size_t curPKToffset = 0;
625 for (
size_t pp=0; pp<plan.getNumSends(); ++pp) {
626 sendPacketOffsets[pp] = curPKToffset;
627 size_t numPackets = 0;
628 for (
size_t j=plan.getStartsTo()[pp]; j<plan.getStartsTo()[pp]+plan.getLengthsTo()[pp]; ++j) {
629 numPackets += numExportPacketsPerLID[j];
631 if (numPackets > maxNumPackets) maxNumPackets = numPackets;
632 packetsPerSend[pp] = numPackets;
633 curPKToffset += numPackets;
638 size_t numBlocks = plan.getNumSends() + plan.hasSelfMessage();
639 size_t procIndex = 0;
640 while ((procIndex < numBlocks) && (plan.getProcsTo()[procIndex] < myProcID)) {
643 if (procIndex == numBlocks) {
648 size_t selfIndex = 0;
649 if (plan.getIndicesTo().is_null()) {
651 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
652 Teuchos::TimeMonitor timeMonSends2 (*timer_doPosts4KV_sends_fast_);
657 for (
size_t i = 0; i < numBlocks; ++i) {
658 size_t p = i + procIndex;
659 if (p > (numBlocks - 1)) {
663 if (plan.getProcsTo()[p] != myProcID && packetsPerSend[p] > 0) {
664 exports_view_type tmpSend =
665 subview_offset(exports, sendPacketOffsets[p], packetsPerSend[p]);
667 if (sendType == Details::DISTRIBUTOR_SEND) {
669 as<int> (tmpSend.size ()),
670 plan.getProcsTo()[p], tag, *plan.getComm());
672 else if (sendType == Details::DISTRIBUTOR_RSEND) {
673 readySend<int> (tmpSend,
674 as<int> (tmpSend.size ()),
675 plan.getProcsTo()[p], tag, *plan.getComm());
677 else if (sendType == Details::DISTRIBUTOR_ISEND) {
678 exports_view_type tmpSendBuf =
679 subview_offset (exports, sendPacketOffsets[p], packetsPerSend[p]);
680 requests_.push_back (isend<int> (tmpSendBuf, plan.getProcsTo()[p],
681 tag, *plan.getComm()));
683 else if (sendType == Details::DISTRIBUTOR_SSEND) {
685 as<int> (tmpSend.size ()),
686 plan.getProcsTo()[p], tag, *plan.getComm());
689 TEUCHOS_TEST_FOR_EXCEPTION(
690 true, std::logic_error,
691 "Tpetra::Distributor::doPosts(4 args, Kokkos): "
692 "Invalid send type. We should never get here. "
693 "Please report this bug to the Tpetra developers.");
701 if (plan.hasSelfMessage()) {
702 deep_copy_offset(imports, exports, selfReceiveOffset,
703 sendPacketOffsets[selfNum], packetsPerSend[selfNum]);
708 #ifdef HAVE_TPETRA_DISTRIBUTOR_TIMINGS
709 Teuchos::TimeMonitor timeMonSends2 (*timer_doPosts4KV_sends_slow_);
713 typedef typename ExpView::non_const_value_type Packet;
714 typedef typename ExpView::array_layout Layout;
715 typedef typename ExpView::device_type Device;
716 typedef typename ExpView::memory_traits Mem;
717 Kokkos::View<Packet*,Layout,Device,Mem> sendArray (
"sendArray", maxNumPackets);
719 TEUCHOS_TEST_FOR_EXCEPTION(
720 sendType == Details::DISTRIBUTOR_ISEND,
722 "Tpetra::Distributor::doPosts(4-arg, Kokkos): "
723 "The \"send buffer\" code path may not necessarily work with nonblocking sends.");
725 Array<size_t> indicesOffsets (numExportPacketsPerLID.size(), 0);
727 for (
int j=0; j<numExportPacketsPerLID.size(); ++j) {
728 indicesOffsets[j] = ioffset;
729 ioffset += numExportPacketsPerLID[j];
732 for (
size_t i = 0; i < numBlocks; ++i) {
733 size_t p = i + procIndex;
734 if (p > (numBlocks - 1)) {
738 if (plan.getProcsTo()[p] != myProcID) {
739 size_t sendArrayOffset = 0;
740 size_t j = plan.getStartsTo()[p];
741 size_t numPacketsTo_p = 0;
742 for (
size_t k = 0; k < plan.getLengthsTo()[p]; ++k, ++j) {
743 numPacketsTo_p += numExportPacketsPerLID[j];
744 deep_copy_offset(sendArray, exports, sendArrayOffset,
745 indicesOffsets[j], numExportPacketsPerLID[j]);
746 sendArrayOffset += numExportPacketsPerLID[j];
748 if (numPacketsTo_p > 0) {
750 subview_offset(sendArray,
size_t(0), numPacketsTo_p);
752 if (sendType == Details::DISTRIBUTOR_RSEND) {
753 readySend<int> (tmpSend,
754 as<int> (tmpSend.size ()),
755 plan.getProcsTo()[p], tag, *plan.getComm());
757 else if (sendType == Details::DISTRIBUTOR_ISEND) {
758 exports_view_type tmpSendBuf =
759 subview_offset (sendArray,
size_t(0), numPacketsTo_p);
760 requests_.push_back (isend<int> (tmpSendBuf, plan.getProcsTo()[p],
761 tag, *plan.getComm()));
763 else if (sendType == Details::DISTRIBUTOR_SSEND) {
765 as<int> (tmpSend.size ()),
766 plan.getProcsTo()[p], tag, *plan.getComm());
770 as<int> (tmpSend.size ()),
771 plan.getProcsTo()[p], tag, *plan.getComm());
777 selfIndex = plan.getStartsTo()[p];
781 if (plan.hasSelfMessage()) {
782 for (
size_t k = 0; k < plan.getLengthsTo()[selfNum]; ++k) {
783 deep_copy_offset(imports, exports, selfReceiveOffset,
784 indicesOffsets[selfIndex],
785 numExportPacketsPerLID[selfIndex]);
786 selfReceiveOffset += numExportPacketsPerLID[selfIndex];
Stand-alone utility functions and macros.
Implementation details of Tpetra.
EDistributorSendType
The type of MPI send that Distributor should use.
Namespace Tpetra contains the class and methods constituting the Tpetra library.