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"
48 #include "Tpetra_Details_copyConvert.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) {
76 ::Tpetra::Details::copyConvert (dst, src);
86 auto src_copy = Kokkos::create_mirror (Kokkos::HostSpace (), src);
88 ::Tpetra::Details::copyConvert (dst, src_copy);
127 template<
class DstViewType,
129 class DstWhichVecsType,
130 class SrcWhichVecsType>
133 const SrcViewType& src,
134 const bool dstConstStride,
135 const bool srcConstStride,
136 const DstWhichVecsType& dstWhichVecs,
137 const SrcWhichVecsType& srcWhichVecs)
140 using Kokkos::subview;
141 using size_type =
typename DstViewType::size_type;
143 if (dstConstStride && srcConstStride) {
144 copyConvertResolvingPossibleAliasing (dst, src);
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);
158 copyConvertResolvingPossibleAliasing (dst_j, src_j);
166 template<
class DstViewType,
170 const SrcViewType& src)
172 return copyConvertResolvingPossibleAliasing (dst, src);
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...
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 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.