Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_getGraphOffRankOffsets_def.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_GETGRAPHOFFRANKOFFSETS_DEF_HPP
11 #define TPETRA_DETAILS_GETGRAPHOFFRANKOFFSETS_DEF_HPP
12 
17 
19 #include "Tpetra_Map.hpp"
20 #include "KokkosSparse_findRelOffset.hpp"
21 
22 namespace Tpetra {
23 namespace Details {
24 namespace Impl {
25 
39 template <class LO,
40  class GO,
41  class Node,
42  class OffsetType>
44  GetGraphOffRankOffsets(const offsets_type& OffRankOffsets,
45  const local_map_type& lclColMap,
46  const local_map_type& lclDomMap,
47  const row_offsets_type& ptr,
48  const lcl_col_inds_type& ind)
49  : OffRankOffsets_(OffRankOffsets)
50  , lclColMap_(lclColMap)
51  , lclDomMap_(lclDomMap)
52  , ptr_(ptr)
53  , ind_(ind) {
54  typedef typename device_type::execution_space execution_space;
55  typedef Kokkos::RangePolicy<execution_space, LO> policy_type;
56 
57  lclNumRows_ = ptr.extent(0) - 1;
58  policy_type range(0, ptr.extent(0));
59  Kokkos::parallel_for(range, *this);
60 }
61 
62 template <class LO,
63  class GO,
64  class Node,
65  class OffsetType>
66 KOKKOS_FUNCTION void
68 operator()(const LO& lclRowInd) const {
69  const LO INVALID =
70  Tpetra::Details::OrdinalTraits<LO>::invalid();
71 
72  if (lclRowInd == lclNumRows_)
73  OffRankOffsets_[lclRowInd] = ptr_[lclRowInd];
74  else {
75  // TODO: use parallel reduce
76  size_t offset = ptr_[lclRowInd + 1];
77  for (size_t j = ptr_[lclRowInd]; j < ptr_[lclRowInd + 1]; j++) {
78  const LO lclColInd = ind_[j];
79  const GO gblColInd = lclColMap_.getGlobalElement(lclColInd);
80  const LO lclDomInd = lclDomMap_.getLocalElement(gblColInd);
81  if ((lclDomInd == INVALID) && (j < offset))
82  offset = j;
83  }
84  OffRankOffsets_[lclRowInd] = offset;
85  }
86 }
87 
88 } // namespace Impl
89 } // namespace Details
90 } // namespace Tpetra
91 
92 // Explicit template instantiation macro for
93 // Tpetra::Details::Impl::GetGraphOffRankOffsets. NOT FOR USERS!!! Must
94 // be used inside the Tpetra namespace.
95 #define TPETRA_DETAILS_IMPL_GETGRAPHOFFRANKOFFSETS_INSTANT(LO, GO, NODE) \
96  template class Details::Impl::GetGraphOffRankOffsets<LO, GO, NODE::device_type>;
97 
98 #endif // TPETRA_DETAILS_GETGRAPHOFFRANKOFFSETS_DEF_HPP
Import KokkosSparse::OrdinalTraits, a traits class for &quot;invalid&quot; (flag) values of integer types...
GetGraphOffRankOffsets(const offsets_type &OffRankOffsets, const local_map_type &lclColMap, const local_map_type &lclDomMap, const row_offsets_type &ptr, const lcl_col_inds_type &ind)
Constructor; also runs the functor.
KOKKOS_FUNCTION void operator()(const LO &lclRowInd) const
Kokkos::parallel_for loop body.