44 #ifndef TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_LOCAL_DEEP_COPY_HPP
45 #define TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_LOCAL_DEEP_COPY_HPP
47 #include "Kokkos_Core.hpp"
55 template<
class DstViewType,
class SrcViewType>
57 copyConvertResolvingPossibleAliasing (
const DstViewType& dst,
58 const SrcViewType& src)
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 ());
72 if (src_beg == dst_beg && src_end == dst_end) {
75 else if (dst_end <= src_beg || src_end <= dst_beg) {
86 auto src_copy = Kokkos::create_mirror (Kokkos::HostSpace (), src);
129 template<
class DstViewType,
131 class DstWhichVecsType,
132 class SrcWhichVecsType>
135 const SrcViewType& src,
136 const bool dstConstStride,
137 const bool srcConstStride,
138 const DstWhichVecsType& dstWhichVecs,
139 const SrcWhichVecsType& srcWhichVecs)
142 using Kokkos::subview;
143 using size_type =
typename DstViewType::size_type;
145 if (dstConstStride && srcConstStride) {
146 copyConvertResolvingPossibleAliasing (dst, src);
149 const size_type numCols = dstConstStride ?
150 static_cast<size_type
> (srcWhichVecs.size ()) :
151 static_cast<size_type> (dstWhichVecs.size ());
152 for (size_type j = 0; j < numCols; ++j) {
153 const size_type dst_col = dstConstStride ? j :
154 static_cast<size_type
> (dstWhichVecs[j]);
155 const auto dst_j = subview (dst, ALL (), dst_col);
156 const size_type src_col = srcConstStride ? j :
157 static_cast<size_type
> (srcWhichVecs[j]);
158 const auto src_j = subview (src, ALL (), src_col);
160 copyConvertResolvingPossibleAliasing (dst_j, src_j);
168 template<
class DstViewType,
172 const SrcViewType& src)
174 return copyConvertResolvingPossibleAliasing (dst, src);
180 #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.