Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_getGraphDiagOffsets_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_GETGRAPHDIAGOFFSETS_DEF_HPP
11 #define TPETRA_DETAILS_GETGRAPHDIAGOFFSETS_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 DiagOffsetType>
44  GetGraphDiagOffsets(const diag_offsets_type& diagOffsets,
45  const local_map_type& lclRowMap,
46  const local_map_type& lclColMap,
47  const row_offsets_type& ptr,
48  const lcl_col_inds_type& ind,
49  const bool isSorted)
50  : diagOffsets_(diagOffsets)
51  , lclRowMap_(lclRowMap)
52  , lclColMap_(lclColMap)
53  , ptr_(ptr)
54  , ind_(ind)
55  , isSorted_(isSorted) {
56  typedef typename device_type::execution_space execution_space;
57  typedef Kokkos::RangePolicy<execution_space, LO> policy_type;
58 
59  const LO lclNumRows = lclRowMap.getLocalNumElements();
60  policy_type range(0, lclNumRows);
61  Kokkos::parallel_for(range, *this);
62 }
63 
64 template <class LO,
65  class GO,
66  class Node,
67  class DiagOffsetType>
68 KOKKOS_FUNCTION void
70 operator()(const LO& lclRowInd) const {
71  const size_t STINV =
72  Tpetra::Details::OrdinalTraits<diag_offset_type>::invalid();
73  const GO gblRowInd = lclRowMap_.getGlobalElement(lclRowInd);
74  const GO gblColInd = gblRowInd;
75  const LO lclColInd = lclColMap_.getLocalElement(gblColInd);
76 
77  if (lclColInd == Tpetra::Details::OrdinalTraits<LO>::invalid()) {
78  diagOffsets_[lclRowInd] = STINV;
79  } else {
80  // Could be empty, but that's OK.
81  const LO numEnt = ptr_[lclRowInd + 1] - ptr_[lclRowInd];
82  // std::pair doesn't have its methods marked as device
83  // functions, so we have to use Kokkos::pair.
84  auto lclColInds =
85  Kokkos::subview(ind_, Kokkos::make_pair(ptr_[lclRowInd],
86  ptr_[lclRowInd + 1]));
87  using ::KokkosSparse::findRelOffset;
88  const LO diagOffset =
89  findRelOffset<LO, lcl_col_inds_type>(lclColInds, numEnt,
90  lclColInd, 0, isSorted_);
91  diagOffsets_[lclRowInd] = (diagOffset == numEnt) ? STINV : static_cast<diag_offset_type>(diagOffset);
92  }
93 }
94 
95 } // namespace Impl
96 } // namespace Details
97 } // namespace Tpetra
98 
99 // Explicit template instantiation macro for
100 // Tpetra::Details::Impl::GetGraphDiagOffsets. NOT FOR USERS!!! Must
101 // be used inside the Tpetra namespace.
102 #define TPETRA_DETAILS_IMPL_GETGRAPHDIAGOFFSETS_INSTANT(LO, GO, NODE) \
103  template class Details::Impl::GetGraphDiagOffsets<LO, GO, NODE::device_type>;
104 
105 #endif // TPETRA_DETAILS_GETGRAPHDIAGOFFSETS_DEF_HPP
KOKKOS_FUNCTION void operator()(const LO &lclRowInd) const
Kokkos::parallel_for loop body.
Import KokkosSparse::OrdinalTraits, a traits class for &quot;invalid&quot; (flag) values of integer types...
KOKKOS_INLINE_FUNCTION LocalOrdinal getLocalNumElements() const
The number of indices that live on the calling process.
GetGraphDiagOffsets(const diag_offsets_type &diagOffsets, const local_map_type &lclRowMap, const local_map_type &lclColMap, const row_offsets_type &ptr, const lcl_col_inds_type &ind, const bool isSorted)
Constructor; also runs the functor.