Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_TpetraMap_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef XPETRA_TPETRAMAP_DEF_HPP
11#define XPETRA_TPETRAMAP_DEF_HPP
12
14
15namespace Xpetra {
16
18
19
20template <class LocalOrdinal, class GlobalOrdinal, class Node>
22 TpetraMap(global_size_t numGlobalElements,
23 GlobalOrdinal indexBase,
24 const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
25 LocalGlobal lg)
26 : map_(Teuchos::rcp(new Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>(numGlobalElements,
27 indexBase, comm,
28 toTpetra(lg)))) {}
29
31
32template <class LocalOrdinal, class GlobalOrdinal, class Node>
34 TpetraMap(global_size_t numGlobalElements,
35 size_t numLocalElements,
36 GlobalOrdinal indexBase,
37 const Teuchos::RCP<const Teuchos::Comm<int> > &comm)
38 : map_(Teuchos::rcp(new Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>(numGlobalElements,
39 numLocalElements,
40 indexBase, comm))) {}
41
43
44template <class LocalOrdinal, class GlobalOrdinal, class Node>
46 TpetraMap(global_size_t numGlobalElements,
48 GlobalOrdinal indexBase,
49 const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
51 : map_(Teuchos::rcp(new Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>(numGlobalElements,
52 elementList,
53 indexBase,
54 comm,
55 params))) {}
56
58template <class LocalOrdinal, class GlobalOrdinal, class Node>
60 TpetraMap(global_size_t numGlobalElements,
61 const Kokkos::View<const GlobalOrdinal *, typename Node::device_type> &indexList,
62 GlobalOrdinal indexBase,
63 const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
65 : map_(Teuchos::rcp(new Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>(numGlobalElements,
66 indexList,
67 indexBase,
68 comm,
69 params))) {}
70
72template <class LocalOrdinal, class GlobalOrdinal, class Node>
74
76
77template <class LocalOrdinal, class GlobalOrdinal, class Node>
79 XPETRA_MONITOR("TpetraMap::getGlobalNumElements");
80 return map_->getGlobalNumElements();
81}
82
83template <class LocalOrdinal, class GlobalOrdinal, class Node>
85 XPETRA_MONITOR("TpetraMap::getLocalNumElements");
86 return map_->getLocalNumElements();
87}
88
89template <class LocalOrdinal, class GlobalOrdinal, class Node>
91 XPETRA_MONITOR("TpetraMap::getIndexBase");
92 return map_->getIndexBase();
93}
94
95template <class LocalOrdinal, class GlobalOrdinal, class Node>
97 XPETRA_MONITOR("TpetraMap::getMinLocalIndex");
98 return map_->getMinLocalIndex();
99}
100
101template <class LocalOrdinal, class GlobalOrdinal, class Node>
103 XPETRA_MONITOR("TpetraMap::getMaxLocalIndex");
104 return map_->getMaxLocalIndex();
105}
106
107template <class LocalOrdinal, class GlobalOrdinal, class Node>
109 XPETRA_MONITOR("TpetraMap::getMinGlobalIndex");
110 return map_->getMinGlobalIndex();
111}
112
113template <class LocalOrdinal, class GlobalOrdinal, class Node>
115 XPETRA_MONITOR("TpetraMap::getMaxGlobalIndex");
116 return map_->getMaxGlobalIndex();
117}
118
119template <class LocalOrdinal, class GlobalOrdinal, class Node>
121 XPETRA_MONITOR("TpetraMap::getMinAllGlobalIndex");
122 return map_->getMinAllGlobalIndex();
123}
124
125template <class LocalOrdinal, class GlobalOrdinal, class Node>
127 XPETRA_MONITOR("TpetraMap::getMaxAllGlobalIndex");
128 return map_->getMaxAllGlobalIndex();
129}
130
131template <class LocalOrdinal, class GlobalOrdinal, class Node>
132LocalOrdinal TpetraMap<LocalOrdinal, GlobalOrdinal, Node>::getLocalElement(GlobalOrdinal globalIndex) const {
133 XPETRA_MONITOR("TpetraMap::getLocalElement");
134 return map_->getLocalElement(globalIndex);
135}
136
137template <class LocalOrdinal, class GlobalOrdinal, class Node>
138GlobalOrdinal TpetraMap<LocalOrdinal, GlobalOrdinal, Node>::getGlobalElement(LocalOrdinal localIndex) const {
139 XPETRA_MONITOR("TpetraMap::getGlobalElement");
140 return map_->getGlobalElement(localIndex);
141}
142
143template <class LocalOrdinal, class GlobalOrdinal, class Node>
145 XPETRA_MONITOR("TpetraMap::getRemoteIndexList");
146 return toXpetra(map_->getRemoteIndexList(GIDList, nodeIDList, LIDList));
147}
148
149template <class LocalOrdinal, class GlobalOrdinal, class Node>
151 XPETRA_MONITOR("TpetraMap::getRemoteIndexList");
152 return toXpetra(map_->getRemoteIndexList(GIDList, nodeIDList));
153}
154
155template <class LocalOrdinal, class GlobalOrdinal, class Node>
157 XPETRA_MONITOR("TpetraMap::getLocalElementList");
158 return map_->getLocalElementList();
159}
160
161template <class LocalOrdinal, class GlobalOrdinal, class Node>
163 XPETRA_MONITOR("TpetraMap::getMyGlobalIndicesDevice");
164 return map_->getMyGlobalIndicesDevice();
165}
166
167template <class LocalOrdinal, class GlobalOrdinal, class Node>
169 XPETRA_MONITOR("TpetraMap::isNodeLocalElement");
170 return map_->isNodeLocalElement(localIndex);
171}
172
173template <class LocalOrdinal, class GlobalOrdinal, class Node>
175 XPETRA_MONITOR("TpetraMap::isNodeGlobalElement");
176 return map_->isNodeGlobalElement(globalIndex);
177}
178
179template <class LocalOrdinal, class GlobalOrdinal, class Node>
181 XPETRA_MONITOR("TpetraMap::isContiguous");
182 return map_->isContiguous();
183}
184
185template <class LocalOrdinal, class GlobalOrdinal, class Node>
187 XPETRA_MONITOR("TpetraMap::isDistributed");
188 return map_->isDistributed();
189}
190
191template <class LocalOrdinal, class GlobalOrdinal, class Node>
193 XPETRA_MONITOR("TpetraMap::isCompatible");
194 return map_->isCompatible(toTpetra(map));
195}
196
197template <class LocalOrdinal, class GlobalOrdinal, class Node>
199 XPETRA_MONITOR("TpetraMap::isSameAs");
200 return map_->isSameAs(toTpetra(map));
201}
202
203template <class LocalOrdinal, class GlobalOrdinal, class Node>
205 XPETRA_MONITOR("TpetraMap::getComm");
206 return map_->getComm();
207}
208
209template <class LocalOrdinal, class GlobalOrdinal, class Node>
211 XPETRA_MONITOR("TpetraMap::description");
212 return map_->description();
213}
214
215template <class LocalOrdinal, class GlobalOrdinal, class Node>
217 XPETRA_MONITOR("TpetraMap::describe");
218 map_->describe(out, verbLevel);
219}
220
221template <class LocalOrdinal, class GlobalOrdinal, class Node>
225
226template <class LocalOrdinal, class GlobalOrdinal, class Node>
230
231template <class LocalOrdinal, class GlobalOrdinal, class Node>
234
235template <class LocalOrdinal, class GlobalOrdinal, class Node>
237
238template <class LocalOrdinal, class GlobalOrdinal, class Node>
240
241template <class LocalOrdinal, class GlobalOrdinal, class Node>
245
246} // namespace Xpetra
247
248// TODO: remove?
250template <class LocalOrdinal, class GlobalOrdinal, class Node>
252 XPETRA_MONITOR("TpetraMap==TpetraMap");
253 return map1.isSameAs(map2);
254}
255
257template <class LocalOrdinal, class GlobalOrdinal, class Node>
259 XPETRA_MONITOR("TpetraMap!=TpetraMap");
260 return !map1.isSameAs(map2);
261}
262
263#endif // XPETRA_TPETRAMAP_DEF_HPP
#define XPETRA_MONITOR(funcName)
bool operator!=(const Xpetra::TpetraMap< LocalOrdinal, GlobalOrdinal, Node > &map1, const Xpetra::TpetraMap< LocalOrdinal, GlobalOrdinal, Node > &map2)
Returns true if map is not identical to this map. Implemented in TpetraMap::isSameAs().
bool operator==(const Xpetra::TpetraMap< LocalOrdinal, GlobalOrdinal, Node > &map1, const Xpetra::TpetraMap< LocalOrdinal, GlobalOrdinal, Node > &map2)
Returns true if map is identical to this map. Implemented in TpetraMap::isSameAs().
Kokkos::View< const global_ordinal_type *, typename Node::device_type > global_indices_array_device_type
bool isDistributed() const
Whether this Map is globally distributed or locally replicated.
UnderlyingLib lib() const
Get the library used by this object (Tpetra or Epetra?)
LocalOrdinal getLocalElement(GlobalOrdinal globalIndex) const
The local index corresponding to the given global index.
LookupStatus getRemoteIndexList(const Teuchos::ArrayView< const GlobalOrdinal > &GIDList, const Teuchos::ArrayView< int > &nodeIDList, const Teuchos::ArrayView< LocalOrdinal > &LIDList) const
Return the process IDs and corresponding local IDs for the given global IDs.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print this object with the given verbosity level to the given FancyOStream.
GlobalOrdinal getMaxAllGlobalIndex() const
The maximum global index over all processes in the communicator.
std::string description() const
Return a simple one-line description of this object.
global_size_t getGlobalNumElements() const
The number of elements in this Map.
Teuchos::ArrayView< const GlobalOrdinal > getLocalElementList() const
Return a view of the global indices owned by this node.
size_t getLocalNumElements() const
The number of elements belonging to the calling node.
bool isContiguous() const
True if this Map is distributed contiguously, else false.
TpetraMap(global_size_t numGlobalElements, GlobalOrdinal indexBase, const Teuchos::RCP< const Teuchos::Comm< int > > &comm, LocalGlobal lg=GloballyDistributed)
Constructor with Tpetra-defined contiguous uniform distribution.
bool isSameAs(const Map< LocalOrdinal, GlobalOrdinal, Node > &map) const
True if and only if map is identical to this Map.
RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_Map() const
Get the underlying Tpetra map.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > removeEmptyProcesses() const
Return a new Map with processes with zero elements removed.
GlobalOrdinal getIndexBase() const
The index base for this Map.
GlobalOrdinal getMinGlobalIndex() const
The minimum global index owned by the calling process.
GlobalOrdinal getMinAllGlobalIndex() const
The minimum global index over all processes in the communicator.
local_map_type getLocalMap() const
Get the local Map for Kokkos kernels.
bool isNodeLocalElement(LocalOrdinal localIndex) const
True if the local index is valid for this Map on this node, else false.
LocalOrdinal getMaxLocalIndex() const
The maximum local index on the calling process.
RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > replaceCommWithSubset(const Teuchos::RCP< const Teuchos::Comm< int > > &newComm) const
Replace this Map's communicator with a subset communicator.
bool isCompatible(const Map< LocalOrdinal, GlobalOrdinal, Node > &map) const
True if and only if map is compatible with this Map.
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Get this Map's Comm object.
GlobalOrdinal getGlobalElement(LocalOrdinal localIndex) const
The global index corresponding to the given local index.
GlobalOrdinal getMaxGlobalIndex() const
The maximum global index owned by the calling process.
LocalOrdinal getMinLocalIndex() const
The minimum local index.
global_indices_array_device_type getMyGlobalIndicesDevice() const
Return a view of the global indices owned by this process.
bool isNodeGlobalElement(GlobalOrdinal globalIndex) const
True if the global index is found in this Map on this node, else false.
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > graph)
size_t global_size_t
Global size_t object.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)