10#ifndef TPETRA_DETAILS_COMPUTEOFFSETS_HPP
11#define TPETRA_DETAILS_COMPUTEOFFSETS_HPP
20#include "TpetraCore_config.h"
43template<
class OffsetType,
46class ComputeOffsetsFromCounts {
48 static_assert (std::is_integral<OffsetType>::value,
49 "OffsetType must be a built-in integer.");
50 static_assert (std::is_integral<CountType>::value,
51 "CountType must be a built-in integer.");
52 static_assert (std::is_integral<SizeType>::value,
53 "SizeType must be a built-in integer.");
55 using offsets_view_type =
56 Kokkos::View<OffsetType*, Kokkos::AnonymousSpace>;
57 using counts_view_type =
58 Kokkos::View<const CountType*, Kokkos::AnonymousSpace>;
65 ComputeOffsetsFromCounts (
const offsets_view_type& offsets,
66 const counts_view_type& counts) :
69 size_ (counts.extent (0))
73 KOKKOS_INLINE_FUNCTION
void
74 operator () (
const SizeType i, OffsetType& update,
75 const bool finalPass)
const
77 const auto curVal = (i < size_) ? counts_[i] : OffsetType ();
81 update += (i < size_) ? curVal : OffsetType ();
84 template<
class ExecutionSpace>
86 run (
const ExecutionSpace& execSpace,
87 const offsets_view_type& offsets,
88 const counts_view_type& counts)
90 const SizeType numCounts (counts.extent (0));
91 using range_type = Kokkos::RangePolicy<ExecutionSpace, SizeType>;
92 range_type range (execSpace, 0, numCounts + SizeType (1));
94 ComputeOffsetsFromCounts<OffsetType, CountType, SizeType>;
95 functor_type functor (offsets, counts);
97 const char funcName[] =
"Tpetra::Details::computeOffsetsFromCounts";
98 Kokkos::parallel_scan (funcName, range, functor, total);
104 offsets_view_type offsets_;
106 counts_view_type counts_;
120template<
class OffsetType,
123class ComputeOffsetsFromConstantCount {
125 static_assert (std::is_integral<OffsetType>::value,
126 "OffsetType must be a built-in integer.");
127 static_assert (std::is_integral<CountType>::value,
128 "CountType must be a built-in integer.");
129 static_assert (std::is_integral<SizeType>::value,
130 "SizeType must be a built-in integer.");
132 using offsets_view_type =
133 Kokkos::View<OffsetType*, Kokkos::AnonymousSpace>;
140 ComputeOffsetsFromConstantCount (
const offsets_view_type& offsets,
141 const CountType count) :
147 KOKKOS_INLINE_FUNCTION
void
148 operator () (
const SizeType i)
const
150 offsets_[i] = count_*i;
153 template<
class ExecutionSpace>
155 run (
const ExecutionSpace& execSpace,
156 const offsets_view_type& offsets,
157 const CountType count)
159 const SizeType numOffsets (offsets.extent (0));
160 if(numOffsets == SizeType(0))
165 using range_type = Kokkos::RangePolicy<ExecutionSpace, SizeType>;
166 range_type range (execSpace, 0, numOffsets);
168 ComputeOffsetsFromConstantCount<OffsetType, CountType, SizeType>;
169 functor_type functor (offsets, count);
170 const OffsetType total = (numOffsets - 1) * count;
171 const char funcName[] =
172 "Tpetra::Details::computeOffsetsFromConstantCount";
173 Kokkos::parallel_for (funcName, range, functor);
179 offsets_view_type offsets_;
211template<
class ExecutionSpace,
212 class OffsetsViewType,
213 class CountsViewType,
214 class SizeType =
typename OffsetsViewType::size_type>
215typename OffsetsViewType::non_const_value_type
220 static_assert (Kokkos::is_execution_space<ExecutionSpace>::value,
221 "ExecutionSpace must be a Kokkos execution space.");
222 static_assert (Kokkos::is_view<OffsetsViewType>::value,
223 "OffsetsViewType (the type of ptr) must be a Kokkos::View.");
224 static_assert (Kokkos::is_view<CountsViewType>::value,
225 "CountsViewType (the type of counts) must be a Kokkos::View.");
226 static_assert (std::is_same<
typename OffsetsViewType::value_type,
227 typename OffsetsViewType::non_const_value_type>::value,
228 "OffsetsViewType (the type of ptr) must be a nonconst Kokkos::View.");
229 static_assert (
static_cast<int> (OffsetsViewType::rank) == 1,
230 "OffsetsViewType (the type of ptr) must be a rank-1 Kokkos::View.");
231 static_assert (
static_cast<int> (CountsViewType::rank) == 1,
232 "CountsViewType (the type of counts) must be a rank-1 Kokkos::View.");
234 using offset_type =
typename OffsetsViewType::non_const_value_type;
235 static_assert (std::is_integral<offset_type>::value,
236 "The entries of ptr must be built-in integers.");
237 using count_type =
typename CountsViewType::non_const_value_type;
238 static_assert (std::is_integral<count_type>::value,
239 "The entries of counts must be built-in integers.");
240 static_assert (std::is_integral<SizeType>::value,
241 "SizeType must be a built-in integer type.");
243 const char funcName[] =
"Tpetra::Details::computeOffsetsFromCounts";
247 offset_type
total (0);
252 ": counts.size() = " <<
numCounts <<
" >= ptr.size() = " <<
255 using Kokkos::AnonymousSpace;
265 typename offsets_device_type::memory_space;
268 Kokkos::SpaceAccessibility<
278 using Kokkos::view_alloc;
279 using Kokkos::WithoutInitializing;
297 class SizeType =
typename OffsetsViewType::size_type>
298typename OffsetsViewType::non_const_value_type
302 using execution_space =
typename OffsetsViewType::execution_space;
330 class SizeType =
typename OffsetsViewType::size_type>
331typename OffsetsViewType::non_const_value_type
335 static_assert (Kokkos::is_view<OffsetsViewType>::value,
336 "ptr must be a Kokkos::View.");
337 static_assert (std::is_same<
typename OffsetsViewType::value_type,
338 typename OffsetsViewType::non_const_value_type>::value,
339 "ptr must be a nonconst Kokkos::View.");
340 static_assert (
static_cast<int> (OffsetsViewType::rank) == 1,
341 "ptr must be a rank-1 Kokkos::View.");
343 using offset_type =
typename OffsetsViewType::non_const_value_type;
344 static_assert (std::is_integral<offset_type>::value,
345 "The type of each entry of ptr must be a "
346 "built-in integer.");
347 static_assert (std::is_integral<CountType>::value,
348 "CountType must be a built-in integer.");
349 static_assert (std::is_integral<SizeType>::value,
350 "SizeType must be a built-in integer.");
352 using device_type =
typename OffsetsViewType::device_type;
353 using execution_space =
typename device_type::execution_space;
355 offset_type
total (0);
356 if (
ptr.extent (0) != 0) {
Declaration and definition of Tpetra::Details::getEntryOnHost.
Struct that holds views of the contents of a CrsMatrix.
Implementation details of Tpetra.
OffsetsViewType::non_const_value_type computeOffsetsFromCounts(const ExecutionSpace &execSpace, const OffsetsViewType &ptr, const CountsViewType &counts)
Compute offsets from counts.
OffsetsViewType::non_const_value_type computeOffsetsFromConstantCount(const OffsetsViewType &ptr, const CountType count)
Compute offsets from a constant count.
Namespace Tpetra contains the class and methods constituting the Tpetra library.