Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_BlockedMap_decl.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 PACKAGES_XPETRA_SUP_BLOCKEDMAP_XPETRA_BLOCKEDMAP_DECL_HPP_
11#define PACKAGES_XPETRA_SUP_BLOCKEDMAP_XPETRA_BLOCKEDMAP_DECL_HPP_
12
13#include "Xpetra_ConfigDefs.hpp"
14
15#include "Xpetra_Import.hpp"
16#include "Xpetra_Map_decl.hpp"
17
18namespace Xpetra {
19
20template <class LocalOrdinal,
21 class GlobalOrdinal,
22 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
23class BlockedMap : public Map<LocalOrdinal, GlobalOrdinal, Node> {
24 public:
25 typedef LocalOrdinal local_ordinal_type;
26 typedef GlobalOrdinal global_ordinal_type;
27 typedef Node node_type;
29
30 private:
31#undef XPETRA_BLOCKEDMAP_SHORT
33
34 public:
36
37
39
41 BlockedMap();
42
60 BlockedMap(const RCP<const Map>& fullmap, const std::vector<RCP<const Map>>& maps, bool bThyraMode = false);
61
63 BlockedMap(const std::vector<RCP<const Map>>& maps, const std::vector<RCP<const Map>>& thyramaps);
64
66 BlockedMap(const BlockedMap& input);
67
69 virtual ~BlockedMap();
70
72
73
75 virtual global_size_t getGlobalNumElements() const;
76
78 virtual size_t getLocalNumElements() const;
79
81 virtual GlobalOrdinal getIndexBase() const;
82
84 virtual LocalOrdinal getMinLocalIndex() const;
85
87 virtual LocalOrdinal getMaxLocalIndex() const;
88
90 virtual GlobalOrdinal getMinGlobalIndex() const;
91
93 virtual GlobalOrdinal getMaxGlobalIndex() const;
94
96 virtual GlobalOrdinal getMinAllGlobalIndex() const;
97
99 virtual GlobalOrdinal getMaxAllGlobalIndex() const;
100
102 virtual LocalOrdinal getLocalElement(GlobalOrdinal globalIndex) const;
103
105 virtual GlobalOrdinal getGlobalElement(LocalOrdinal localIndex) const;
106
109 const Teuchos::ArrayView<int>& /* nodeIDList */,
110 const Teuchos::ArrayView<LocalOrdinal>& /* LIDList */) const;
111
114 const Teuchos::ArrayView<int>& /* nodeIDList */) const;
115
118
121
123
125
126
128 virtual bool isNodeLocalElement(LocalOrdinal localIndex) const;
129
131 virtual bool isNodeGlobalElement(GlobalOrdinal globalIndex) const;
132
134 virtual bool isContiguous() const;
135
137 virtual bool isDistributed() const;
138
141
143 virtual bool isSameAs(const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>& map) const;
144
146
148
149
152
154
163
165
167
168
169 /*virtual size_t getLocalLength() const {
170 throw Xpetra::Exceptions::RuntimeError("BlockedMap::getLocalLength: routine not implemented.");
171 return 0;
172 }*/
173
175 /*virtual global_size_t getGlobalLength() const {
176 throw Xpetra::Exceptions::RuntimeError("BlockedMap::getGlobalLength: routine not implemented.");
177 return 0;
178 }*/
179
181 virtual bool getThyraMode() const;
182
184
186
187
190
193 replaceCommWithSubset(const Teuchos::RCP<const Teuchos::Comm<int>>& /* newComm */) const;
194
196
198
199
201 virtual UnderlyingLib lib() const;
202
203 // TODO: find a better solution for this hack
204 // The problem is that EpetraMap, TpetraMap and StridedMap all inherit Map. To have proper toEpetra() we
205 // need to understand the type of underlying matrix. But in src/Map we have no knowledge of StridedMaps, so
206 // we cannot check for it by casting. This function allows us to avoid the restriction, as StridedMap redefines
207 // it to return the base map.
209 getMap() const;
210
212
214 size_t getNumMaps() const;
215
221 getMap(size_t i, bool bThyraMode = false) const;
222
225 getImporter(size_t i) const;
226
229 getFullMap() const;
230
232 size_t getMapIndexForGID(GlobalOrdinal gid) const;
233
234#ifdef HAVE_XPETRA_TPETRA
236
238 local_map_type getLocalMap() const { return fullmap_->getLocalMap(); }
239
240#else // HAVE_XPETRA_TPETRA
241#ifdef __GNUC__
242#warning "Xpetra Kokkos interface for CrsMatrix is enabled (HAVE_XPETRA_KOKKOS_REFACTOR) but Tpetra is disabled. The Kokkos interface needs Tpetra to be enabled, too."
243#endif // __GNUC__
244#endif // #else !HAVE_XPETRA_TPETRA
245
247
249
250
252 virtual std::string description() const;
253
256
258
259 protected:
262 virtual void assign(const BlockedMap& input);
263
280
281 private:
282 bool CheckConsistency() const;
283
284 private:
286 std::vector<RCP<const Map>> maps_;
287 std::vector<RCP<Import>> importers_;
288 bool bThyraMode_; //< boolean flag: use Thyra numbering for local sub-block maps. default = false (for Xpetra mode)
289 std::vector<RCP<const Map>> thyraMaps_; //< store Thyra-style numbering maps here in Thyra mode. In Xpetra mode this vector is empty.
290}; // BlockedMap class
291
292} // namespace Xpetra
293
294#define XPETRA_BLOCKEDMAP_SHORT
295
296#endif /* PACKAGES_XPETRA_SUP_BLOCKEDMAP_XPETRA_BLOCKEDMAP_DECL_HPP_ */
static const EVerbosityLevel verbLevel_default
virtual GlobalOrdinal getMinAllGlobalIndex() const
The minimum global index over all processes in the communicator.
virtual global_size_t getGlobalNumElements() const
The number of elements in this Map.
virtual GlobalOrdinal getMaxGlobalIndex() const
The maximum global index owned by the calling process.
BlockedMap< LocalOrdinal, GlobalOrdinal, Node > & operator=(const BlockedMap &rhs)
Assignment operator: Does a deep copy.
virtual GlobalOrdinal getMaxAllGlobalIndex() const
The maximum global index over all processes in the communicator.
virtual GlobalOrdinal getIndexBase() const
The index base for this Map.
virtual bool isNodeGlobalElement(GlobalOrdinal globalIndex) const
Whether the given global index is valid for this Map on this process.
virtual RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
const RCP< Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > getImporter(size_t i) const
get the importer between full map and partial map
static Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > concatenateMaps(const std::vector< Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > > &subMaps)
Helper function to concatenate several maps.
virtual Teuchos::ArrayView< const GlobalOrdinal > getLocalElementList() const
Return a view of the global indices owned by this process.
virtual GlobalOrdinal getGlobalElement(LocalOrdinal localIndex) const
The global index corresponding to the given local index.
virtual bool isDistributed() const
Whether this Map is globally distributed or locally replicated.
typename Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node >::local_map_type local_map_type
size_t getMapIndexForGID(GlobalOrdinal gid) const
returns map index in map extractor which contains GID
virtual bool isSameAs(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map) const
True if and only if map is identical to this Map.
virtual bool isCompatible(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map) const
True if and only if map is compatible with this Map.
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getFullMap() const
the full map
virtual global_indices_array_device_type getMyGlobalIndicesDevice() const
Return a view of the global indices owned by this process.
std::vector< RCP< const Map > > thyraMaps_
size_t getNumMaps() const
number of partial maps
virtual GlobalOrdinal getMinGlobalIndex() const
The minimum global index owned by the calling process.
virtual LocalOrdinal getMaxLocalIndex() const
The maximum local index on the calling process.
std::vector< RCP< const Map > > maps_
virtual bool isNodeLocalElement(LocalOrdinal localIndex) const
Whether the given local index is valid for this Map on this process.
Map< LocalOrdinal, GlobalOrdinal, Node >::global_indices_array_device_type global_indices_array_device_type
std::vector< RCP< Import > > importers_
local_map_type getLocalMap() const
Get the local Map for Kokkos kernels.
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Get this Map's Comm object.
virtual ~BlockedMap()
Destructor.
virtual bool getThyraMode() const
Local number of rows on the calling process.
virtual RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > replaceCommWithSubset(const Teuchos::RCP< const Teuchos::Comm< int > > &) const
Replace this Map's communicator with a subset communicator.
virtual bool isContiguous() const
True if this Map is distributed contiguously, else false.
virtual UnderlyingLib lib() const
Get the library used by this object (Tpetra or Epetra?)
virtual LocalOrdinal getMinLocalIndex() const
The minimum local index.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
virtual LocalOrdinal getLocalElement(GlobalOrdinal globalIndex) const
The local index corresponding to the given global index.
virtual std::string description() const
A simple one-line description of this object.
virtual RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > removeEmptyProcesses() const
Return a new Map with processes with zero elements removed.
virtual size_t getLocalNumElements() const
The number of elements belonging to the calling process.
virtual void assign(const BlockedMap &input)
Implementation of the assignment operator (operator=); does a deep copy.
virtual LookupStatus getRemoteIndexList(const Teuchos::ArrayView< const GlobalOrdinal > &, const Teuchos::ArrayView< int > &, const Teuchos::ArrayView< LocalOrdinal > &) const
Return the process ranks and corresponding local indices for the given global indices.
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
size_t global_size_t
Global size_t object.