Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_LocalMap.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_LOCALMAP_HPP
11 #define TPETRA_DETAILS_LOCALMAP_HPP
12 
16 
17 #include "Tpetra_Details_FixedHashTable.hpp"
18 // #include "Tpetra_Details_OrdinalTraits.hpp" // comes in above
19 // #include "Kokkos_Core.hpp" // comes in above
21 
22 namespace Tpetra {
23 namespace Details {
24 
38 template<class LocalOrdinal, class GlobalOrdinal, class DeviceType>
39 class LocalMap {
40 public:
42  using local_ordinal_type = LocalOrdinal;
43 
45  using global_ordinal_type = GlobalOrdinal;
46 
48  using device_type = DeviceType;
49 
51  using execution_space = typename device_type::execution_space;
52 
54  using memory_space = typename device_type::memory_space;
55 
57  KOKKOS_DEFAULTED_FUNCTION LocalMap() = default;
58 
59  LocalMap (const ::Tpetra::Details::FixedHashTable<GlobalOrdinal, LocalOrdinal, device_type>& glMap,
60  const ::Kokkos::View<const GlobalOrdinal*, ::Kokkos::LayoutLeft, device_type>& lgMap,
61  const GlobalOrdinal indexBase,
62  const GlobalOrdinal myMinGid,
63  const GlobalOrdinal myMaxGid,
64  const GlobalOrdinal firstContiguousGid,
65  const GlobalOrdinal lastContiguousGid,
66  const LocalOrdinal numLocalElements,
67  const bool contiguous) :
68  glMap_ (glMap),
69  lgMap_ (lgMap),
70  indexBase_ (indexBase),
71  myMinGid_ (myMinGid),
72  myMaxGid_ (myMaxGid),
73  firstContiguousGid_ (firstContiguousGid),
74  lastContiguousGid_ (lastContiguousGid),
75  numLocalElements_ (numLocalElements),
76  contiguous_ (contiguous)
77  {}
78 
80  KOKKOS_INLINE_FUNCTION LocalOrdinal getLocalNumElements () const {
81  return numLocalElements_;
82  }
83 
85  KOKKOS_INLINE_FUNCTION GlobalOrdinal getIndexBase () const {
86  return indexBase_;
87  }
88 
93  KOKKOS_INLINE_FUNCTION bool isContiguous () const {
94  return contiguous_;
95  }
96 
98  KOKKOS_INLINE_FUNCTION LocalOrdinal getMinLocalIndex () const {
99  return 0;
100  }
101 
103  KOKKOS_INLINE_FUNCTION LocalOrdinal
105  {
106  if (numLocalElements_ == 0) {
107  return ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid ();
108  } else { // Local indices are always zero-based.
109  return static_cast<LocalOrdinal> (numLocalElements_ - 1);
110  }
111  }
112 
114  KOKKOS_INLINE_FUNCTION GlobalOrdinal getMinGlobalIndex () const {
115  return myMinGid_;
116  }
117 
119  KOKKOS_INLINE_FUNCTION GlobalOrdinal getMaxGlobalIndex () const {
120  return myMaxGid_;
121  }
122 
124  KOKKOS_INLINE_FUNCTION LocalOrdinal
125  getLocalElement (const GlobalOrdinal globalIndex) const
126  {
127  if (contiguous_) {
128  if (globalIndex < myMinGid_ || globalIndex > myMaxGid_) {
129  return ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid ();
130  }
131  return static_cast<LocalOrdinal> (globalIndex - myMinGid_);
132  }
133  else if (globalIndex >= firstContiguousGid_ &&
134  globalIndex <= lastContiguousGid_) {
135  return static_cast<LocalOrdinal> (globalIndex - firstContiguousGid_);
136  }
137  else {
138  // If the given global index is not in the table, this returns
139  // the same value as OrdinalTraits<LocalOrdinal>::invalid().
140  return glMap_.get (globalIndex);
141  }
142  }
143 
145  KOKKOS_INLINE_FUNCTION GlobalOrdinal
146  getGlobalElement (const LocalOrdinal localIndex) const
147  {
148  if (localIndex < getMinLocalIndex () || localIndex > getMaxLocalIndex ()) {
149  return ::Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ();
150  }
151  if (isContiguous ()) {
152  return getMinGlobalIndex () + localIndex;
153  }
154  else {
155  return lgMap_(localIndex);
156  }
157  }
158 
159 private:
176  ::Kokkos::View<const GlobalOrdinal*, ::Kokkos::LayoutLeft, device_type> lgMap_;
177 
178  GlobalOrdinal indexBase_ = 0;
179  GlobalOrdinal myMinGid_ = Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
180  GlobalOrdinal myMaxGid_ = Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
181  GlobalOrdinal firstContiguousGid_ = Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
182  GlobalOrdinal lastContiguousGid_ = Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
183  LocalOrdinal numLocalElements_ = 0;
184  bool contiguous_ = false;
185 };
186 
187 } // namespace Details
188 } // namespace Tpetra
189 
190 #endif // TPETRA_DETAILS_LOCALMAP_HPP
191 
KOKKOS_INLINE_FUNCTION GlobalOrdinal getMinGlobalIndex() const
The minimum global index on the calling process.
GO global_ordinal_type
The type of global indices.
KOKKOS_INLINE_FUNCTION ValueType get(const KeyType &key) const
Get the value corresponding to the given key.
KOKKOS_DEFAULTED_FUNCTION LocalMap()=default
Default constructor.
KOKKOS_INLINE_FUNCTION LocalOrdinal getLocalElement(const GlobalOrdinal globalIndex) const
Get the local index corresponding to the given global index. (device only)
typename device_type::execution_space execution_space
The Kokkos execution space.
&quot;Local&quot; part of Map suitable for Kokkos kernels.
KOKKOS_INLINE_FUNCTION GlobalOrdinal getGlobalElement(const LocalOrdinal localIndex) const
Get the global index corresponding to the given local index. (device only)
KOKKOS_INLINE_FUNCTION bool isContiguous() const
Whether the Map is (locally) contiguous.
KOKKOS_INLINE_FUNCTION LocalOrdinal getMinLocalIndex() const
The minimum local index.
KOKKOS_INLINE_FUNCTION GlobalOrdinal getMaxGlobalIndex() const
The maximum global index on the calling process.
Forward declaration of Tpetra::Details::LocalMap.
KOKKOS_INLINE_FUNCTION LocalOrdinal getLocalNumElements() const
The number of indices that live on the calling process.
typename device_type::memory_space memory_space
The Kokkos memory space.
KOKKOS_INLINE_FUNCTION GlobalOrdinal getIndexBase() const
The (global) index base.
KOKKOS_INLINE_FUNCTION LocalOrdinal getMaxLocalIndex() const
The maximum local index.
LO local_ordinal_type
The type of local indices.