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 {
55  typedef typename device_type::execution_space execution_space;
56  typedef Kokkos::RangePolicy<execution_space, LO> policy_type;
57 
58  lclNumRows_ = ptr.extent(0)-1;
59  policy_type range (0, ptr.extent(0));
60  Kokkos::parallel_for (range, *this);
61 }
62 
63 template<class LO,
64  class GO,
65  class Node,
66  class OffsetType>
67 KOKKOS_FUNCTION void
69 operator() (const LO& lclRowInd) const
70 {
71  const LO INVALID =
72  Tpetra::Details::OrdinalTraits<LO>::invalid ();
73 
74  if (lclRowInd == lclNumRows_)
75  OffRankOffsets_[lclRowInd] = ptr_[lclRowInd];
76  else {
77  // TODO: use parallel reduce
78  size_t offset = ptr_[lclRowInd+1];
79  for (size_t j = ptr_[lclRowInd]; j < ptr_[lclRowInd+1]; j++) {
80  const LO lclColInd = ind_[j];
81  const GO gblColInd = lclColMap_.getGlobalElement (lclColInd);
82  const LO lclDomInd = lclDomMap_.getLocalElement (gblColInd);
83  if ((lclDomInd == INVALID) && (j < offset))
84  offset = j;
85  }
86  OffRankOffsets_[lclRowInd] = offset;
87  }
88 }
89 
90 } // namespace Impl
91 } // namespace Details
92 } // namespace Tpetra
93 
94 // Explicit template instantiation macro for
95 // Tpetra::Details::Impl::GetGraphOffRankOffsets. NOT FOR USERS!!! Must
96 // be used inside the Tpetra namespace.
97 #define TPETRA_DETAILS_IMPL_GETGRAPHOFFRANKOFFSETS_INSTANT( LO, GO, NODE ) \
98  template class Details::Impl::GetGraphOffRankOffsets< LO, GO, NODE::device_type >;
99 
100 #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.