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)
50 : map_(Teuchos::rcp(new Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>(numGlobalElements,
51 elementList,
52 indexBase,
53 comm))) {}
54
56template <class LocalOrdinal, class GlobalOrdinal, class Node>
58 TpetraMap(global_size_t numGlobalElements,
59 const Kokkos::View<const GlobalOrdinal *, typename Node::device_type> &indexList,
60 GlobalOrdinal indexBase,
61 const Teuchos::RCP<const Teuchos::Comm<int> > &comm)
62 : map_(Teuchos::rcp(new Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>(numGlobalElements,
63 indexList,
64 indexBase,
65 comm))) {}
66
68template <class LocalOrdinal, class GlobalOrdinal, class Node>
70
72
73template <class LocalOrdinal, class GlobalOrdinal, class Node>
75 XPETRA_MONITOR("TpetraMap::getGlobalNumElements");
76 return map_->getGlobalNumElements();
77}
78
79template <class LocalOrdinal, class GlobalOrdinal, class Node>
81 XPETRA_MONITOR("TpetraMap::getLocalNumElements");
82 return map_->getLocalNumElements();
83}
84
85template <class LocalOrdinal, class GlobalOrdinal, class Node>
87 XPETRA_MONITOR("TpetraMap::getIndexBase");
88 return map_->getIndexBase();
89}
90
91template <class LocalOrdinal, class GlobalOrdinal, class Node>
93 XPETRA_MONITOR("TpetraMap::getMinLocalIndex");
94 return map_->getMinLocalIndex();
95}
96
97template <class LocalOrdinal, class GlobalOrdinal, class Node>
99 XPETRA_MONITOR("TpetraMap::getMaxLocalIndex");
100 return map_->getMaxLocalIndex();
101}
102
103template <class LocalOrdinal, class GlobalOrdinal, class Node>
105 XPETRA_MONITOR("TpetraMap::getMinGlobalIndex");
106 return map_->getMinGlobalIndex();
107}
108
109template <class LocalOrdinal, class GlobalOrdinal, class Node>
111 XPETRA_MONITOR("TpetraMap::getMaxGlobalIndex");
112 return map_->getMaxGlobalIndex();
113}
114
115template <class LocalOrdinal, class GlobalOrdinal, class Node>
117 XPETRA_MONITOR("TpetraMap::getMinAllGlobalIndex");
118 return map_->getMinAllGlobalIndex();
119}
120
121template <class LocalOrdinal, class GlobalOrdinal, class Node>
123 XPETRA_MONITOR("TpetraMap::getMaxAllGlobalIndex");
124 return map_->getMaxAllGlobalIndex();
125}
126
127template <class LocalOrdinal, class GlobalOrdinal, class Node>
128LocalOrdinal TpetraMap<LocalOrdinal, GlobalOrdinal, Node>::getLocalElement(GlobalOrdinal globalIndex) const {
129 XPETRA_MONITOR("TpetraMap::getLocalElement");
130 return map_->getLocalElement(globalIndex);
131}
132
133template <class LocalOrdinal, class GlobalOrdinal, class Node>
134GlobalOrdinal TpetraMap<LocalOrdinal, GlobalOrdinal, Node>::getGlobalElement(LocalOrdinal localIndex) const {
135 XPETRA_MONITOR("TpetraMap::getGlobalElement");
136 return map_->getGlobalElement(localIndex);
137}
138
139template <class LocalOrdinal, class GlobalOrdinal, class Node>
141 XPETRA_MONITOR("TpetraMap::getRemoteIndexList");
142 return toXpetra(map_->getRemoteIndexList(GIDList, nodeIDList, LIDList));
143}
144
145template <class LocalOrdinal, class GlobalOrdinal, class Node>
147 XPETRA_MONITOR("TpetraMap::getRemoteIndexList");
148 return toXpetra(map_->getRemoteIndexList(GIDList, nodeIDList));
149}
150
151template <class LocalOrdinal, class GlobalOrdinal, class Node>
153 XPETRA_MONITOR("TpetraMap::getLocalElementList");
154 return map_->getLocalElementList();
155}
156
157template <class LocalOrdinal, class GlobalOrdinal, class Node>
159 XPETRA_MONITOR("TpetraMap::getMyGlobalIndicesDevice");
160 return map_->getMyGlobalIndicesDevice();
161}
162
163template <class LocalOrdinal, class GlobalOrdinal, class Node>
165 XPETRA_MONITOR("TpetraMap::isNodeLocalElement");
166 return map_->isNodeLocalElement(localIndex);
167}
168
169template <class LocalOrdinal, class GlobalOrdinal, class Node>
171 XPETRA_MONITOR("TpetraMap::isNodeGlobalElement");
172 return map_->isNodeGlobalElement(globalIndex);
173}
174
175template <class LocalOrdinal, class GlobalOrdinal, class Node>
177 XPETRA_MONITOR("TpetraMap::isContiguous");
178 return map_->isContiguous();
179}
180
181template <class LocalOrdinal, class GlobalOrdinal, class Node>
183 XPETRA_MONITOR("TpetraMap::isDistributed");
184 return map_->isDistributed();
185}
186
187template <class LocalOrdinal, class GlobalOrdinal, class Node>
189 XPETRA_MONITOR("TpetraMap::isCompatible");
190 return map_->isCompatible(toTpetra(map));
191}
192
193template <class LocalOrdinal, class GlobalOrdinal, class Node>
195 XPETRA_MONITOR("TpetraMap::isSameAs");
196 return map_->isSameAs(toTpetra(map));
197}
198
199template <class LocalOrdinal, class GlobalOrdinal, class Node>
201 XPETRA_MONITOR("TpetraMap::getComm");
202 return map_->getComm();
203}
204
205template <class LocalOrdinal, class GlobalOrdinal, class Node>
207 XPETRA_MONITOR("TpetraMap::description");
208 return map_->description();
209}
210
211template <class LocalOrdinal, class GlobalOrdinal, class Node>
213 XPETRA_MONITOR("TpetraMap::describe");
214 map_->describe(out, verbLevel);
215}
216
217template <class LocalOrdinal, class GlobalOrdinal, class Node>
221
222template <class LocalOrdinal, class GlobalOrdinal, class Node>
226
227template <class LocalOrdinal, class GlobalOrdinal, class Node>
230
231template <class LocalOrdinal, class GlobalOrdinal, class Node>
233
234template <class LocalOrdinal, class GlobalOrdinal, class Node>
236
237template <class LocalOrdinal, class GlobalOrdinal, class Node>
241
242} // namespace Xpetra
243
244// TODO: remove?
246template <class LocalOrdinal, class GlobalOrdinal, class Node>
248 XPETRA_MONITOR("TpetraMap==TpetraMap");
249 return map1.isSameAs(map2);
250}
251
253template <class LocalOrdinal, class GlobalOrdinal, class Node>
255 XPETRA_MONITOR("TpetraMap!=TpetraMap");
256 return !map1.isSameAs(map2);
257}
258
259#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)