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 /*
2 // @HEADER
3 // ***********************************************************************
4 //
5 // Tpetra: Templated Linear Algebra Services Package
6 // Copyright (2008) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ************************************************************************
41 // @HEADER
42 */
43 
44 #ifndef TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_LOCAL_DEEP_COPY_HPP
45 #define TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_LOCAL_DEEP_COPY_HPP
46 
47 #include "Kokkos_Core.hpp"
49 
50 namespace Tpetra {
51 namespace Details {
52 
53 namespace { // (anonymous)
54 
55 template<class DstViewType, class SrcViewType>
56 void
57 copyConvertResolvingPossibleAliasing (const DstViewType& dst,
58  const SrcViewType& src)
59 {
60  // NOTE: It's important to do the addition _inside_ the
61  // reinterpret-cast. If you reinterpret_cast the separate results,
62  // you may get the wrong answer (e.g., because ptrdiff_t is signed,
63  // and pointers may have arbitrary 64-bit virtual addresses). I'm
64  // speaking from experience here.
65  const ptrdiff_t dst_beg = reinterpret_cast<ptrdiff_t> (dst.data ());
66  const ptrdiff_t dst_end =
67  reinterpret_cast<ptrdiff_t> (dst.data () + dst.span ());
68  const ptrdiff_t src_beg = reinterpret_cast<ptrdiff_t> (src.data ());
69  const ptrdiff_t src_end =
70  reinterpret_cast<ptrdiff_t> (src.data () + src.span ());
71 
72  if (src_beg == dst_beg && src_end == dst_end) {
73  // Do nothing; there's no need to copy
74  }
75  else if (dst_end <= src_beg || src_end <= dst_beg) { // no aliasing
77  }
78  else {
79  // dst and src alias each other, so we can't call
80  // Kokkos::deep_copy(dst,src) directly (Kokkos detects this and
81  // throws, at least in debug mode). Instead, we make temporary
82  // host storage (create_mirror always makes a new allocation,
83  // unlike create_mirror_view). Use host because it's cheaper to
84  // allocate. Hopefully users aren't doing aliased copies in a
85  // tight loop.
86  auto src_copy = Kokkos::create_mirror (Kokkos::HostSpace (), src);
87  Kokkos::deep_copy (src_copy, src);
88  ::Tpetra::Details::copyConvert (dst, src_copy);
89  }
90 }
91 
92 } // namespace (anonymous)
93 
127 template<class DstViewType,
128  class SrcViewType,
129  class DstWhichVecsType,
130  class SrcWhichVecsType>
131 void
132 localDeepCopy (const DstViewType& dst,
133  const SrcViewType& src,
134  const bool dstConstStride,
135  const bool srcConstStride,
136  const DstWhichVecsType& dstWhichVecs,
137  const SrcWhichVecsType& srcWhichVecs)
138 {
139  using Kokkos::ALL;
140  using Kokkos::subview;
141  using size_type = typename DstViewType::size_type;
142 
143  if (dstConstStride && srcConstStride) {
144  copyConvertResolvingPossibleAliasing (dst, src);
145  }
146  else {
147  const size_type numCols = dstConstStride ?
148  static_cast<size_type> (srcWhichVecs.size ()) :
149  static_cast<size_type> (dstWhichVecs.size ());
150  for (size_type j = 0; j < numCols; ++j) {
151  const size_type dst_col = dstConstStride ? j :
152  static_cast<size_type> (dstWhichVecs[j]);
153  const auto dst_j = subview (dst, ALL (), dst_col);
154  const size_type src_col = srcConstStride ? j :
155  static_cast<size_type> (srcWhichVecs[j]);
156  const auto src_j = subview (src, ALL (), src_col);
157 
158  copyConvertResolvingPossibleAliasing (dst_j, src_j);
159  }
160  }
161 }
162 
166 template<class DstViewType,
167  class SrcViewType>
168 void
169 localDeepCopyConstStride (const DstViewType& dst,
170  const SrcViewType& src)
171 {
172  return copyConvertResolvingPossibleAliasing (dst, src);
173 }
174 
175 } // Details namespace
176 } // Tpetra namespace
177 
178 #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.