Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_computeOffsets.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef TPETRA_DETAILS_COMPUTEOFFSETS_HPP
11 #define TPETRA_DETAILS_COMPUTEOFFSETS_HPP
12 
19 
20 #include "TpetraCore_config.h"
22 #include <limits>
23 #include <type_traits>
24 
25 namespace Tpetra {
26 namespace Details {
27 
28 //
29 // Implementation details for computeOffsetsFromCounts (see below).
30 // Users should skip over this anonymous namespace.
31 //
32 namespace { // (anonymous)
33 
43 template <class OffsetType,
44  class CountType,
45  class SizeType>
46 class ComputeOffsetsFromCounts {
47  public:
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.");
54 
55  using offsets_view_type =
56  Kokkos::View<OffsetType*, Kokkos::AnonymousSpace>;
57  using counts_view_type =
58  Kokkos::View<const CountType*, Kokkos::AnonymousSpace>;
59 
65  ComputeOffsetsFromCounts(const offsets_view_type& offsets,
66  const counts_view_type& counts)
67  : offsets_(offsets)
68  , counts_(counts)
69  , size_(counts.extent(0)) {}
70 
72  KOKKOS_INLINE_FUNCTION void
73  operator()(const SizeType i, OffsetType& update,
74  const bool finalPass) const {
75  const auto curVal = (i < size_) ? counts_[i] : OffsetType();
76  if (finalPass) {
77  offsets_[i] = update;
78  }
79  update += (i < size_) ? curVal : OffsetType();
80  }
81 
82  template <class ExecutionSpace>
83  static OffsetType
84  run(const ExecutionSpace& execSpace,
85  const offsets_view_type& offsets,
86  const counts_view_type& counts) {
87  const SizeType numCounts(counts.extent(0));
88  using range_type = Kokkos::RangePolicy<ExecutionSpace, SizeType>;
89  range_type range(execSpace, 0, numCounts + SizeType(1));
90  using functor_type =
91  ComputeOffsetsFromCounts<OffsetType, CountType, SizeType>;
92  functor_type functor(offsets, counts);
93  OffsetType total(0);
94  const char funcName[] = "Tpetra::Details::computeOffsetsFromCounts";
95  Kokkos::parallel_scan(funcName, range, functor, total);
96  return total;
97  }
98 
99  private:
101  offsets_view_type offsets_;
103  counts_view_type counts_;
105  SizeType size_;
106 };
107 
117 template <class OffsetType,
118  class CountType,
119  class SizeType>
120 class ComputeOffsetsFromConstantCount {
121  public:
122  static_assert(std::is_integral<OffsetType>::value,
123  "OffsetType must be a built-in integer.");
124  static_assert(std::is_integral<CountType>::value,
125  "CountType must be a built-in integer.");
126  static_assert(std::is_integral<SizeType>::value,
127  "SizeType must be a built-in integer.");
128 
129  using offsets_view_type =
130  Kokkos::View<OffsetType*, Kokkos::AnonymousSpace>;
131 
137  ComputeOffsetsFromConstantCount(const offsets_view_type& offsets,
138  const CountType count)
139  : offsets_(offsets)
140  , count_(count) {}
141 
143  KOKKOS_INLINE_FUNCTION void
144  operator()(const SizeType i) const {
145  offsets_[i] = count_ * i;
146  }
147 
148  template <class ExecutionSpace>
149  static OffsetType
150  run(const ExecutionSpace& execSpace,
151  const offsets_view_type& offsets,
152  const CountType count) {
153  const SizeType numOffsets(offsets.extent(0));
154  if (numOffsets == SizeType(0)) {
155  // Special case that is possible with zero rows
156  return 0;
157  }
158  using range_type = Kokkos::RangePolicy<ExecutionSpace, SizeType>;
159  range_type range(execSpace, 0, numOffsets);
160  using functor_type =
161  ComputeOffsetsFromConstantCount<OffsetType, CountType, SizeType>;
162  functor_type functor(offsets, count);
163  const OffsetType total = (numOffsets - 1) * count;
164  const char funcName[] =
165  "Tpetra::Details::computeOffsetsFromConstantCount";
166  Kokkos::parallel_for(funcName, range, functor);
167  return total;
168  }
169 
170  private:
172  offsets_view_type offsets_;
174  CountType count_;
175 };
176 
177 } // namespace
178 
204 template <class ExecutionSpace,
205  class OffsetsViewType,
206  class CountsViewType,
207  class SizeType = typename OffsetsViewType::size_type>
208 typename OffsetsViewType::non_const_value_type
209 computeOffsetsFromCounts(const ExecutionSpace& execSpace,
210  const OffsetsViewType& ptr,
211  const CountsViewType& counts) {
212  static_assert(Kokkos::is_execution_space<ExecutionSpace>::value,
213  "ExecutionSpace must be a Kokkos execution space.");
214  static_assert(Kokkos::is_view<OffsetsViewType>::value,
215  "OffsetsViewType (the type of ptr) must be a Kokkos::View.");
216  static_assert(Kokkos::is_view<CountsViewType>::value,
217  "CountsViewType (the type of counts) must be a Kokkos::View.");
218  static_assert(std::is_same<typename OffsetsViewType::value_type,
219  typename OffsetsViewType::non_const_value_type>::value,
220  "OffsetsViewType (the type of ptr) must be a nonconst Kokkos::View.");
221  static_assert(static_cast<int>(OffsetsViewType::rank) == 1,
222  "OffsetsViewType (the type of ptr) must be a rank-1 Kokkos::View.");
223  static_assert(static_cast<int>(CountsViewType::rank) == 1,
224  "CountsViewType (the type of counts) must be a rank-1 Kokkos::View.");
225 
226  using offset_type = typename OffsetsViewType::non_const_value_type;
227  static_assert(std::is_integral<offset_type>::value,
228  "The entries of ptr must be built-in integers.");
229  using count_type = typename CountsViewType::non_const_value_type;
230  static_assert(std::is_integral<count_type>::value,
231  "The entries of counts must be built-in integers.");
232  static_assert(std::is_integral<SizeType>::value,
233  "SizeType must be a built-in integer type.");
234 
235  const char funcName[] = "Tpetra::Details::computeOffsetsFromCounts";
236 
237  const auto numOffsets = ptr.size();
238  const auto numCounts = counts.size();
239  offset_type total(0);
240 
241  if (numOffsets != 0) {
242  TEUCHOS_TEST_FOR_EXCEPTION(numCounts >= numOffsets, std::invalid_argument, funcName << ": counts.size() = " << numCounts << " >= ptr.size() = " << numOffsets << ".");
243 
244  using Kokkos::AnonymousSpace;
245  using Kokkos::View;
246  View<offset_type*, AnonymousSpace> ptr_a = ptr;
247  View<const count_type*, AnonymousSpace> counts_a;
248 
249  using offsets_device_type = typename OffsetsViewType::device_type;
250  using counts_copy_type = View<count_type*, offsets_device_type>;
251  counts_copy_type counts_copy;
252 
253  using offsets_memory_space =
254  typename offsets_device_type::memory_space;
255  using counts_memory_space = typename CountsViewType::memory_space;
256  constexpr bool countsAccessibleFromOffsetsExecSpace =
257  Kokkos::SpaceAccessibility<
258  offsets_memory_space, counts_memory_space>::accessible;
259  if (countsAccessibleFromOffsetsExecSpace) {
260  // NOTE (mfh 21 Aug 2019) Some compilers have trouble deducing
261  // that operator= works if more than one template argument
262  // differ. If that should happen, introduce an intermediate
263  // type here.
264  counts_a = counts;
265  } else {
266  using Kokkos::view_alloc;
267  using Kokkos::WithoutInitializing;
268  counts_copy = counts_copy_type(view_alloc("counts_copy", WithoutInitializing), numCounts);
269  Kokkos::deep_copy(execSpace, counts_copy, counts);
270  counts_a = counts_copy;
271  }
272 
273  using functor_type =
274  ComputeOffsetsFromCounts<offset_type, count_type, SizeType>;
275  total = functor_type::run(execSpace, ptr_a, counts_a);
276  }
277 
278  return total;
279 }
280 
282 template <class OffsetsViewType,
283  class CountsViewType,
284  class SizeType = typename OffsetsViewType::size_type>
285 typename OffsetsViewType::non_const_value_type
286 computeOffsetsFromCounts(const OffsetsViewType& ptr,
287  const CountsViewType& counts) {
288  using execution_space = typename OffsetsViewType::execution_space;
289  return computeOffsetsFromCounts(execution_space(), ptr, counts);
290 }
291 
314 template <class OffsetsViewType,
315  class CountType,
316  class SizeType = typename OffsetsViewType::size_type>
317 typename OffsetsViewType::non_const_value_type
318 computeOffsetsFromConstantCount(const OffsetsViewType& ptr,
319  const CountType count) {
320  static_assert(Kokkos::is_view<OffsetsViewType>::value,
321  "ptr must be a Kokkos::View.");
322  static_assert(std::is_same<typename OffsetsViewType::value_type,
323  typename OffsetsViewType::non_const_value_type>::value,
324  "ptr must be a nonconst Kokkos::View.");
325  static_assert(static_cast<int>(OffsetsViewType::rank) == 1,
326  "ptr must be a rank-1 Kokkos::View.");
327 
328  using offset_type = typename OffsetsViewType::non_const_value_type;
329  static_assert(std::is_integral<offset_type>::value,
330  "The type of each entry of ptr must be a "
331  "built-in integer.");
332  static_assert(std::is_integral<CountType>::value,
333  "CountType must be a built-in integer.");
334  static_assert(std::is_integral<SizeType>::value,
335  "SizeType must be a built-in integer.");
336 
337  using device_type = typename OffsetsViewType::device_type;
338  using execution_space = typename device_type::execution_space;
339 
340  offset_type total(0);
341  if (ptr.extent(0) != 0) {
342  using CT = CountType;
343  using functor_type =
344  ComputeOffsetsFromConstantCount<offset_type, CT, SizeType>;
345  execution_space execSpace;
346  total = functor_type::run(execSpace, ptr, count);
347  }
348  return total;
349 }
350 
351 } // namespace Details
352 } // namespace Tpetra
353 
354 #endif // TPETRA_DETAILS_COMPUTEOFFSETS_HPP
counts_view_type counts_
Bucket counts (input argument).
OffsetsViewType::non_const_value_type computeOffsetsFromConstantCount(const OffsetsViewType &ptr, const CountType count)
Compute offsets from a constant count.
offsets_view_type offsets_
Offsets (output argument)
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.
CountType count_
&quot;Count&quot; input argument
SizeType size_
Number of entries in counts_.
OffsetsViewType::non_const_value_type computeOffsetsFromCounts(const ExecutionSpace &execSpace, const OffsetsViewType &ptr, const CountsViewType &counts)
Compute offsets from counts.
Declaration and definition of Tpetra::Details::getEntryOnHost.