25 #ifndef TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_DIST_OBJECT_KERNELS_HPP
26 #define TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_DIST_OBJECT_KERNELS_HPP
28 #include "Kokkos_Core.hpp"
29 #include "Kokkos_ArithTraits.hpp"
34 namespace KokkosRefactor {
50 template<
class IntegerType,
51 const bool isSigned = std::numeric_limits<IntegerType>::is_signed>
53 static KOKKOS_INLINE_FUNCTION
bool
54 test (
const IntegerType x,
55 const IntegerType exclusiveUpperBound);
59 template<
class IntegerType>
61 static KOKKOS_INLINE_FUNCTION
bool
62 test (
const IntegerType x,
63 const IntegerType exclusiveUpperBound)
65 return x < static_cast<IntegerType> (0) || x >= exclusiveUpperBound;
70 template<
class IntegerType>
71 struct OutOfBounds<IntegerType, false> {
72 static KOKKOS_INLINE_FUNCTION
bool
73 test (
const IntegerType x,
74 const IntegerType exclusiveUpperBound)
76 return x >= exclusiveUpperBound;
82 template<
class IntegerType>
83 KOKKOS_INLINE_FUNCTION
bool
84 outOfBounds (
const IntegerType x,
const IntegerType exclusiveUpperBound)
94 template <
typename DstView,
typename SrcView,
typename IdxView,
95 typename Enabled =
void>
96 struct PackArraySingleColumn {
97 typedef typename DstView::execution_space execution_space;
98 typedef typename execution_space::size_type size_type;
105 PackArraySingleColumn (
const DstView& dst_,
109 dst(dst_), src(src_), idx(idx_), col(col_) {}
111 KOKKOS_INLINE_FUNCTION
void
112 operator() (
const size_type k)
const {
113 dst(k) = src(idx(k), col);
117 pack (
const DstView& dst,
121 const execution_space &space)
123 typedef Kokkos::RangePolicy<execution_space, size_type> range_type;
125 (
"Tpetra::MultiVector pack one col",
126 range_type (space, 0, idx.size ()),
127 PackArraySingleColumn (dst, src, idx, col));
131 template <
typename DstView,
134 typename SizeType =
typename DstView::execution_space::size_type,
135 typename Enabled =
void>
136 class PackArraySingleColumnWithBoundsCheck {
138 static_assert (Kokkos::is_view<DstView>::value,
139 "DstView must be a Kokkos::View.");
140 static_assert (Kokkos::is_view<SrcView>::value,
141 "SrcView must be a Kokkos::View.");
142 static_assert (Kokkos::is_view<IdxView>::value,
143 "IdxView must be a Kokkos::View.");
144 static_assert (static_cast<int> (DstView::rank) == 1,
145 "DstView must be a rank-1 Kokkos::View.");
146 static_assert (static_cast<int> (SrcView::rank) == 2,
147 "SrcView must be a rank-2 Kokkos::View.");
148 static_assert (static_cast<int> (IdxView::rank) == 1,
149 "IdxView must be a rank-1 Kokkos::View.");
150 static_assert (std::is_integral<SizeType>::value,
151 "SizeType must be a built-in integer type.");
153 using execution_space =
typename DstView::execution_space;
156 typedef SizeType size_type;
157 using value_type = size_t;
164 execution_space space;
167 PackArraySingleColumnWithBoundsCheck (
const DstView& dst_,
170 const size_type col_) :
171 dst (dst_), src (src_), idx (idx_), col (col_) {}
173 KOKKOS_INLINE_FUNCTION
void
174 operator() (
const size_type k, value_type& lclErrCount)
const {
175 using index_type =
typename IdxView::non_const_value_type;
177 const index_type lclRow = idx(k);
178 if (lclRow < static_cast<index_type> (0) ||
179 lclRow >= static_cast<index_type> (src.extent (0))) {
183 dst(k) = src(lclRow, col);
187 KOKKOS_INLINE_FUNCTION
188 void init (value_type& initialErrorCount)
const {
189 initialErrorCount = 0;
192 KOKKOS_INLINE_FUNCTION
void
193 join (value_type& dstErrorCount,
194 const value_type& srcErrorCount)
const
196 dstErrorCount += srcErrorCount;
200 pack (
const DstView& dst,
204 const execution_space &space)
206 typedef Kokkos::RangePolicy<execution_space, size_type> range_type;
207 typedef typename IdxView::non_const_value_type index_type;
209 size_t errorCount = 0;
210 Kokkos::parallel_reduce
211 (
"Tpetra::MultiVector pack one col debug only",
212 range_type (space, 0, idx.size ()),
213 PackArraySingleColumnWithBoundsCheck (dst, src, idx, col),
216 if (errorCount != 0) {
220 auto idx_h = Kokkos::create_mirror_view (idx);
225 std::vector<index_type> badIndices;
226 const size_type numInds = idx_h.extent (0);
227 for (size_type k = 0; k < numInds; ++k) {
228 if (idx_h(k) < static_cast<index_type> (0) ||
229 idx_h(k) >= static_cast<index_type> (src.extent (0))) {
230 badIndices.push_back (idx_h(k));
234 TEUCHOS_TEST_FOR_EXCEPTION
235 (errorCount != badIndices.size (), std::logic_error,
236 "PackArraySingleColumnWithBoundsCheck: errorCount = " << errorCount
237 <<
" != badIndices.size() = " << badIndices.size () <<
". This sho"
238 "uld never happen. Please report this to the Tpetra developers.");
240 std::ostringstream os;
241 os <<
"MultiVector single-column pack kernel had "
242 << badIndices.size () <<
" out-of bounds index/ices. "
244 for (
size_t k = 0; k < badIndices.size (); ++k) {
246 if (k + 1 < badIndices.size ()) {
251 throw std::runtime_error (os.str ());
257 template <
typename DstView,
typename SrcView,
typename IdxView>
259 pack_array_single_column (
const DstView& dst,
264 const typename DstView::execution_space &space)
266 static_assert (Kokkos::is_view<DstView>::value,
267 "DstView must be a Kokkos::View.");
268 static_assert (Kokkos::is_view<SrcView>::value,
269 "SrcView must be a Kokkos::View.");
270 static_assert (Kokkos::is_view<IdxView>::value,
271 "IdxView must be a Kokkos::View.");
272 static_assert (static_cast<int> (DstView::rank) == 1,
273 "DstView must be a rank-1 Kokkos::View.");
274 static_assert (static_cast<int> (SrcView::rank) == 2,
275 "SrcView must be a rank-2 Kokkos::View.");
276 static_assert (static_cast<int> (IdxView::rank) == 1,
277 "IdxView must be a rank-1 Kokkos::View.");
279 using execution_space =
typename DstView::execution_space;
281 static_assert (Kokkos::SpaceAccessibility<execution_space,
282 typename DstView::memory_space>::accessible,
283 "DstView not accessible from execution space");
284 static_assert (Kokkos::SpaceAccessibility<execution_space,
285 typename SrcView::memory_space>::accessible,
286 "SrcView not accessible from execution space");
287 static_assert (Kokkos::SpaceAccessibility<execution_space,
288 typename IdxView::memory_space>::accessible,
289 "IdxView not accessible from execution space");
292 typedef PackArraySingleColumnWithBoundsCheck<DstView,SrcView,IdxView> impl_type;
293 impl_type::pack (dst, src, idx, col, space);
296 typedef PackArraySingleColumn<DstView,SrcView,IdxView> impl_type;
297 impl_type::pack (dst, src, idx, col, space);
303 template <
typename DstView,
typename SrcView,
typename IdxView>
305 pack_array_single_column (
const DstView& dst,
309 const bool debug =
true)
311 pack_array_single_column(dst, src, idx, col, debug,
typename DstView::execution_space());
314 template <
typename DstView,
typename SrcView,
typename IdxView,
315 typename Enabled =
void>
316 struct PackArrayMultiColumn {
317 using execution_space =
typename DstView::execution_space;
318 typedef typename execution_space::size_type size_type;
325 PackArrayMultiColumn (
const DstView& dst_,
328 const size_t numCols_) :
329 dst(dst_), src(src_), idx(idx_), numCols(numCols_) {}
331 KOKKOS_INLINE_FUNCTION
void
332 operator() (
const size_type k)
const {
333 const typename IdxView::value_type localRow = idx(k);
334 const size_t offset = k*numCols;
335 for (
size_t j = 0; j < numCols; ++j) {
336 dst(offset + j) = src(localRow, j);
340 static void pack(
const DstView& dst,
344 const execution_space &space) {
345 typedef Kokkos::RangePolicy<execution_space, size_type> range_type;
347 (
"Tpetra::MultiVector pack multicol const stride",
348 range_type (space, 0, idx.size ()),
349 PackArrayMultiColumn (dst, src, idx, numCols));
353 template <
typename DstView,
356 typename SizeType =
typename DstView::execution_space::size_type,
357 typename Enabled =
void>
358 class PackArrayMultiColumnWithBoundsCheck {
360 using size_type = SizeType;
361 using value_type = size_t;
362 using execution_space =
typename DstView::execution_space;
371 PackArrayMultiColumnWithBoundsCheck (
const DstView& dst_,
374 const size_type numCols_) :
375 dst (dst_), src (src_), idx (idx_), numCols (numCols_) {}
377 KOKKOS_INLINE_FUNCTION
void
378 operator() (
const size_type k, value_type& lclErrorCount)
const {
379 typedef typename IdxView::non_const_value_type index_type;
381 const index_type lclRow = idx(k);
382 if (lclRow < static_cast<index_type> (0) ||
383 lclRow >= static_cast<index_type> (src.extent (0))) {
387 const size_type offset = k*numCols;
388 for (size_type j = 0; j < numCols; ++j) {
389 dst(offset + j) = src(lclRow, j);
394 KOKKOS_INLINE_FUNCTION
395 void init (value_type& initialErrorCount)
const {
396 initialErrorCount = 0;
399 KOKKOS_INLINE_FUNCTION
void
400 join (value_type& dstErrorCount,
401 const value_type& srcErrorCount)
const
403 dstErrorCount += srcErrorCount;
407 pack (
const DstView& dst,
410 const size_type numCols,
411 const execution_space &space)
413 typedef Kokkos::RangePolicy<execution_space, size_type> range_type;
414 typedef typename IdxView::non_const_value_type index_type;
416 size_t errorCount = 0;
417 Kokkos::parallel_reduce
418 (
"Tpetra::MultiVector pack multicol const stride debug only",
419 range_type (space, 0, idx.size ()),
420 PackArrayMultiColumnWithBoundsCheck (dst, src, idx, numCols),
422 if (errorCount != 0) {
426 auto idx_h = Kokkos::create_mirror_view (idx);
431 std::vector<index_type> badIndices;
432 const size_type numInds = idx_h.extent (0);
433 for (size_type k = 0; k < numInds; ++k) {
434 if (idx_h(k) < static_cast<index_type> (0) ||
435 idx_h(k) >= static_cast<index_type> (src.extent (0))) {
436 badIndices.push_back (idx_h(k));
440 TEUCHOS_TEST_FOR_EXCEPTION
441 (errorCount != badIndices.size (), std::logic_error,
442 "PackArrayMultiColumnWithBoundsCheck: errorCount = " << errorCount
443 <<
" != badIndices.size() = " << badIndices.size () <<
". This sho"
444 "uld never happen. Please report this to the Tpetra developers.");
446 std::ostringstream os;
447 os <<
"Tpetra::MultiVector multiple-column pack kernel had "
448 << badIndices.size () <<
" out-of bounds index/ices (errorCount = "
449 << errorCount <<
"): [";
450 for (
size_t k = 0; k < badIndices.size (); ++k) {
452 if (k + 1 < badIndices.size ()) {
457 throw std::runtime_error (os.str ());
463 template <
typename DstView,
467 pack_array_multi_column (
const DstView& dst,
470 const size_t numCols,
472 const typename DstView::execution_space &space)
474 static_assert (Kokkos::is_view<DstView>::value,
475 "DstView must be a Kokkos::View.");
476 static_assert (Kokkos::is_view<SrcView>::value,
477 "SrcView must be a Kokkos::View.");
478 static_assert (Kokkos::is_view<IdxView>::value,
479 "IdxView must be a Kokkos::View.");
480 static_assert (static_cast<int> (DstView::rank) == 1,
481 "DstView must be a rank-1 Kokkos::View.");
482 static_assert (static_cast<int> (SrcView::rank) == 2,
483 "SrcView must be a rank-2 Kokkos::View.");
484 static_assert (static_cast<int> (IdxView::rank) == 1,
485 "IdxView must be a rank-1 Kokkos::View.");
487 using execution_space =
typename DstView::execution_space;
489 static_assert (Kokkos::SpaceAccessibility<execution_space,
490 typename DstView::memory_space>::accessible,
491 "DstView not accessible from execution space");
492 static_assert (Kokkos::SpaceAccessibility<execution_space,
493 typename SrcView::memory_space>::accessible,
494 "SrcView not accessible from execution space");
495 static_assert (Kokkos::SpaceAccessibility<execution_space,
496 typename IdxView::memory_space>::accessible,
497 "IdxView not accessible from execution space");
500 typedef PackArrayMultiColumnWithBoundsCheck<DstView,
501 SrcView, IdxView> impl_type;
502 impl_type::pack (dst, src, idx, numCols, space);
505 typedef PackArrayMultiColumn<DstView, SrcView, IdxView> impl_type;
506 impl_type::pack (dst, src, idx, numCols, space);
510 template <
typename DstView,
514 pack_array_multi_column (
const DstView& dst,
517 const size_t numCols,
518 const bool debug =
true) {
519 pack_array_multi_column(dst, src, idx, numCols, debug,
typename DstView::execution_space());
522 template <
typename DstView,
typename SrcView,
typename IdxView,
523 typename ColView,
typename Enabled =
void>
524 struct PackArrayMultiColumnVariableStride {
525 using execution_space =
typename DstView::execution_space;
526 typedef typename execution_space::size_type size_type;
534 PackArrayMultiColumnVariableStride (
const DstView& dst_,
538 const size_t numCols_) :
539 dst(dst_), src(src_), idx(idx_), col(col_), numCols(numCols_) {}
541 KOKKOS_INLINE_FUNCTION
542 void operator() (
const size_type k)
const {
543 const typename IdxView::value_type localRow = idx(k);
544 const size_t offset = k*numCols;
545 for (
size_t j = 0; j < numCols; ++j) {
546 dst(offset + j) = src(localRow, col(j));
550 static void pack(
const DstView& dst,
555 const execution_space &space) {
556 typedef Kokkos::RangePolicy<execution_space, size_type> range_type;
558 (
"Tpetra::MultiVector pack multicol var stride",
559 range_type (space, 0, idx.size ()),
560 PackArrayMultiColumnVariableStride (dst, src, idx, col, numCols));
564 template <
typename DstView,
568 typename SizeType =
typename DstView::execution_space::size_type,
569 typename Enabled =
void>
570 class PackArrayMultiColumnVariableStrideWithBoundsCheck {
572 using size_type = SizeType;
573 using value_type = size_t;
574 using execution_space =
typename DstView::execution_space;
584 PackArrayMultiColumnVariableStrideWithBoundsCheck (
const DstView& dst_,
588 const size_type numCols_) :
589 dst (dst_), src (src_), idx (idx_), col (col_), numCols (numCols_) {}
591 KOKKOS_INLINE_FUNCTION
void
592 operator() (
const size_type k, value_type& lclErrorCount)
const {
593 typedef typename IdxView::non_const_value_type row_index_type;
594 typedef typename ColView::non_const_value_type col_index_type;
596 const row_index_type lclRow = idx(k);
597 if (lclRow < static_cast<row_index_type> (0) ||
598 lclRow >= static_cast<row_index_type> (src.extent (0))) {
602 const size_type offset = k*numCols;
603 for (size_type j = 0; j < numCols; ++j) {
604 const col_index_type lclCol = col(j);
605 if (Impl::outOfBounds<col_index_type> (lclCol, src.extent (1))) {
609 dst(offset + j) = src(lclRow, lclCol);
615 KOKKOS_INLINE_FUNCTION
void
616 init (value_type& initialErrorCount)
const {
617 initialErrorCount = 0;
620 KOKKOS_INLINE_FUNCTION
void
621 join (value_type& dstErrorCount,
622 const value_type& srcErrorCount)
const
624 dstErrorCount += srcErrorCount;
628 pack (
const DstView& dst,
632 const size_type numCols,
633 const execution_space &space)
635 using range_type = Kokkos::RangePolicy<execution_space, size_type>;
636 using row_index_type =
typename IdxView::non_const_value_type;
637 using col_index_type =
typename ColView::non_const_value_type;
639 size_t errorCount = 0;
640 Kokkos::parallel_reduce
641 (
"Tpetra::MultiVector pack multicol var stride debug only",
642 range_type (space, 0, idx.size ()),
643 PackArrayMultiColumnVariableStrideWithBoundsCheck (dst, src, idx,
646 if (errorCount != 0) {
647 constexpr
size_t maxNumBadIndicesToPrint = 100;
649 std::ostringstream os;
650 os <<
"Tpetra::MultiVector multicolumn variable stride pack kernel "
651 "found " << errorCount
652 <<
" error" << (errorCount != size_t (1) ?
"s" :
"") <<
". ";
657 auto idx_h = Kokkos::create_mirror_view (idx);
662 std::vector<row_index_type> badRows;
663 const size_type numRowInds = idx_h.extent (0);
664 for (size_type k = 0; k < numRowInds; ++k) {
665 if (Impl::outOfBounds<row_index_type> (idx_h(k), src.extent (0))) {
666 badRows.push_back (idx_h(k));
670 if (badRows.size () != 0) {
671 os << badRows.size () <<
" out-of-bounds row ind"
672 << (badRows.size () != size_t (1) ?
"ices" :
"ex");
673 if (badRows.size () <= maxNumBadIndicesToPrint) {
675 for (
size_t k = 0; k < badRows.size (); ++k) {
677 if (k + 1 < badRows.size ()) {
688 os <<
"No out-of-bounds row indices. ";
693 auto col_h = Kokkos::create_mirror_view (col);
698 std::vector<col_index_type> badCols;
699 const size_type numColInds = col_h.extent (0);
700 for (size_type k = 0; k < numColInds; ++k) {
701 if (Impl::outOfBounds<col_index_type> (col_h(k), src.extent (1))) {
702 badCols.push_back (col_h(k));
706 if (badCols.size () != 0) {
707 os << badCols.size () <<
" out-of-bounds column ind"
708 << (badCols.size () != size_t (1) ?
"ices" :
"ex");
709 if (badCols.size () <= maxNumBadIndicesToPrint) {
711 for (
size_t k = 0; k < badCols.size (); ++k) {
713 if (k + 1 < badCols.size ()) {
724 os <<
"No out-of-bounds column indices. ";
727 TEUCHOS_TEST_FOR_EXCEPTION
728 (errorCount != 0 && badRows.size () == 0 && badCols.size () == 0,
729 std::logic_error,
"Tpetra::MultiVector variable stride pack "
730 "kernel reports errorCount=" << errorCount <<
", but we failed "
731 "to find any bad rows or columns. This should never happen. "
732 "Please report this to the Tpetra developers.");
734 throw std::runtime_error (os.str ());
739 template <
typename DstView,
744 pack_array_multi_column_variable_stride (
const DstView& dst,
748 const size_t numCols,
750 const typename DstView::execution_space &space)
752 static_assert (Kokkos::is_view<DstView>::value,
753 "DstView must be a Kokkos::View.");
754 static_assert (Kokkos::is_view<SrcView>::value,
755 "SrcView must be a Kokkos::View.");
756 static_assert (Kokkos::is_view<IdxView>::value,
757 "IdxView must be a Kokkos::View.");
758 static_assert (Kokkos::is_view<ColView>::value,
759 "ColView must be a Kokkos::View.");
760 static_assert (static_cast<int> (DstView::rank) == 1,
761 "DstView must be a rank-1 Kokkos::View.");
762 static_assert (static_cast<int> (SrcView::rank) == 2,
763 "SrcView must be a rank-2 Kokkos::View.");
764 static_assert (static_cast<int> (IdxView::rank) == 1,
765 "IdxView must be a rank-1 Kokkos::View.");
766 static_assert (static_cast<int> (ColView::rank) == 1,
767 "ColView must be a rank-1 Kokkos::View.");
769 using execution_space =
typename DstView::execution_space;
771 static_assert (Kokkos::SpaceAccessibility<execution_space,
772 typename DstView::memory_space>::accessible,
773 "DstView not accessible from execution space");
774 static_assert (Kokkos::SpaceAccessibility<execution_space,
775 typename SrcView::memory_space>::accessible,
776 "SrcView not accessible from execution space");
777 static_assert (Kokkos::SpaceAccessibility<execution_space,
778 typename IdxView::memory_space>::accessible,
779 "IdxView not accessible from execution space");
782 typedef PackArrayMultiColumnVariableStrideWithBoundsCheck<DstView,
783 SrcView, IdxView, ColView> impl_type;
784 impl_type::pack (dst, src, idx, col, numCols, space);
787 typedef PackArrayMultiColumnVariableStride<DstView,
788 SrcView, IdxView, ColView> impl_type;
789 impl_type::pack (dst, src, idx, col, numCols, space);
793 template <
typename DstView,
798 pack_array_multi_column_variable_stride (
const DstView& dst,
802 const size_t numCols,
803 const bool debug =
true) {
804 pack_array_multi_column_variable_stride(dst, src, idx, col, numCols, debug,
805 typename DstView::execution_space());
810 struct atomic_tag {};
811 struct nonatomic_tag {};
815 KOKKOS_INLINE_FUNCTION
816 void operator() (atomic_tag, SC& dest,
const SC& src)
const {
817 Kokkos::atomic_add (&dest, src);
821 KOKKOS_INLINE_FUNCTION
822 void operator() (nonatomic_tag, SC& dest,
const SC& src)
const {
833 KOKKOS_INLINE_FUNCTION
834 void operator() (atomic_tag, SC& dest,
const SC& src)
const {
839 KOKKOS_INLINE_FUNCTION
840 void operator() (nonatomic_tag, SC& dest,
const SC& src)
const {
846 template <
class Scalar>
850 KOKKOS_FUNCTION AbsMaxHelper& operator+=(AbsMaxHelper
const& rhs) {
851 auto lhs_abs_value = Kokkos::ArithTraits<Scalar>::abs(value);
852 auto rhs_abs_value = Kokkos::ArithTraits<Scalar>::abs(rhs.value);
853 value = lhs_abs_value > rhs_abs_value ? lhs_abs_value : rhs_abs_value;
857 KOKKOS_FUNCTION AbsMaxHelper operator+(AbsMaxHelper
const& rhs)
const {
858 AbsMaxHelper ret = *
this;
864 template <
typename SC>
865 KOKKOS_INLINE_FUNCTION
866 void operator() (atomic_tag, SC& dst,
const SC& src)
const {
867 Kokkos::atomic_add(
reinterpret_cast<AbsMaxHelper<SC>*
>(&dst), AbsMaxHelper<SC>{src});
870 template <
typename SC>
871 KOKKOS_INLINE_FUNCTION
872 void operator() (nonatomic_tag, SC& dst,
const SC& src)
const {
873 auto dst_abs_value = Kokkos::ArithTraits<SC>::abs(dst);
874 auto src_abs_value = Kokkos::ArithTraits<SC>::abs(src);
875 dst = dst_abs_value > src_abs_value ? dst_abs_value : src_abs_value;
879 template <
typename ExecutionSpace,
884 typename Enabled =
void>
885 class UnpackArrayMultiColumn {
887 static_assert (Kokkos::is_view<DstView>::value,
888 "DstView must be a Kokkos::View.");
889 static_assert (Kokkos::is_view<SrcView>::value,
890 "SrcView must be a Kokkos::View.");
891 static_assert (Kokkos::is_view<IdxView>::value,
892 "IdxView must be a Kokkos::View.");
893 static_assert (static_cast<int> (DstView::rank) == 2,
894 "DstView must be a rank-2 Kokkos::View.");
895 static_assert (static_cast<int> (SrcView::rank) == 1,
896 "SrcView must be a rank-1 Kokkos::View.");
897 static_assert (static_cast<int> (IdxView::rank) == 1,
898 "IdxView must be a rank-1 Kokkos::View.");
901 typedef typename ExecutionSpace::execution_space execution_space;
902 typedef typename execution_space::size_type size_type;
912 UnpackArrayMultiColumn (
const ExecutionSpace& ,
917 const size_t numCols_) :
925 template<
class TagType>
926 KOKKOS_INLINE_FUNCTION
void
927 operator() (TagType tag,
const size_type k)
const
930 (std::is_same<TagType, atomic_tag>::value ||
931 std::is_same<TagType, nonatomic_tag>::value,
932 "TagType must be atomic_tag or nonatomic_tag.");
934 const typename IdxView::value_type localRow = idx(k);
935 const size_t offset = k*numCols;
936 for (
size_t j = 0; j < numCols; ++j) {
937 op (tag, dst(localRow, j), src(offset+j));
942 unpack (
const ExecutionSpace& execSpace,
947 const size_t numCols,
948 const bool use_atomic_updates)
950 if (use_atomic_updates) {
952 Kokkos::RangePolicy<atomic_tag, execution_space, size_type>;
954 (
"Tpetra::MultiVector unpack const stride atomic",
955 range_type (0, idx.size ()),
956 UnpackArrayMultiColumn (execSpace, dst, src, idx, op, numCols));
960 Kokkos::RangePolicy<nonatomic_tag, execution_space, size_type>;
962 (
"Tpetra::MultiVector unpack const stride nonatomic",
963 range_type (0, idx.size ()),
964 UnpackArrayMultiColumn (execSpace, dst, src, idx, op, numCols));
969 template <
typename ExecutionSpace,
974 typename SizeType =
typename ExecutionSpace::execution_space::size_type,
975 typename Enabled =
void>
976 class UnpackArrayMultiColumnWithBoundsCheck {
978 static_assert (Kokkos::is_view<DstView>::value,
979 "DstView must be a Kokkos::View.");
980 static_assert (Kokkos::is_view<SrcView>::value,
981 "SrcView must be a Kokkos::View.");
982 static_assert (Kokkos::is_view<IdxView>::value,
983 "IdxView must be a Kokkos::View.");
984 static_assert (static_cast<int> (DstView::rank) == 2,
985 "DstView must be a rank-2 Kokkos::View.");
986 static_assert (static_cast<int> (SrcView::rank) == 1,
987 "SrcView must be a rank-1 Kokkos::View.");
988 static_assert (static_cast<int> (IdxView::rank) == 1,
989 "IdxView must be a rank-1 Kokkos::View.");
990 static_assert (std::is_integral<SizeType>::value,
991 "SizeType must be a built-in integer type.");
994 using execution_space =
typename ExecutionSpace::execution_space;
995 using size_type = SizeType;
996 using value_type = size_t;
1006 UnpackArrayMultiColumnWithBoundsCheck (
const ExecutionSpace& ,
1007 const DstView& dst_,
1008 const SrcView& src_,
1009 const IdxView& idx_,
1011 const size_type numCols_) :
1019 template<
class TagType>
1020 KOKKOS_INLINE_FUNCTION
void
1021 operator() (TagType tag,
1023 size_t& lclErrCount)
const
1026 (std::is_same<TagType, atomic_tag>::value ||
1027 std::is_same<TagType, nonatomic_tag>::value,
1028 "TagType must be atomic_tag or nonatomic_tag.");
1029 using index_type =
typename IdxView::non_const_value_type;
1031 const index_type lclRow = idx(k);
1032 if (lclRow < static_cast<index_type> (0) ||
1033 lclRow >= static_cast<index_type> (dst.extent (0))) {
1037 const size_type offset = k*numCols;
1038 for (size_type j = 0; j < numCols; ++j) {
1039 op (tag, dst(lclRow,j), src(offset+j));
1044 template<
class TagType>
1045 KOKKOS_INLINE_FUNCTION
void
1046 init (TagType,
size_t& initialErrorCount)
const {
1047 initialErrorCount = 0;
1050 template<
class TagType>
1051 KOKKOS_INLINE_FUNCTION
void
1053 size_t& dstErrorCount,
1054 const size_t& srcErrorCount)
const
1056 dstErrorCount += srcErrorCount;
1060 unpack (
const ExecutionSpace& execSpace,
1065 const size_type numCols,
1066 const bool use_atomic_updates)
1068 using index_type =
typename IdxView::non_const_value_type;
1070 size_t errorCount = 0;
1071 if (use_atomic_updates) {
1073 Kokkos::RangePolicy<atomic_tag, execution_space, size_type>;
1074 Kokkos::parallel_reduce
1075 (
"Tpetra::MultiVector unpack multicol const stride atomic debug only",
1076 range_type (0, idx.size ()),
1077 UnpackArrayMultiColumnWithBoundsCheck (execSpace, dst, src,
1083 Kokkos::RangePolicy<nonatomic_tag, execution_space, size_type>;
1084 Kokkos::parallel_reduce
1085 (
"Tpetra::MultiVector unpack multicol const stride nonatomic debug only",
1086 range_type (0, idx.size ()),
1087 UnpackArrayMultiColumnWithBoundsCheck (execSpace, dst, src,
1092 if (errorCount != 0) {
1096 auto idx_h = Kokkos::create_mirror_view (idx);
1101 std::vector<index_type> badIndices;
1102 const size_type numInds = idx_h.extent (0);
1103 for (size_type k = 0; k < numInds; ++k) {
1104 if (idx_h(k) < static_cast<index_type> (0) ||
1105 idx_h(k) >= static_cast<index_type> (dst.extent (0))) {
1106 badIndices.push_back (idx_h(k));
1110 if (errorCount != badIndices.size ()) {
1111 std::ostringstream os;
1112 os <<
"MultiVector unpack kernel: errorCount = " << errorCount
1113 <<
" != badIndices.size() = " << badIndices.size ()
1114 <<
". This should never happen. "
1115 "Please report this to the Tpetra developers.";
1116 throw std::logic_error (os.str ());
1119 std::ostringstream os;
1120 os <<
"MultiVector unpack kernel had " << badIndices.size ()
1121 <<
" out-of bounds index/ices. Here they are: [";
1122 for (
size_t k = 0; k < badIndices.size (); ++k) {
1123 os << badIndices[k];
1124 if (k + 1 < badIndices.size ()) {
1129 throw std::runtime_error (os.str ());
1134 template <
typename ExecutionSpace,
1140 unpack_array_multi_column (
const ExecutionSpace& execSpace,
1145 const size_t numCols,
1146 const bool use_atomic_updates,
1149 static_assert (Kokkos::is_view<DstView>::value,
1150 "DstView must be a Kokkos::View.");
1151 static_assert (Kokkos::is_view<SrcView>::value,
1152 "SrcView must be a Kokkos::View.");
1153 static_assert (Kokkos::is_view<IdxView>::value,
1154 "IdxView must be a Kokkos::View.");
1155 static_assert (static_cast<int> (DstView::rank) == 2,
1156 "DstView must be a rank-2 Kokkos::View.");
1157 static_assert (static_cast<int> (SrcView::rank) == 1,
1158 "SrcView must be a rank-1 Kokkos::View.");
1159 static_assert (static_cast<int> (IdxView::rank) == 1,
1160 "IdxView must be a rank-1 Kokkos::View.");
1163 typedef UnpackArrayMultiColumnWithBoundsCheck<ExecutionSpace,
1164 DstView, SrcView, IdxView, Op> impl_type;
1165 impl_type::unpack (execSpace, dst, src, idx, op, numCols,
1166 use_atomic_updates);
1169 typedef UnpackArrayMultiColumn<ExecutionSpace,
1170 DstView, SrcView, IdxView, Op> impl_type;
1171 impl_type::unpack (execSpace, dst, src, idx, op, numCols,
1172 use_atomic_updates);
1176 template <
typename ExecutionSpace,
1182 typename Enabled =
void>
1183 class UnpackArrayMultiColumnVariableStride {
1185 static_assert (Kokkos::is_view<DstView>::value,
1186 "DstView must be a Kokkos::View.");
1187 static_assert (Kokkos::is_view<SrcView>::value,
1188 "SrcView must be a Kokkos::View.");
1189 static_assert (Kokkos::is_view<IdxView>::value,
1190 "IdxView must be a Kokkos::View.");
1191 static_assert (Kokkos::is_view<ColView>::value,
1192 "ColView must be a Kokkos::View.");
1193 static_assert (static_cast<int> (DstView::rank) == 2,
1194 "DstView must be a rank-2 Kokkos::View.");
1195 static_assert (static_cast<int> (SrcView::rank) == 1,
1196 "SrcView must be a rank-1 Kokkos::View.");
1197 static_assert (static_cast<int> (IdxView::rank) == 1,
1198 "IdxView must be a rank-1 Kokkos::View.");
1199 static_assert (static_cast<int> (ColView::rank) == 1,
1200 "ColView must be a rank-1 Kokkos::View.");
1203 using execution_space =
typename ExecutionSpace::execution_space;
1204 using size_type =
typename execution_space::size_type;
1215 UnpackArrayMultiColumnVariableStride (
const ExecutionSpace& ,
1216 const DstView& dst_,
1217 const SrcView& src_,
1218 const IdxView& idx_,
1219 const ColView& col_,
1221 const size_t numCols_) :
1230 template<
class TagType>
1231 KOKKOS_INLINE_FUNCTION
void
1232 operator() (TagType tag,
const size_type k)
const
1235 (std::is_same<TagType, atomic_tag>::value ||
1236 std::is_same<TagType, nonatomic_tag>::value,
1237 "TagType must be atomic_tag or nonatomic_tag.");
1239 const typename IdxView::value_type localRow = idx(k);
1240 const size_t offset = k*numCols;
1241 for (
size_t j = 0; j < numCols; ++j) {
1242 op (tag, dst(localRow, col(j)), src(offset+j));
1247 unpack (
const ExecutionSpace& execSpace,
1253 const size_t numCols,
1254 const bool use_atomic_updates)
1256 if (use_atomic_updates) {
1258 Kokkos::RangePolicy<atomic_tag, execution_space, size_type>;
1259 Kokkos::parallel_for
1260 (
"Tpetra::MultiVector unpack var stride atomic",
1261 range_type (0, idx.size ()),
1262 UnpackArrayMultiColumnVariableStride (execSpace, dst, src,
1263 idx, col, op, numCols));
1267 Kokkos::RangePolicy<nonatomic_tag, execution_space, size_type>;
1268 Kokkos::parallel_for
1269 (
"Tpetra::MultiVector unpack var stride nonatomic",
1270 range_type (0, idx.size ()),
1271 UnpackArrayMultiColumnVariableStride (execSpace, dst, src,
1272 idx, col, op, numCols));
1277 template <
typename ExecutionSpace,
1283 typename SizeType =
typename ExecutionSpace::execution_space::size_type,
1284 typename Enabled =
void>
1285 class UnpackArrayMultiColumnVariableStrideWithBoundsCheck {
1287 static_assert (Kokkos::is_view<DstView>::value,
1288 "DstView must be a Kokkos::View.");
1289 static_assert (Kokkos::is_view<SrcView>::value,
1290 "SrcView must be a Kokkos::View.");
1291 static_assert (Kokkos::is_view<IdxView>::value,
1292 "IdxView must be a Kokkos::View.");
1293 static_assert (Kokkos::is_view<ColView>::value,
1294 "ColView must be a Kokkos::View.");
1295 static_assert (static_cast<int> (DstView::rank) == 2,
1296 "DstView must be a rank-2 Kokkos::View.");
1297 static_assert (static_cast<int> (SrcView::rank) == 1,
1298 "SrcView must be a rank-1 Kokkos::View.");
1299 static_assert (static_cast<int> (IdxView::rank) == 1,
1300 "IdxView must be a rank-1 Kokkos::View.");
1301 static_assert (static_cast<int> (ColView::rank) == 1,
1302 "ColView must be a rank-1 Kokkos::View.");
1303 static_assert (std::is_integral<SizeType>::value,
1304 "SizeType must be a built-in integer type.");
1307 using execution_space =
typename ExecutionSpace::execution_space;
1308 using size_type = SizeType;
1309 using value_type = size_t;
1320 UnpackArrayMultiColumnVariableStrideWithBoundsCheck
1321 (
const ExecutionSpace& ,
1322 const DstView& dst_,
1323 const SrcView& src_,
1324 const IdxView& idx_,
1325 const ColView& col_,
1327 const size_t numCols_) :
1336 template<
class TagType>
1337 KOKKOS_INLINE_FUNCTION
void
1338 operator() (TagType tag,
1340 value_type& lclErrorCount)
const
1343 (std::is_same<TagType, atomic_tag>::value ||
1344 std::is_same<TagType, nonatomic_tag>::value,
1345 "TagType must be atomic_tag or nonatomic_tag.");
1346 using row_index_type =
typename IdxView::non_const_value_type;
1347 using col_index_type =
typename ColView::non_const_value_type;
1349 const row_index_type lclRow = idx(k);
1350 if (lclRow < static_cast<row_index_type> (0) ||
1351 lclRow >= static_cast<row_index_type> (dst.extent (0))) {
1355 const size_type offset = k * numCols;
1356 for (size_type j = 0; j < numCols; ++j) {
1357 const col_index_type lclCol = col(j);
1358 if (Impl::outOfBounds<col_index_type> (lclCol, dst.extent (1))) {
1362 op (tag, dst(lclRow, col(j)), src(offset+j));
1368 KOKKOS_INLINE_FUNCTION
void
1369 init (value_type& initialErrorCount)
const {
1370 initialErrorCount = 0;
1373 KOKKOS_INLINE_FUNCTION
void
1374 join (value_type& dstErrorCount,
1375 const value_type& srcErrorCount)
const
1377 dstErrorCount += srcErrorCount;
1381 unpack (
const ExecutionSpace& execSpace,
1387 const size_type numCols,
1388 const bool use_atomic_updates)
1390 using row_index_type =
typename IdxView::non_const_value_type;
1391 using col_index_type =
typename ColView::non_const_value_type;
1393 size_t errorCount = 0;
1394 if (use_atomic_updates) {
1396 Kokkos::RangePolicy<atomic_tag, execution_space, size_type>;
1397 Kokkos::parallel_reduce
1398 (
"Tpetra::MultiVector unpack var stride atomic debug only",
1399 range_type (0, idx.size ()),
1400 UnpackArrayMultiColumnVariableStrideWithBoundsCheck
1401 (execSpace, dst, src, idx, col, op, numCols),
1406 Kokkos::RangePolicy<nonatomic_tag, execution_space, size_type>;
1407 Kokkos::parallel_reduce
1408 (
"Tpetra::MultiVector unpack var stride nonatomic debug only",
1409 range_type (0, idx.size ()),
1410 UnpackArrayMultiColumnVariableStrideWithBoundsCheck
1411 (execSpace, dst, src, idx, col, op, numCols),
1415 if (errorCount != 0) {
1416 constexpr
size_t maxNumBadIndicesToPrint = 100;
1418 std::ostringstream os;
1419 os <<
"Tpetra::MultiVector multicolumn variable stride unpack kernel "
1420 "found " << errorCount
1421 <<
" error" << (errorCount != size_t (1) ?
"s" :
"") <<
". ";
1427 auto idx_h = Kokkos::create_mirror_view (idx);
1432 std::vector<row_index_type> badRows;
1433 const size_type numRowInds = idx_h.extent (0);
1434 for (size_type k = 0; k < numRowInds; ++k) {
1435 if (idx_h(k) < static_cast<row_index_type> (0) ||
1436 idx_h(k) >= static_cast<row_index_type> (dst.extent (0))) {
1437 badRows.push_back (idx_h(k));
1441 if (badRows.size () != 0) {
1442 os << badRows.size () <<
" out-of-bounds row ind"
1443 << (badRows.size () != size_t (1) ?
"ices" :
"ex");
1444 if (badRows.size () <= maxNumBadIndicesToPrint) {
1446 for (
size_t k = 0; k < badRows.size (); ++k) {
1448 if (k + 1 < badRows.size ()) {
1459 os <<
"No out-of-bounds row indices. ";
1464 auto col_h = Kokkos::create_mirror_view (col);
1469 std::vector<col_index_type> badCols;
1470 const size_type numColInds = col_h.extent (0);
1471 for (size_type k = 0; k < numColInds; ++k) {
1472 if (Impl::outOfBounds<col_index_type> (col_h(k), dst.extent (1))) {
1473 badCols.push_back (col_h(k));
1477 if (badCols.size () != 0) {
1478 os << badCols.size () <<
" out-of-bounds column ind"
1479 << (badCols.size () != size_t (1) ?
"ices" :
"ex");
1480 if (badCols.size () <= maxNumBadIndicesToPrint) {
1481 for (
size_t k = 0; k < badCols.size (); ++k) {
1484 if (k + 1 < badCols.size ()) {
1495 os <<
"No out-of-bounds column indices. ";
1498 TEUCHOS_TEST_FOR_EXCEPTION
1499 (errorCount != 0 && badRows.size () == 0 && badCols.size () == 0,
1500 std::logic_error,
"Tpetra::MultiVector variable stride unpack "
1501 "kernel reports errorCount=" << errorCount <<
", but we failed "
1502 "to find any bad rows or columns. This should never happen. "
1503 "Please report this to the Tpetra developers.");
1505 throw std::runtime_error (os.str ());
1510 template <
typename ExecutionSpace,
1517 unpack_array_multi_column_variable_stride (
const ExecutionSpace& execSpace,
1523 const size_t numCols,
1524 const bool use_atomic_updates,
1527 static_assert (Kokkos::is_view<DstView>::value,
1528 "DstView must be a Kokkos::View.");
1529 static_assert (Kokkos::is_view<SrcView>::value,
1530 "SrcView must be a Kokkos::View.");
1531 static_assert (Kokkos::is_view<IdxView>::value,
1532 "IdxView must be a Kokkos::View.");
1533 static_assert (Kokkos::is_view<ColView>::value,
1534 "ColView must be a Kokkos::View.");
1535 static_assert (static_cast<int> (DstView::rank) == 2,
1536 "DstView must be a rank-2 Kokkos::View.");
1537 static_assert (static_cast<int> (SrcView::rank) == 1,
1538 "SrcView must be a rank-1 Kokkos::View.");
1539 static_assert (static_cast<int> (IdxView::rank) == 1,
1540 "IdxView must be a rank-1 Kokkos::View.");
1541 static_assert (static_cast<int> (ColView::rank) == 1,
1542 "ColView must be a rank-1 Kokkos::View.");
1546 UnpackArrayMultiColumnVariableStrideWithBoundsCheck<ExecutionSpace,
1547 DstView, SrcView, IdxView, ColView, Op>;
1548 impl_type::unpack (execSpace, dst, src, idx, col, op, numCols,
1549 use_atomic_updates);
1552 using impl_type = UnpackArrayMultiColumnVariableStride<ExecutionSpace,
1553 DstView, SrcView, IdxView, ColView, Op>;
1554 impl_type::unpack (execSpace, dst, src, idx, col, op, numCols,
1555 use_atomic_updates);
1559 template <
typename DstView,
typename SrcView,
1560 typename DstIdxView,
typename SrcIdxView,
typename Op,
1561 typename Enabled =
void>
1562 struct PermuteArrayMultiColumn {
1563 using size_type =
typename DstView::size_type;
1572 PermuteArrayMultiColumn (
const DstView& dst_,
1573 const SrcView& src_,
1574 const DstIdxView& dst_idx_,
1575 const SrcIdxView& src_idx_,
1576 const size_t numCols_,
1578 dst(dst_), src(src_), dst_idx(dst_idx_), src_idx(src_idx_),
1579 numCols(numCols_), op(op_) {}
1581 KOKKOS_INLINE_FUNCTION
void
1582 operator() (
const size_type k)
const {
1583 const typename DstIdxView::value_type toRow = dst_idx(k);
1584 const typename SrcIdxView::value_type fromRow = src_idx(k);
1586 for (
size_t j = 0; j < numCols; ++j) {
1587 op(tag, dst(toRow, j), src(fromRow, j));
1591 template <
typename ExecutionSpace>
1593 permute (
const ExecutionSpace &space,
1596 const DstIdxView& dst_idx,
1597 const SrcIdxView& src_idx,
1598 const size_t numCols,
1602 Kokkos::RangePolicy<ExecutionSpace, size_type>;
1604 const size_type n = std::min (dst_idx.size (), src_idx.size ());
1605 Kokkos::parallel_for
1606 (
"Tpetra::MultiVector permute multicol const stride",
1607 range_type (space, 0, n),
1608 PermuteArrayMultiColumn (dst, src, dst_idx, src_idx, numCols, op));
1615 template <
typename ExecutionSpace,
typename DstView,
typename SrcView,
1616 typename DstIdxView,
typename SrcIdxView,
typename Op>
1617 void permute_array_multi_column(
const ExecutionSpace &space,
const DstView &dst,
1618 const SrcView &src,
const DstIdxView &dst_idx,
1619 const SrcIdxView &src_idx,
size_t numCols,
1621 PermuteArrayMultiColumn<DstView, SrcView, DstIdxView, SrcIdxView,
1622 Op>::permute(space, dst, src, dst_idx, src_idx,
1629 template <
typename DstView,
typename SrcView,
1630 typename DstIdxView,
typename SrcIdxView,
typename Op>
1631 void permute_array_multi_column(
const DstView& dst,
1633 const DstIdxView& dst_idx,
1634 const SrcIdxView& src_idx,
1637 using execution_space =
typename DstView::execution_space;
1638 PermuteArrayMultiColumn<DstView,SrcView,DstIdxView,SrcIdxView,Op>::permute(
1639 execution_space(), dst, src, dst_idx, src_idx, numCols, op);
1642 template <
typename DstView,
typename SrcView,
1643 typename DstIdxView,
typename SrcIdxView,
1644 typename DstColView,
typename SrcColView,
typename Op,
1645 typename Enabled =
void>
1646 struct PermuteArrayMultiColumnVariableStride {
1647 using size_type =
typename DstView::size_type;
1658 PermuteArrayMultiColumnVariableStride(
const DstView& dst_,
1659 const SrcView& src_,
1660 const DstIdxView& dst_idx_,
1661 const SrcIdxView& src_idx_,
1662 const DstColView& dst_col_,
1663 const SrcColView& src_col_,
1664 const size_t numCols_,
1666 dst(dst_), src(src_), dst_idx(dst_idx_), src_idx(src_idx_),
1667 dst_col(dst_col_), src_col(src_col_),
1668 numCols(numCols_), op(op_) {}
1670 KOKKOS_INLINE_FUNCTION
void
1671 operator() (
const size_type k)
const {
1672 const typename DstIdxView::value_type toRow = dst_idx(k);
1673 const typename SrcIdxView::value_type fromRow = src_idx(k);
1674 const nonatomic_tag tag;
1675 for (
size_t j = 0; j < numCols; ++j) {
1676 op(tag, dst(toRow, dst_col(j)), src(fromRow, src_col(j)));
1680 template <
typename ExecutionSpace>
1682 permute (
const ExecutionSpace &space,
1685 const DstIdxView& dst_idx,
1686 const SrcIdxView& src_idx,
1687 const DstColView& dst_col,
1688 const SrcColView& src_col,
1689 const size_t numCols,
1693 static_assert(Kokkos::SpaceAccessibility<
1694 ExecutionSpace,
typename DstView::memory_space>::accessible,
1695 "ExecutionSpace must be able to access DstView");
1697 using range_type = Kokkos::RangePolicy<ExecutionSpace, size_type>;
1698 const size_type n = std::min (dst_idx.size (), src_idx.size ());
1699 Kokkos::parallel_for
1700 (
"Tpetra::MultiVector permute multicol var stride",
1701 range_type (space, 0, n),
1702 PermuteArrayMultiColumnVariableStride (dst, src, dst_idx, src_idx,
1703 dst_col, src_col, numCols, op));
1710 template <
typename ExecutionSpace,
typename DstView,
typename SrcView,
1711 typename DstIdxView,
typename SrcIdxView,
typename DstColView,
1712 typename SrcColView,
typename Op>
1713 void permute_array_multi_column_variable_stride(
1714 const ExecutionSpace &space,
const DstView &dst,
const SrcView &src,
1715 const DstIdxView &dst_idx,
const SrcIdxView &src_idx,
1716 const DstColView &dst_col,
const SrcColView &src_col,
size_t numCols,
1718 PermuteArrayMultiColumnVariableStride<DstView, SrcView, DstIdxView,
1719 SrcIdxView, DstColView, SrcColView,
1720 Op>::permute(space, dst, src, dst_idx,
1721 src_idx, dst_col, src_col,
1728 template <
typename DstView,
typename SrcView,
1729 typename DstIdxView,
typename SrcIdxView,
1730 typename DstColView,
typename SrcColView,
typename Op>
1731 void permute_array_multi_column_variable_stride(
const DstView& dst,
1733 const DstIdxView& dst_idx,
1734 const SrcIdxView& src_idx,
1735 const DstColView& dst_col,
1736 const SrcColView& src_col,
1739 using execution_space =
typename DstView::execution_space;
1740 PermuteArrayMultiColumnVariableStride<DstView,SrcView,
1741 DstIdxView,SrcIdxView,DstColView,SrcColView,Op>::permute(
1742 execution_space(), dst, src, dst_idx, src_idx, dst_col, src_col, numCols, op);
1749 #endif // TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_DIST_OBJECT_KERNELS_HPP
KOKKOS_INLINE_FUNCTION bool outOfBounds(const IntegerType x, const IntegerType exclusiveUpperBound)
Is x out of bounds? That is, is x less than zero, or greater than or equal to the given exclusive upp...
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst.
Is x out of bounds? That is, is x less than zero, or greater than or equal to the given exclusive upp...