10#ifndef XPETRA_TPETRAMAP_DEF_HPP
11#define XPETRA_TPETRAMAP_DEF_HPP
20template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
23 GlobalOrdinal indexBase,
26 : map_(
Teuchos::rcp(new
Tpetra::
Map<LocalOrdinal, GlobalOrdinal, Node>(numGlobalElements,
32template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
35 size_t numLocalElements,
36 GlobalOrdinal indexBase,
38 : map_(
Teuchos::rcp(new
Tpetra::
Map<LocalOrdinal, GlobalOrdinal, Node>(numGlobalElements,
44template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
48 GlobalOrdinal indexBase,
51 : map_(
Teuchos::rcp(new
Tpetra::
Map<LocalOrdinal, GlobalOrdinal, Node>(numGlobalElements,
58template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
61 const Kokkos::View<const GlobalOrdinal *, typename Node::device_type> &indexList,
62 GlobalOrdinal indexBase,
65 : map_(
Teuchos::rcp(new
Tpetra::
Map<LocalOrdinal, GlobalOrdinal, Node>(numGlobalElements,
72template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
77template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
80 return map_->getGlobalNumElements();
83template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
86 return map_->getLocalNumElements();
89template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
92 return map_->getIndexBase();
95template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
98 return map_->getMinLocalIndex();
101template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
104 return map_->getMaxLocalIndex();
107template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
110 return map_->getMinGlobalIndex();
113template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
116 return map_->getMaxGlobalIndex();
119template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
122 return map_->getMinAllGlobalIndex();
125template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
128 return map_->getMaxAllGlobalIndex();
131template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
134 return map_->getLocalElement(globalIndex);
137template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
140 return map_->getGlobalElement(localIndex);
143template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
146 return toXpetra(map_->getRemoteIndexList(GIDList, nodeIDList, LIDList));
149template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
152 return toXpetra(map_->getRemoteIndexList(GIDList, nodeIDList));
155template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
158 return map_->getLocalElementList();
161template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
164 return map_->getMyGlobalIndicesDevice();
167template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
170 return map_->isNodeLocalElement(localIndex);
173template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
176 return map_->isNodeGlobalElement(globalIndex);
179template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
182 return map_->isContiguous();
185template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
188 return map_->isDistributed();
191template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
194 return map_->isCompatible(
toTpetra(map));
197template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
200 return map_->isSameAs(
toTpetra(map));
203template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
206 return map_->getComm();
209template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
212 return map_->description();
215template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
218 map_->describe(out, verbLevel);
221template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
223 return toXpetra(map_->removeEmptyProcesses());
226template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
228 return toXpetra(map_->replaceCommWithSubset(newComm));
231template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
235template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
238template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
241template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
243 return map_->getLocalMap();
250template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
257template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
#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)