42 #ifndef TPETRA_DETAILS_LOCALMAP_HPP
43 #define TPETRA_DETAILS_LOCALMAP_HPP
49 #include "Tpetra_Details_FixedHashTable.hpp"
70 template<
class LocalOrdinal,
class GlobalOrdinal,
class DeviceType>
73 typedef LocalOrdinal local_ordinal_type;
74 typedef GlobalOrdinal global_ordinal_type;
75 typedef DeviceType device_type;
79 myMinGid_ (Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ()),
80 myMaxGid_ (Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ()),
81 firstContiguousGid_ (Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ()),
82 lastContiguousGid_ (Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ()),
83 numLocalElements_ (0),
86 LocalMap (const ::Tpetra::Details::FixedHashTable<GlobalOrdinal, LocalOrdinal, DeviceType>& glMap,
87 const ::Kokkos::View<const GlobalOrdinal*, ::Kokkos::LayoutLeft, DeviceType>& lgMap,
88 const GlobalOrdinal indexBase,
89 const GlobalOrdinal myMinGid,
90 const GlobalOrdinal myMaxGid,
91 const GlobalOrdinal firstContiguousGid,
92 const GlobalOrdinal lastContiguousGid,
93 const LocalOrdinal numLocalElements,
94 const bool contiguous) :
97 indexBase_ (indexBase),
100 firstContiguousGid_ (firstContiguousGid),
101 lastContiguousGid_ (lastContiguousGid),
102 numLocalElements_ (numLocalElements),
103 contiguous_ (contiguous)
108 return numLocalElements_;
130 KOKKOS_INLINE_FUNCTION LocalOrdinal
133 if (numLocalElements_ == 0) {
134 return ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid ();
136 return static_cast<LocalOrdinal
> (numLocalElements_ - 1);
151 KOKKOS_INLINE_FUNCTION LocalOrdinal
155 if (globalIndex < myMinGid_ || globalIndex > myMaxGid_) {
156 return ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid ();
158 return static_cast<LocalOrdinal
> (globalIndex - myMinGid_);
160 else if (globalIndex >= firstContiguousGid_ &&
161 globalIndex <= lastContiguousGid_) {
162 return static_cast<LocalOrdinal
> (globalIndex - firstContiguousGid_);
167 return glMap_.
get (globalIndex);
172 KOKKOS_INLINE_FUNCTION GlobalOrdinal
176 return ::Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ();
182 return lgMap_(localIndex);
203 ::Kokkos::View<const GlobalOrdinal*, ::Kokkos::LayoutLeft, DeviceType> lgMap_;
204 GlobalOrdinal indexBase_;
205 GlobalOrdinal myMinGid_;
206 GlobalOrdinal myMaxGid_;
207 GlobalOrdinal firstContiguousGid_;
208 GlobalOrdinal lastContiguousGid_;
209 LocalOrdinal numLocalElements_;
216 #endif // TPETRA_DETAILS_LOCALMAP_HPP
KOKKOS_INLINE_FUNCTION GlobalOrdinal getMinGlobalIndex() const
The minimum global index on the calling process.
KOKKOS_INLINE_FUNCTION ValueType get(const KeyType &key) const
Get the value corresponding to the given key.
KOKKOS_INLINE_FUNCTION LocalOrdinal getLocalElement(const GlobalOrdinal globalIndex) const
Get the local index corresponding to the given global index.
"Local" 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.
KOKKOS_INLINE_FUNCTION bool isContiguous() const
Whether the Map is (locally) contiguous.
KOKKOS_INLINE_FUNCTION LocalOrdinal getNodeNumElements() const
The number of indices that live on the calling process.
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 GlobalOrdinal getIndexBase() const
The (global) index base.
KOKKOS_INLINE_FUNCTION LocalOrdinal getMaxLocalIndex() const
The maximum local index.