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 {
40public:
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 {}
78
81 return numLocalElements_;
82 }
83
86 return indexBase_;
87 }
88
94 return contiguous_;
95 }
96
101
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
115 return myMinGid_;
116 }
117
120 return myMaxGid_;
121 }
122
126 {
127 if (contiguous_) {
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
147 {
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
159private:
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
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.