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)) {}
 
   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();
 
   79    update += (i < size_) ? curVal : OffsetType();
 
   82  template <
class ExecutionSpace>
 
   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));
 
   91        ComputeOffsetsFromCounts<OffsetType, CountType, SizeType>;
 
   92    functor_type functor(offsets, counts);
 
   94    const char funcName[] = 
"Tpetra::Details::computeOffsetsFromCounts";
 
   95    Kokkos::parallel_scan(funcName, range, functor, total);
 
  101  offsets_view_type offsets_;
 
  103  counts_view_type counts_;
 
  117template <
class OffsetType,
 
  120class ComputeOffsetsFromConstantCount {
 
  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.");
 
  129  using offsets_view_type =
 
  130      Kokkos::View<OffsetType*, Kokkos::AnonymousSpace>;
 
  137  ComputeOffsetsFromConstantCount(
const offsets_view_type& offsets,
 
  138                                  const CountType count)
 
  143  KOKKOS_INLINE_FUNCTION 
void 
  144  operator()(
const SizeType i)
 const {
 
  145    offsets_[i] = count_ * i;
 
  148  template <
class ExecutionSpace>
 
  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)) {
 
  158    using range_type = Kokkos::RangePolicy<ExecutionSpace, SizeType>;
 
  159    range_type range(execSpace, 0, numOffsets);
 
  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);
 
  172  offsets_view_type offsets_;
 
  204template <
class ExecutionSpace,
 
  205          class OffsetsViewType,
 
  206          class CountsViewType,
 
  207          class SizeType = 
typename OffsetsViewType::size_type>
 
  208typename OffsetsViewType::non_const_value_type
 
  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.");
 
  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.");
 
  235  const char funcName[] = 
"Tpetra::Details::computeOffsetsFromCounts";
 
  239  offset_type 
total(0);
 
  244    using Kokkos::AnonymousSpace;
 
  254        typename offsets_device_type::memory_space;
 
  257        Kokkos::SpaceAccessibility<
 
  266      using Kokkos::view_alloc;
 
  267      using Kokkos::WithoutInitializing;
 
 
  284          class SizeType = 
typename OffsetsViewType::size_type>
 
  285typename OffsetsViewType::non_const_value_type
 
  288  using execution_space = 
typename OffsetsViewType::execution_space;
 
 
  316          class SizeType = 
typename OffsetsViewType::size_type>
 
  317typename OffsetsViewType::non_const_value_type
 
  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.");
 
  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.");
 
  337  using device_type     = 
typename OffsetsViewType::device_type;
 
  338  using execution_space = 
typename device_type::execution_space;
 
  340  offset_type 
total(0);
 
  341  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.