42 #ifndef TPETRA_DETAILS_LOCALMAP_HPP 43 #define TPETRA_DETAILS_LOCALMAP_HPP 49 #include "Tpetra_Details_FixedHashTable.hpp" 69 template<
class LocalOrdinal,
class GlobalOrdinal,
class DeviceType>
72 typedef LocalOrdinal local_ordinal_type;
73 typedef GlobalOrdinal global_ordinal_type;
74 typedef DeviceType device_type;
76 LocalMap (const ::Tpetra::Details::FixedHashTable<GlobalOrdinal, LocalOrdinal, DeviceType>& glMap,
77 const ::Kokkos::View<const GlobalOrdinal*, ::Kokkos::LayoutLeft, DeviceType>& lgMap,
78 const GlobalOrdinal indexBase,
79 const GlobalOrdinal myMinGid,
80 const GlobalOrdinal myMaxGid,
81 const GlobalOrdinal firstContiguousGid,
82 const GlobalOrdinal lastContiguousGid,
83 const LocalOrdinal numLocalElements,
84 const bool contiguous) :
87 indexBase_ (indexBase),
90 firstContiguousGid_ (firstContiguousGid),
91 lastContiguousGid_ (lastContiguousGid),
92 numLocalElements_ (numLocalElements),
93 contiguous_ (contiguous)
98 return numLocalElements_;
120 KOKKOS_INLINE_FUNCTION LocalOrdinal
123 if (numLocalElements_ == 0) {
124 return ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid ();
126 return static_cast<LocalOrdinal
> (numLocalElements_ - 1);
141 KOKKOS_INLINE_FUNCTION LocalOrdinal
145 if (globalIndex < myMinGid_ || globalIndex > myMaxGid_) {
146 return ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid ();
148 return static_cast<LocalOrdinal
> (globalIndex - myMinGid_);
150 else if (globalIndex >= firstContiguousGid_ &&
151 globalIndex <= lastContiguousGid_) {
152 return static_cast<LocalOrdinal
> (globalIndex - firstContiguousGid_);
157 return glMap_.
get (globalIndex);
162 KOKKOS_INLINE_FUNCTION GlobalOrdinal
166 return ::Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ();
172 return lgMap_(localIndex);
193 ::Kokkos::View<const GlobalOrdinal*, ::Kokkos::LayoutLeft, DeviceType> lgMap_;
194 GlobalOrdinal indexBase_;
195 GlobalOrdinal myMinGid_;
196 GlobalOrdinal myMaxGid_;
197 GlobalOrdinal firstContiguousGid_;
198 GlobalOrdinal lastContiguousGid_;
199 LocalOrdinal numLocalElements_;
206 #endif // TPETRA_DETAILS_LOCALMAP_HPP KOKKOS_INLINE_FUNCTION LocalOrdinal getMinLocalIndex() const
The minimum local index.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
KOKKOS_INLINE_FUNCTION GlobalOrdinal getIndexBase() const
The (global) index base.
"Local" part of Map suitable for Kokkos kernels.
Implementation details of Tpetra.
KOKKOS_INLINE_FUNCTION LocalOrdinal getMaxLocalIndex() const
The maximum local index.
KOKKOS_INLINE_FUNCTION ValueType get(const KeyType &key) const
Get the value corresponding to the given key.
KOKKOS_INLINE_FUNCTION GlobalOrdinal getMinGlobalIndex() const
The minimum global index on the calling process.
KOKKOS_INLINE_FUNCTION bool isContiguous() const
Whether the Map is (locally) contiguous.
KOKKOS_INLINE_FUNCTION LocalOrdinal getLocalElement(const GlobalOrdinal globalIndex) const
Get the local index corresponding to the given global index.
KOKKOS_INLINE_FUNCTION GlobalOrdinal getGlobalElement(const LocalOrdinal localIndex) const
Get the global index corresponding to the given local index.
KOKKOS_INLINE_FUNCTION GlobalOrdinal getMaxGlobalIndex() const
The maximum global index on the calling process.
KOKKOS_INLINE_FUNCTION LocalOrdinal getNodeNumElements() const
The number of indices that live on the calling process.