10 #ifndef TPETRA_DETAILS_LOCALMAP_HPP
11 #define TPETRA_DETAILS_LOCALMAP_HPP
17 #include "Tpetra_Details_FixedHashTable.hpp"
38 template <
class LocalOrdinal,
class GlobalOrdinal,
class DeviceType>
57 KOKKOS_DEFAULTED_FUNCTION
LocalMap() =
default;
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)
70 , indexBase_(indexBase)
73 , firstContiguousGid_(firstContiguousGid)
74 , lastContiguousGid_(lastContiguousGid)
75 , numLocalElements_(numLocalElements)
76 , contiguous_(contiguous) {}
80 return numLocalElements_;
102 KOKKOS_INLINE_FUNCTION LocalOrdinal
104 if (numLocalElements_ == 0) {
105 return ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid();
107 return static_cast<LocalOrdinal
>(numLocalElements_ - 1);
122 KOKKOS_INLINE_FUNCTION LocalOrdinal
125 if (globalIndex < myMinGid_ || globalIndex > myMaxGid_) {
126 return ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid();
128 return static_cast<LocalOrdinal
>(globalIndex - myMinGid_);
129 }
else if (globalIndex >= firstContiguousGid_ &&
130 globalIndex <= lastContiguousGid_) {
131 return static_cast<LocalOrdinal
>(globalIndex - firstContiguousGid_);
135 return glMap_.
get(globalIndex);
140 KOKKOS_INLINE_FUNCTION GlobalOrdinal
143 return ::Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
148 return lgMap_(localIndex);
169 ::Kokkos::View<const GlobalOrdinal*, ::Kokkos::LayoutLeft, device_type> lgMap_;
171 GlobalOrdinal indexBase_ = 0;
172 GlobalOrdinal myMinGid_ = Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
173 GlobalOrdinal myMaxGid_ = Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
174 GlobalOrdinal firstContiguousGid_ = Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
175 GlobalOrdinal lastContiguousGid_ = Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
176 LocalOrdinal numLocalElements_ = 0;
177 bool contiguous_ =
false;
183 #endif // TPETRA_DETAILS_LOCALMAP_HPP
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.
"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. (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.
DT device_type
The device type.
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.