Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
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
22namespace Tpetra {
23namespace Details {
24
38template <class LocalOrdinal, class GlobalOrdinal, class DeviceType>
39class LocalMap {
40 public:
43
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
58
59 LocalMap(const ::Tpetra::Details::FixedHashTable<GlobalOrdinal, LocalOrdinal, device_type>& glMap,
60 const ::Kokkos::View<const GlobalOrdinal*, ::Kokkos::LayoutLeft, device_type>& lgMap,
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
80 return numLocalElements_;
81 }
82
85 return indexBase_;
86 }
87
93 return contiguous_;
94 }
95
100
104 if (numLocalElements_ == 0) {
105 return ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid();
106 } else { // Local indices are always zero-based.
107 return static_cast<LocalOrdinal>(numLocalElements_ - 1);
108 }
109 }
110
113 return myMinGid_;
114 }
115
118 return myMaxGid_;
119 }
120
124 if (contiguous_) {
126 return ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid();
127 }
128 return static_cast<LocalOrdinal>(globalIndex - myMinGid_);
129 } else if (globalIndex >= firstContiguousGid_ &&
130 globalIndex <= lastContiguousGid_) {
131 return static_cast<LocalOrdinal>(globalIndex - firstContiguousGid_);
132 } else {
133 // If the given global index is not in the table, this returns
134 // the same value as OrdinalTraits<LocalOrdinal>::invalid().
135 return glMap_.get(globalIndex);
136 }
137 }
138
143 return ::Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
144 }
145 if (isContiguous()) {
146 return getMinGlobalIndex() + localIndex;
147 } else {
148 return lgMap_(localIndex);
149 }
150 }
151
152 private:
169 ::Kokkos::View<const GlobalOrdinal*, ::Kokkos::LayoutLeft, device_type> lgMap_;
170
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;
178};
179
180} // namespace Details
181} // namespace Tpetra
182
183#endif // TPETRA_DETAILS_LOCALMAP_HPP
Forward declaration of Tpetra::Details::LocalMap.
Struct that holds views of the contents of a CrsMatrix.
KOKKOS_INLINE_FUNCTION ValueType get(const KeyType &key) const
Get the value corresponding to the given key.
"Local" part of Map suitable for Kokkos kernels.
KOKKOS_INLINE_FUNCTION LocalOrdinal getLocalNumElements() const
The number of indices that live on the calling process.
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.
KOKKOS_INLINE_FUNCTION GlobalOrdinal getMaxGlobalIndex() const
The maximum global index on the calling process.
KOKKOS_INLINE_FUNCTION bool isContiguous() const
Whether the Map is (locally) contiguous.
KOKKOS_INLINE_FUNCTION LocalOrdinal getMinLocalIndex() const
The minimum local index.
KOKKOS_DEFAULTED_FUNCTION LocalMap()=default
Default constructor.
typename device_type::memory_space memory_space
The Kokkos memory space.
KOKKOS_INLINE_FUNCTION LocalOrdinal getMaxLocalIndex() const
The maximum local index.
KOKKOS_INLINE_FUNCTION GlobalOrdinal getIndexBase() const
The (global) index base.
KOKKOS_INLINE_FUNCTION GlobalOrdinal getGlobalElement(const LocalOrdinal localIndex) const
Get the global index corresponding to the given local index. (device only)
KOKKOS_INLINE_FUNCTION GlobalOrdinal getMinGlobalIndex() const
The minimum global index on the calling process.
DeviceType device_type
The device type.
Implementation details of Tpetra.
Namespace Tpetra contains the class and methods constituting the Tpetra library.