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
23 copyConvertResolvingPossibleAliasing (const DstViewType& dst,
24  const SrcViewType& src)
25 {
26  // NOTE: It's important to do the addition _inside_ the
27  // reinterpret-cast. If you reinterpret_cast the separate results,
28  // you may get the wrong answer (e.g., because ptrdiff_t is signed,
29  // and pointers may have arbitrary 64-bit virtual addresses). I'm
30  // speaking from experience here.
31  const ptrdiff_t dst_beg = reinterpret_cast<ptrdiff_t> (dst.data ());
32  const ptrdiff_t dst_end =
33  reinterpret_cast<ptrdiff_t> (dst.data () + dst.span ());
34  const ptrdiff_t src_beg = reinterpret_cast<ptrdiff_t> (src.data ());
35  const ptrdiff_t src_end =
36  reinterpret_cast<ptrdiff_t> (src.data () + src.span ());
37 
38  if (src_beg == dst_beg && src_end == dst_end) {
39  // Do nothing; there's no need to copy
40  }
41  else if (dst_end <= src_beg || src_end <= dst_beg) { // no aliasing
43  }
44  else {
45  // dst and src alias each other, so we can't call
46  // Kokkos::deep_copy(dst,src) directly (Kokkos detects this and
47  // throws, at least in debug mode). Instead, we make temporary
48  // host storage (create_mirror always makes a new allocation,
49  // unlike create_mirror_view). Use host because it's cheaper to
50  // allocate. Hopefully users aren't doing aliased copies in a
51  // tight loop.
52  auto src_copy = Kokkos::create_mirror (Kokkos::HostSpace (), src);
53 
54  // DEEP_COPY REVIEW - NOT TESTED
55  Kokkos::deep_copy (src_copy, src);
56  ::Tpetra::Details::copyConvert (dst, src_copy);
57  }
58 }
59 
60 } // namespace (anonymous)
61 
95 template<class DstViewType,
96  class SrcViewType,
97  class DstWhichVecsType,
98  class SrcWhichVecsType>
99 void
100 localDeepCopy (const DstViewType& dst,
101  const SrcViewType& src,
102  const bool dstConstStride,
103  const bool srcConstStride,
104  const DstWhichVecsType& dstWhichVecs,
105  const SrcWhichVecsType& srcWhichVecs)
106 {
107  using Kokkos::ALL;
108  using Kokkos::subview;
109  using size_type = typename DstViewType::size_type;
110 
111  if (dstConstStride && srcConstStride) {
112  copyConvertResolvingPossibleAliasing (dst, src);
113  }
114  else {
115  const size_type numCols = dstConstStride ?
116  static_cast<size_type> (srcWhichVecs.size ()) :
117  static_cast<size_type> (dstWhichVecs.size ());
118  for (size_type j = 0; j < numCols; ++j) {
119  const size_type dst_col = dstConstStride ? j :
120  static_cast<size_type> (dstWhichVecs[j]);
121  const auto dst_j = subview (dst, ALL (), dst_col);
122  const size_type src_col = srcConstStride ? j :
123  static_cast<size_type> (srcWhichVecs[j]);
124  const auto src_j = subview (src, ALL (), src_col);
125 
126  copyConvertResolvingPossibleAliasing (dst_j, src_j);
127  }
128  }
129 }
130 
134 template<class DstViewType,
135  class SrcViewType>
136 void
137 localDeepCopyConstStride (const DstViewType& dst,
138  const SrcViewType& src)
139 {
140  return copyConvertResolvingPossibleAliasing (dst, src);
141 }
142 
143 } // Details namespace
144 } // Tpetra namespace
145 
146 #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.