12#ifndef XPETRA_STRIDEDMAP_DECL_HPP
13#define XPETRA_STRIDEDMAP_DECL_HPP
15#include <Tpetra_KokkosCompat_DefaultNode.hpp>
16#include <Teuchos_Describable.hpp>
61template <
class LocalOrdinal,
63 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
64class StridedMap :
public virtual Map<LocalOrdinal, GlobalOrdinal, Node> {
72#undef XPETRA_STRIDEDMAP_SHORT
101 GlobalOrdinal indexBase,
102 std::vector<size_t>& stridingInfo,
104 LocalOrdinal stridedBlockId = -1,
105 GlobalOrdinal offset = 0,
132 size_t numLocalElements,
133 GlobalOrdinal indexBase,
134 std::vector<size_t>& stridingInfo,
136 LocalOrdinal stridedBlockId = -1,
137 GlobalOrdinal offset = 0);
153 GlobalOrdinal indexBase,
154 std::vector<size_t>& stridingInfo,
156 LocalOrdinal stridedBlockId = -1);
159 std::vector<size_t>& stridingInfo,
161 LocalOrdinal stridedBlockId = -1,
162 GlobalOrdinal offset = 0);
201#ifdef HAVE_XPETRA_TPETRA
206 return map_->getLocalMap();
211 "Xpetra Kokkos interface for CrsMatrix is enabled (HAVE_XPETRA_KOKKOS_REFACTOR) but Tpetra is disabled. The Kokkos interface needs Tpetra to be enabled, too."
366#define XPETRA_STRIDEDMAP_SHORT
static const EVerbosityLevel verbLevel_default
Kokkos::View< const global_ordinal_type *, typename Node::device_type > global_indices_array_device_type
Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node >::local_map_type local_map_type
Class that stores a strided map.
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
LocalOrdinal getMaxLocalIndex() const
Returns maximum local index.
std::vector< size_t > stridingInfo_
Vector with size of strided blocks (dofs)
GlobalOrdinal getMinAllGlobalIndex() const
Return the minimum global index over all nodes.
LocalOrdinal local_ordinal_type
virtual ~StridedMap()
Destructor.
global_indices_array_device_type getMyGlobalIndicesDevice() const
Return a view of the global indices owned by this process on the Map's device.
bool isSameAs(const Map &map) const
Returns true if map is identical to this Map.
bool isNodeGlobalElement(GlobalOrdinal globalIndex) const
Returns true if the global index is found in this Map on this node; returns false if it isn't.
size_t GID2StridingBlockId(GlobalOrdinal gid) const
GlobalOrdinal indexBase_
Index base for the strided map (default = 0)
Map< LocalOrdinal, GlobalOrdinal, Node >::global_indices_array_device_type global_indices_array_device_type
bool isNodeLocalElement(LocalOrdinal localIndex) const
Returns true if the local index is valid for this Map on this node; returns false if it isn't.
LocalOrdinal stridedBlockId_
Member variable denoting which dofs are stored in map.
GlobalOrdinal getGlobalElement(LocalOrdinal localIndex) const
Return the global index for a given local index.
UnderlyingLib lib() const
Get the library used by this object (Tpetra or Epetra?)
global_size_t getGlobalNumElements() const
Returns the number of elements in this Map.
bool isDistributed() const
Returns true if this Map is distributed across more than one node; returns false otherwise.
LocalOrdinal getStridedBlockId() const
void setStridingData(std::vector< size_t > stridingInfo)
LookupStatus getRemoteIndexList(const Teuchos::ArrayView< const GlobalOrdinal > &GIDList, const Teuchos::ArrayView< int > &nodeIDList, const Teuchos::ArrayView< LocalOrdinal > &LIDList) const
Returns the node IDs and corresponding local indices for a given list of global indices.
size_t getLocalNumElements() const
Returns the number of elements belonging to the calling node.
GlobalOrdinal getOffset() const
bool isContiguous() const
Returns true if this Map is distributed contiguously; returns false otherwise.
GlobalOrdinal global_ordinal_type
typename Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node >::local_map_type local_map_type
GlobalOrdinal getMaxAllGlobalIndex() const
Return the maximum global index over all nodes.
RCP< const Map > replaceCommWithSubset(const Teuchos::RCP< const Teuchos::Comm< int > > &newComm) const
Replace this Map's communicator with a subset communicator.
void setOffset(GlobalOrdinal offset)
GlobalOrdinal offset_
Offset for gids in map (default = 0)
GlobalOrdinal getMinGlobalIndex() const
Returns minimum global index owned by this node.
std::vector< size_t > getStridingData() const
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Get the Comm object for this Map.
bool isCompatible(const Map &map) const
Returns true if map is compatible with this Map.
size_t getFixedBlockSize() const
LocalOrdinal getLocalElement(GlobalOrdinal globalIndex) const
Return the local index for a given global index.
bool isStrided() const
returns true, if this is a strided map (i.e. more than 1 strided blocks)
virtual bool CheckConsistency()
Teuchos::ArrayView< const GlobalOrdinal > getLocalElementList() const
Return a list of the global indices owned by this node.
std::string description() const
Return a simple one-line description of this object.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to a FancyOStream object.
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > map_
LocalOrdinal getMinLocalIndex() const
Returns minimum local index.
local_map_type getLocalMap() const
Get the local Map for Kokkos kernels.
GlobalOrdinal getMaxGlobalIndex() const
Returns maximum global index owned by this node.
RCP< const Map > removeEmptyProcesses() const
Return a new Map with processes with zero elements removed.
GlobalOrdinal getIndexBase() const
Returns the index base for this Map.
size_t global_size_t
Global size_t object.