10#ifndef TPETRA_DETAILS_COPYCONVERT_HPP
11#define TPETRA_DETAILS_COPYCONVERT_HPP
18#include "TpetraCore_config.h"
19#include "Kokkos_Core.hpp"
20#include "KokkosKernels_ArithTraits.hpp"
37template <
class OutputValueType,
39 const bool outputIsComplex =
40 KokkosKernels::ArithTraits<OutputValueType>::is_complex,
41 const bool inputIsComplex =
42 KokkosKernels::ArithTraits<InputValueType>::is_complex>
44 static KOKKOS_INLINE_FUNCTION
void
45 convert(OutputValueType& dst,
const InputValueType& src) {
50 dst = OutputValueType(src);
54template <
class OutputRealType,
class InputComplexType>
55struct ConvertValue<OutputRealType, InputComplexType, false, true> {
56 static KOKKOS_INLINE_FUNCTION
void
57 convert(OutputRealType& dst,
58 const InputComplexType& src) {
61 using KAI = KokkosKernels::ArithTraits<InputComplexType>;
62 dst = OutputRealType(KAI::real(src));
66template <
class OutputComplexType,
class InputRealType>
67struct ConvertValue<OutputComplexType, InputRealType, true, false> {
68 static KOKKOS_INLINE_FUNCTION
void
69 convert(OutputComplexType& dst,
70 const InputRealType& src) {
73 using output_mag_type =
74 typename KokkosKernels::ArithTraits<OutputComplexType>::mag_type;
75 using KAM = KokkosKernels::ArithTraits<output_mag_type>;
76 dst = OutputComplexType(src, KAM::zero());
80template <
class OutputValueType,
82KOKKOS_INLINE_FUNCTION
void
83convertValue(OutputValueType& dst,
const InputValueType& src) {
84 ConvertValue<OutputValueType, InputValueType>::convert(dst, src);
91template <
class OutputViewType,
93 const int rank =
static_cast<int>(OutputViewType::rank)>
94class CopyConvertFunctor {};
96template <
class OutputViewType,
98class CopyConvertFunctor<OutputViewType, InputViewType, 1> {
100 static_assert(
static_cast<int>(OutputViewType::rank) == 1 &&
101 static_cast<int>(InputViewType::rank) == 1,
102 "CopyConvertFunctor (implements Tpetra::Details::copyConvert): "
103 "OutputViewType and InputViewType must both have rank 1.");
108 using index_type =
typename OutputViewType::size_type;
110 CopyConvertFunctor(
const OutputViewType& dst,
111 const InputViewType& src)
115 KOKKOS_INLINE_FUNCTION
void
116 operator()(
const index_type i)
const {
117 convertValue(dst_(i), src_(i));
121template <
class OutputViewType,
123class CopyConvertFunctor<OutputViewType, InputViewType, 2> {
125 using index_type =
typename OutputViewType::size_type;
128 static_assert(
static_cast<int>(OutputViewType::rank) == 2 &&
129 static_cast<int>(InputViewType::rank) == 2,
130 "CopyConvertFunctor (implements Tpetra::Details::copyConvert): "
131 "OutputViewType and InputViewType must both have rank 2.");
137 CopyConvertFunctor(
const OutputViewType& dst,
138 const InputViewType& src)
141 , numCols_(dst.extent(1)) {}
143 KOKKOS_INLINE_FUNCTION
void
144 operator()(
const index_type i)
const {
145 const index_type numCols = numCols_;
146 for (index_type j = 0; j < numCols; ++j) {
147 convertValue(dst_(i, j), src_(i, j));
153template <
class OutputViewType,
class InputViewType>
154class CanUseKokkosDeepCopy {
156 static constexpr bool sameValueType =
157 std::is_same<
typename OutputViewType::non_const_value_type,
158 typename InputViewType::non_const_value_type>::value;
159 static constexpr bool sameMemorySpace =
160 std::is_same<
typename OutputViewType::memory_space,
161 typename InputViewType::memory_space>::value;
162 static constexpr bool sameLayout =
163 std::is_same<
typename OutputViewType::array_layout,
164 typename InputViewType::array_layout>::value;
167 static constexpr bool value =
168 sameValueType && (sameMemorySpace || sameLayout);
189template <
class OutputViewType,
191 const bool canUseKokkosDeepCopy =
192 CanUseKokkosDeepCopy<OutputViewType, InputViewType>::value,
193 const bool outputExecSpaceCanAccessInputMemSpace =
194 Kokkos::SpaceAccessibility<
195 typename OutputViewType::memory_space,
196 typename InputViewType::memory_space>::accessible>
197struct CopyConvertImpl {
199 run(
const OutputViewType& dst,
200 const InputViewType& src);
204template <
class OutputViewType,
206 const bool outputExecSpaceCanAccessInputMemSpace>
207struct CopyConvertImpl<OutputViewType, InputViewType,
208 true, outputExecSpaceCanAccessInputMemSpace> {
210 run(
const OutputViewType& dst,
211 const InputViewType& src) {
217 const ptrdiff_t dst_beg =
reinterpret_cast<ptrdiff_t
>(dst.data());
218 const ptrdiff_t dst_end =
219 reinterpret_cast<ptrdiff_t
>(dst.data() + dst.span());
220 const ptrdiff_t src_beg =
reinterpret_cast<ptrdiff_t
>(src.data());
221 const ptrdiff_t src_end =
222 reinterpret_cast<ptrdiff_t
>(src.data() + src.span());
224 if (dst_end > src_beg && src_end > dst_beg) {
232 auto src_copy = Kokkos::create_mirror(Kokkos::HostSpace(), src);
234 Kokkos::deep_copy(src_copy, src);
236 Kokkos::deep_copy(dst, src_copy);
239 using execution_space =
typename OutputViewType::execution_space;
240 Kokkos::deep_copy(execution_space(), dst, src);
247template <
class OutputViewType,
249struct CopyConvertImpl<OutputViewType,
254 run(
const OutputViewType& dst,
255 const InputViewType& src) {
256 using functor_type = CopyConvertFunctor<OutputViewType, InputViewType>;
257 using execution_space =
typename OutputViewType::execution_space;
258 using index_type =
typename OutputViewType::size_type;
259 using range_type = Kokkos::RangePolicy<execution_space, index_type>;
260 Kokkos::parallel_for(
"Tpetra::Details::copyConvert",
261 range_type(0, dst.extent(0)),
262 functor_type(dst, src));
272template <
class OutputViewType,
274struct CopyConvertImpl<OutputViewType, InputViewType, false, false> {
276 run(
const OutputViewType& dst,
277 const InputViewType& src) {
278 using output_memory_space =
typename OutputViewType::memory_space;
279 using output_execution_space =
typename OutputViewType::execution_space;
280 auto src_outputSpaceCopy =
281 Kokkos::create_mirror_view(output_memory_space(), src);
283 Kokkos::deep_copy(output_execution_space(), src_outputSpaceCopy, src);
287 using output_space_copy_type =
decltype(src_outputSpaceCopy);
289 CopyConvertFunctor<OutputViewType, output_space_copy_type>;
290 using execution_space =
typename OutputViewType::execution_space;
291 using index_type =
typename OutputViewType::size_type;
292 using range_type = Kokkos::RangePolicy<execution_space, index_type>;
293 Kokkos::parallel_for(
"Tpetra::Details::copyConvert",
294 range_type(0, dst.extent(0)),
295 functor_type(dst, src_outputSpaceCopy));
308template <
class OutputViewType,
312 static_assert(Kokkos::is_view<OutputViewType>::value,
313 "OutputViewType must be a Kokkos::View.");
314 static_assert(Kokkos::is_view<InputViewType>::value,
315 "InputViewType must be a Kokkos::View.");
316 static_assert(std::is_same<
typename OutputViewType::value_type,
317 typename OutputViewType::non_const_value_type>::value,
318 "OutputViewType must be a nonconst Kokkos::View.");
319 static_assert(
static_cast<int>(OutputViewType::rank) ==
320 static_cast<int>(InputViewType::rank),
321 "src and dst must have the same rank.");
323 if (dst.extent(0) != src.extent(0)) {
324 std::ostringstream
os;
325 os <<
"Tpetra::Details::copyConvert: "
326 <<
"dst.extent(0) = " << dst.extent(0)
327 <<
" != src.extent(0) = " << src.extent(0)
329 throw std::invalid_argument(
os.str());
331 if (
static_cast<int>(OutputViewType::rank) > 1 &&
332 dst.extent(1) != src.extent(1)) {
333 std::ostringstream
os;
334 os <<
"Tpetra::Details::copyConvert: "
335 <<
"dst.extent(1) = " << dst.extent(1)
336 <<
" != src.extent(1) = " << src.extent(1)
338 throw std::invalid_argument(
os.str());
342 using output_view_type =
343 Kokkos::View<
typename OutputViewType::non_const_data_type,
344 typename OutputViewType::array_layout,
345 typename OutputViewType::device_type>;
346 using input_view_type =
347 Kokkos::View<
typename InputViewType::const_data_type,
348 typename InputViewType::array_layout,
349 typename InputViewType::device_type>;
350 CopyConvertImpl<output_view_type, input_view_type>::run(dst, src);
Struct that holds views of the contents of a CrsMatrix.
Implementation details of Tpetra.
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....
Namespace Tpetra contains the class and methods constituting the Tpetra library.