Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
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
25namespace Tpetra {
26namespace Details {
27
28//
29// Implementation details for computeOffsetsFromCounts (see below).
30// Users should skip over this anonymous namespace.
31//
32namespace { // (anonymous)
33
43template<class OffsetType,
44 class CountType,
45 class SizeType>
46class ComputeOffsetsFromCounts {
47public:
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 {}
71
73 KOKKOS_INLINE_FUNCTION void
74 operator () (const SizeType i, OffsetType& update,
75 const bool finalPass) const
76 {
77 const auto curVal = (i < size_) ? counts_[i] : OffsetType ();
78 if (finalPass) {
79 offsets_[i] = update;
80 }
81 update += (i < size_) ? curVal : OffsetType ();
82 }
83
84 template<class ExecutionSpace>
85 static OffsetType
86 run (const ExecutionSpace& execSpace,
87 const offsets_view_type& offsets,
88 const counts_view_type& counts)
89 {
90 const SizeType numCounts (counts.extent (0));
91 using range_type = Kokkos::RangePolicy<ExecutionSpace, SizeType>;
92 range_type range (execSpace, 0, numCounts + SizeType (1));
93 using functor_type =
94 ComputeOffsetsFromCounts<OffsetType, CountType, SizeType>;
95 functor_type functor (offsets, counts);
96 OffsetType total (0);
97 const char funcName[] = "Tpetra::Details::computeOffsetsFromCounts";
98 Kokkos::parallel_scan (funcName, range, functor, total);
99 return total;
100 }
101
102private:
104 offsets_view_type offsets_;
106 counts_view_type counts_;
108 SizeType size_;
109};
110
120template<class OffsetType,
121 class CountType,
122 class SizeType>
123class ComputeOffsetsFromConstantCount {
124public:
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.");
131
132 using offsets_view_type =
133 Kokkos::View<OffsetType*, Kokkos::AnonymousSpace>;
134
140 ComputeOffsetsFromConstantCount (const offsets_view_type& offsets,
141 const CountType count) :
142 offsets_ (offsets),
143 count_ (count)
144 {}
145
147 KOKKOS_INLINE_FUNCTION void
148 operator () (const SizeType i) const
149 {
150 offsets_[i] = count_*i;
151 }
152
153 template<class ExecutionSpace>
154 static OffsetType
155 run (const ExecutionSpace& execSpace,
156 const offsets_view_type& offsets,
157 const CountType count)
158 {
159 const SizeType numOffsets (offsets.extent (0));
160 if(numOffsets == SizeType(0))
161 {
162 // Special case that is possible with zero rows
163 return 0;
164 }
165 using range_type = Kokkos::RangePolicy<ExecutionSpace, SizeType>;
166 range_type range (execSpace, 0, numOffsets);
167 using functor_type =
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);
174 return total;
175 }
176
177private:
179 offsets_view_type offsets_;
181 CountType count_;
182};
183
184} // namespace (anonymous)
185
211template<class ExecutionSpace,
212 class OffsetsViewType,
213 class CountsViewType,
214 class SizeType = typename OffsetsViewType::size_type>
215typename OffsetsViewType::non_const_value_type
216computeOffsetsFromCounts (const ExecutionSpace& execSpace,
217 const OffsetsViewType& ptr,
218 const CountsViewType& counts)
219{
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.");
233
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.");
242
243 const char funcName[] = "Tpetra::Details::computeOffsetsFromCounts";
244
245 const auto numOffsets = ptr.size ();
246 const auto numCounts = counts.size ();
247 offset_type total (0);
248
249 if (numOffsets != 0) {
250 TEUCHOS_TEST_FOR_EXCEPTION
251 (numCounts >= numOffsets, std::invalid_argument, funcName <<
252 ": counts.size() = " << numCounts << " >= ptr.size() = " <<
253 numOffsets << ".");
254
255 using Kokkos::AnonymousSpace;
256 using Kokkos::View;
257 View<offset_type*, AnonymousSpace> ptr_a = ptr;
258 View<const count_type*, AnonymousSpace> counts_a;
259
260 using offsets_device_type = typename OffsetsViewType::device_type;
261 using counts_copy_type = View<count_type*, offsets_device_type>;
262 counts_copy_type counts_copy;
263
264 using offsets_memory_space =
265 typename offsets_device_type::memory_space;
266 using counts_memory_space = typename CountsViewType::memory_space;
267 constexpr bool countsAccessibleFromOffsetsExecSpace =
268 Kokkos::SpaceAccessibility<
269 offsets_memory_space, counts_memory_space>::accessible;
270 if (countsAccessibleFromOffsetsExecSpace) {
271 // NOTE (mfh 21 Aug 2019) Some compilers have trouble deducing
272 // that operator= works if more than one template argument
273 // differ. If that should happen, introduce an intermediate
274 // type here.
275 counts_a = counts;
276 }
277 else {
278 using Kokkos::view_alloc;
279 using Kokkos::WithoutInitializing;
280 counts_copy = counts_copy_type
281 (view_alloc ("counts_copy", WithoutInitializing), numCounts);
282 Kokkos::deep_copy (execSpace, counts_copy, counts);
283 counts_a = counts_copy;
284 }
285
286 using functor_type =
287 ComputeOffsetsFromCounts<offset_type, count_type, SizeType>;
288 total = functor_type::run (execSpace, ptr_a, counts_a);
289 }
290
291 return total;
292}
293
295template<class OffsetsViewType,
296 class CountsViewType,
297 class SizeType = typename OffsetsViewType::size_type>
298typename OffsetsViewType::non_const_value_type
299computeOffsetsFromCounts (const OffsetsViewType& ptr,
300 const CountsViewType& counts)
301{
302 using execution_space = typename OffsetsViewType::execution_space;
303 return computeOffsetsFromCounts (execution_space (), ptr, counts);
304}
305
328template<class OffsetsViewType,
329 class CountType,
330 class SizeType = typename OffsetsViewType::size_type>
331typename OffsetsViewType::non_const_value_type
332computeOffsetsFromConstantCount (const OffsetsViewType& ptr,
333 const CountType count)
334{
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.");
342
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.");
351
352 using device_type = typename OffsetsViewType::device_type;
353 using execution_space = typename device_type::execution_space;
354
355 offset_type total (0);
356 if (ptr.extent (0) != 0) {
357 using CT = CountType;
358 using functor_type =
359 ComputeOffsetsFromConstantCount<offset_type, CT, SizeType>;
360 execution_space execSpace;
361 total = functor_type::run (execSpace, ptr, count);
362 }
363 return total;
364}
365
366} // namespace Details
367} // namespace Tpetra
368
369#endif // TPETRA_DETAILS_COMPUTEOFFSETS_HPP
Declaration and definition of Tpetra::Details::getEntryOnHost.
Nonmember function that computes a residual Computes R = B - A * X.
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.