Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_KokkosRefactor_Details_MultiVectorLocalDeepCopy.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_LOCAL_DEEP_COPY_HPP
11 #define TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_LOCAL_DEEP_COPY_HPP
12 
13 #include "Kokkos_Core.hpp"
15 
16 namespace Tpetra {
17 namespace Details {
18 
19 namespace { // (anonymous)
20 
21 template <class DstViewType, class SrcViewType>
22 void copyConvertResolvingPossibleAliasing(const DstViewType& dst,
23  const SrcViewType& src) {
24  // NOTE: It's important to do the addition _inside_ the
25  // reinterpret-cast. If you reinterpret_cast the separate results,
26  // you may get the wrong answer (e.g., because ptrdiff_t is signed,
27  // and pointers may have arbitrary 64-bit virtual addresses). I'm
28  // speaking from experience here.
29  const ptrdiff_t dst_beg = reinterpret_cast<ptrdiff_t>(dst.data());
30  const ptrdiff_t dst_end =
31  reinterpret_cast<ptrdiff_t>(dst.data() + dst.span());
32  const ptrdiff_t src_beg = reinterpret_cast<ptrdiff_t>(src.data());
33  const ptrdiff_t src_end =
34  reinterpret_cast<ptrdiff_t>(src.data() + src.span());
35 
36  if (src_beg == dst_beg && src_end == dst_end) {
37  // Do nothing; there's no need to copy
38  } else if (dst_end <= src_beg || src_end <= dst_beg) { // no aliasing
40  } else {
41  // dst and src alias each other, so we can't call
42  // Kokkos::deep_copy(dst,src) directly (Kokkos detects this and
43  // throws, at least in debug mode). Instead, we make temporary
44  // host storage (create_mirror always makes a new allocation,
45  // unlike create_mirror_view). Use host because it's cheaper to
46  // allocate. Hopefully users aren't doing aliased copies in a
47  // tight loop.
48  auto src_copy = Kokkos::create_mirror(Kokkos::HostSpace(), src);
49 
50  // DEEP_COPY REVIEW - NOT TESTED
51  Kokkos::deep_copy(src_copy, src);
52  ::Tpetra::Details::copyConvert(dst, src_copy);
53  }
54 }
55 
56 } // namespace
57 
91 template <class DstViewType,
92  class SrcViewType,
93  class DstWhichVecsType,
94  class SrcWhichVecsType>
95 void localDeepCopy(const DstViewType& dst,
96  const SrcViewType& src,
97  const bool dstConstStride,
98  const bool srcConstStride,
99  const DstWhichVecsType& dstWhichVecs,
100  const SrcWhichVecsType& srcWhichVecs) {
101  using Kokkos::ALL;
102  using Kokkos::subview;
103  using size_type = typename DstViewType::size_type;
104 
105  if (dstConstStride && srcConstStride) {
106  copyConvertResolvingPossibleAliasing(dst, src);
107  } else {
108  const size_type numCols = dstConstStride ? static_cast<size_type>(srcWhichVecs.size()) : static_cast<size_type>(dstWhichVecs.size());
109  for (size_type j = 0; j < numCols; ++j) {
110  const size_type dst_col = dstConstStride ? j : static_cast<size_type>(dstWhichVecs[j]);
111  const auto dst_j = subview(dst, ALL(), dst_col);
112  const size_type src_col = srcConstStride ? j : static_cast<size_type>(srcWhichVecs[j]);
113  const auto src_j = subview(src, ALL(), src_col);
114 
115  copyConvertResolvingPossibleAliasing(dst_j, src_j);
116  }
117  }
118 }
119 
123 template <class DstViewType,
124  class SrcViewType>
125 void localDeepCopyConstStride(const DstViewType& dst,
126  const SrcViewType& src) {
127  return copyConvertResolvingPossibleAliasing(dst, src);
128 }
129 
130 } // namespace Details
131 } // namespace Tpetra
132 
133 #endif // TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_LOCAL_DEEP_COPY_HPP
void localDeepCopyConstStride(const DstViewType &dst, const SrcViewType &src)
Implementation of Tpetra::MultiVector deep copy of local data, for when both the source and destinati...
Declare and define Tpetra::Details::copyConvert, an implementation detail of Tpetra (in particular...
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.
void copyConvert(const OutputViewType &dst, const InputViewType &src)
Copy values from the 1-D Kokkos::View src, to the 1-D Kokkos::View dst, of the same length...
void localDeepCopy(const DstViewType &dst, const SrcViewType &src, const bool dstConstStride, const bool srcConstStride, const DstWhichVecsType &dstWhichVecs, const SrcWhichVecsType &srcWhichVecs)
Implementation of Tpetra::MultiVector deep copy of local data.