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 "Kokkos_ArithTraits.hpp"
71 template<
class DstViewType,
73 class DstWhichVecsType,
74 class SrcWhichVecsType,
75 const bool DstConstStride,
76 const bool SrcConstStride>
77 struct LocalDeepCopyFunctor {
78 typedef typename DstViewType::execution_space execution_space;
79 typedef typename DstViewType::size_type index_type;
83 DstWhichVecsType dstWhichVecs_;
84 SrcWhichVecsType srcWhichVecs_;
85 const index_type numVecs_;
90 LocalDeepCopyFunctor (
const DstViewType& dst,
91 const SrcViewType& src,
92 const DstWhichVecsType& dstWhichVecs,
93 const SrcWhichVecsType& srcWhichVecs) :
96 dstWhichVecs_ (dstWhichVecs),
97 srcWhichVecs_ (srcWhichVecs),
98 numVecs_ (DstConstStride ? dst.extent (1) : dstWhichVecs.extent (0))
100 TEUCHOS_TEST_FOR_EXCEPTION(
101 ! DstConstStride && ! SrcConstStride &&
102 dstWhichVecs.extent (0) != srcWhichVecs.extent (0),
103 std::invalid_argument,
"LocalDeepCopyFunctor (4-arg constructor): "
104 "Neither src nor dst have constant stride, but "
105 "dstWhichVecs.extent(0) = " << dstWhichVecs.extent (0)
106 <<
" != srcWhichVecs.extent(0) = " << srcWhichVecs.extent (0)
108 TEUCHOS_TEST_FOR_EXCEPTION(
109 DstConstStride && ! SrcConstStride &&
110 srcWhichVecs.extent (0) != dst.extent (1),
111 std::invalid_argument,
"LocalDeepCopyFunctor (4-arg constructor): "
112 "src does not have constant stride, but srcWhichVecs.extent(0) = "
113 << srcWhichVecs.extent (0) <<
" != dst.extent(1) = "
114 << dst.extent (1) <<
".");
115 TEUCHOS_TEST_FOR_EXCEPTION(
116 ! DstConstStride && SrcConstStride &&
117 dstWhichVecs.extent (0) != src.extent (1),
118 std::invalid_argument,
"LocalDeepCopyFunctor (4-arg constructor): "
119 "dst does not have constant stride, but dstWhichVecs.extent(0) = "
120 << dstWhichVecs.extent (0) <<
" != src.extent(1) = "
121 << src.extent (1) <<
".");
126 LocalDeepCopyFunctor (
const DstViewType& dst,
const SrcViewType& src) :
129 numVecs_ (dst.extent (1))
131 TEUCHOS_TEST_FOR_EXCEPTION(
132 ! DstConstStride || ! SrcConstStride, std::logic_error,
133 "Tpetra::LocalDeepCopyFunctor: You may not use the constant-stride "
134 "constructor if either of the Boolean template parameters is false.");
137 void KOKKOS_INLINE_FUNCTION operator () (
const index_type i)
const {
138 if (DstConstStride) {
139 if (SrcConstStride) {
140 for (index_type j = 0; j < numVecs_; ++j) {
141 dst_(i,j) = src_(i,j);
144 for (index_type j = 0; j < numVecs_; ++j) {
145 dst_(i,j) = src_(i,srcWhichVecs_(j));
149 if (SrcConstStride) {
150 for (index_type j = 0; j < numVecs_; ++j) {
151 dst_(i,dstWhichVecs_(j)) = src_(i,j);
154 for (index_type j = 0; j < numVecs_; ++j) {
155 dst_(i,dstWhichVecs_(j)) = src_(i,srcWhichVecs_(j));
162 template<
class DstViewType,
163 class SrcViewType = DstViewType,
164 class DstWhichVecsType = Kokkos::View<const typename DstViewType::size_type*, typename DstViewType::execution_space>,
165 class SrcWhichVecsType = Kokkos::View<const typename SrcViewType::size_type*, typename SrcViewType::execution_space> >
166 struct LocalDeepCopy {
167 typedef typename DstViewType::size_type index_type;
170 run (
const DstViewType& dst,
const SrcViewType& src,
171 const bool dstConstStride,
const bool srcConstStride)
173 TEUCHOS_TEST_FOR_EXCEPTION(
174 ! dstConstStride || ! srcConstStride, std::invalid_argument,
175 "LocalDeepCopy::run: You may only call the 4-argument version of this "
176 "function if dstConstStride and srcConstStride are both true.");
189 typedef LocalDeepCopyFunctor<DstViewType, SrcViewType,
190 DstWhichVecsType, SrcWhichVecsType,
true,
true> functor_type;
191 functor_type f (dst, src);
192 typedef typename DstViewType::execution_space execution_space;
193 typedef decltype (dst.extent (0)) size_type;
194 typedef Kokkos::RangePolicy<execution_space, size_type> range_type;
196 Kokkos::parallel_for (
"Tpetra::Details::LocalDeepCopy(x,y,stride)",range_type (0, dst.extent (0)), f);
200 run (
const DstViewType& dst,
const SrcViewType& src,
201 const bool dstConstStride,
const bool srcConstStride,
202 const DstWhichVecsType& dstWhichVecs,
203 const SrcWhichVecsType& srcWhichVecs)
214 typedef typename DstViewType::execution_space execution_space;
215 typedef decltype (dst.extent (0)) size_type;
216 typedef Kokkos::RangePolicy<execution_space, size_type> range_type;
218 if (dstConstStride) {
219 if (srcConstStride) {
222 typedef LocalDeepCopyFunctor<DstViewType, SrcViewType,
223 DstWhichVecsType, SrcWhichVecsType,
true,
true> functor_type;
224 functor_type f (dst, src);
225 Kokkos::parallel_for (
"Tpetra::Details::LocalDeepCopy(x,y,stride,whichVecs,0)",range_type (0, dst.extent (0)), f);
228 typedef LocalDeepCopyFunctor<DstViewType, SrcViewType,
229 DstWhichVecsType, SrcWhichVecsType,
true,
false> functor_type;
230 functor_type f (dst, src, srcWhichVecs, srcWhichVecs);
231 Kokkos::parallel_for (
"Tpetra::Details::LocalDeepCopy(x,y,stride,whichVecs,1)",range_type (0, dst.extent (0)), f);
235 if (srcConstStride) {
236 typedef LocalDeepCopyFunctor<DstViewType, SrcViewType,
237 DstWhichVecsType, SrcWhichVecsType,
false,
true> functor_type;
238 functor_type f (dst, src, dstWhichVecs, dstWhichVecs);
239 Kokkos::parallel_for (
"Tpetra::Details::LocalDeepCopy(x,y,stride,whichVecs,2)",range_type (0, dst.extent (0)), f);
242 typedef LocalDeepCopyFunctor<DstViewType, SrcViewType,
243 DstWhichVecsType, SrcWhichVecsType,
false,
false> functor_type;
244 functor_type f (dst, src, dstWhichVecs, srcWhichVecs);
245 Kokkos::parallel_for (
"Tpetra::Details::LocalDeepCopy(x,y,stride,whichVecs,3)",range_type (0, dst.extent (0)), f);
252 template<
class DstViewType,
254 class DstWhichVecsType,
255 class SrcWhichVecsType>
257 localDeepCopy (
const DstViewType& dst,
const SrcViewType& src,
258 const bool dstConstStride,
const bool srcConstStride,
259 const DstWhichVecsType& dstWhichVecs,
260 const SrcWhichVecsType& srcWhichVecs)
262 typedef LocalDeepCopy<DstViewType, SrcViewType,
263 DstWhichVecsType, SrcWhichVecsType> impl_type;
264 impl_type::run (dst, src, dstConstStride, srcConstStride,
265 dstWhichVecs, srcWhichVecs);
268 template<
class DstViewType,
class SrcViewType>
270 localDeepCopyConstStride (
const DstViewType& dst,
const SrcViewType& src)
272 typedef LocalDeepCopy<DstViewType, SrcViewType> impl_type;
273 impl_type::run (dst, src,
true,
true);
279 #endif // TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_LOCAL_DEEP_COPY_HPP