Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Map_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef TPETRA_MAP_DECL_HPP
11 #define TPETRA_MAP_DECL_HPP
12 
16 
17 #include "Tpetra_ConfigDefs.hpp"
18 #include "Tpetra_Map_fwd.hpp"
19 #include "Tpetra_Directory_fwd.hpp"
20 #include "Tpetra_TieBreak_fwd.hpp"
22 #include "Tpetra_KokkosCompat_DefaultNode.hpp"
23 #include "Kokkos_DualView.hpp"
24 #include "Teuchos_Array.hpp"
25 #include "Teuchos_Comm.hpp"
26 #include "Teuchos_Describable.hpp"
27 
28 
29 namespace Tpetra {
30 
195  template <class LocalOrdinal,
196  class GlobalOrdinal,
197  class Node>
198  class Map : public Teuchos::Describable {
199  public:
201 
202 
204  using local_ordinal_type = LocalOrdinal;
205 
207  using global_ordinal_type = GlobalOrdinal;
208 
214  using device_type = typename Node::device_type;
215 
217  using execution_space = typename device_type::execution_space;
218 
220  using memory_space = typename device_type::memory_space;
221 
223  using node_type = Node;
224 
238  using local_map_type =
242 
244 
246 
296  Map (const global_size_t numGlobalElements,
297  const global_ordinal_type indexBase,
298  const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
299  const LocalGlobal lg=GloballyDistributed);
300 
301 
337  Map (const global_size_t numGlobalElements,
338  const size_t numLocalElements,
339  const global_ordinal_type indexBase,
340  const Teuchos::RCP<const Teuchos::Comm<int> > &comm);
341 
342 
383  Map (const global_size_t numGlobalElements,
384  const Kokkos::View<const global_ordinal_type*, device_type>& indexList,
385  const global_ordinal_type indexBase,
386  const Teuchos::RCP<const Teuchos::Comm<int> >& comm);
387 
429  Map (const global_size_t numGlobalElements,
430  const global_ordinal_type indexList[],
431  const local_ordinal_type indexListSize,
432  const global_ordinal_type indexBase,
433  const Teuchos::RCP<const Teuchos::Comm<int> >& comm);
434 
476  Map (const global_size_t numGlobalElements,
477  const Teuchos::ArrayView<const global_ordinal_type>& indexList,
478  const global_ordinal_type indexBase,
479  const Teuchos::RCP<const Teuchos::Comm<int> >& comm);
480 
481 
495  Map ();
496 
499 
502 
504  Map&
506 
508  Map&
510 
520  virtual ~Map ();
521 
523 
525 
530  bool isOneToOne () const;
531 
538  return numGlobalElements_;
539  }
540 
546  size_t getLocalNumElements () const {
547  return numLocalElements_;
548  }
549 
556  return indexBase_;
557  }
558 
565  return static_cast<local_ordinal_type> (0);
566  }
567 
579  if (this->getLocalNumElements () == 0) {
580  return Tpetra::Details::OrdinalTraits<local_ordinal_type>::invalid ();
581  } else { // Local indices are always zero-based.
582  return static_cast<local_ordinal_type> (this->getLocalNumElements () - 1);
583  }
584  }
585 
592  return minMyGID_;
593  }
594 
601  return maxMyGID_;
602  }
603 
610  return minAllGID_;
611  }
612 
619  return maxAllGID_;
620  }
621 
635 
645 
650  local_map_type getLocalMap () const;
651 
682  getRemoteIndexList (const Teuchos::ArrayView<const global_ordinal_type>& GIDList,
683  const Teuchos::ArrayView< int>& nodeIDList,
684  const Teuchos::ArrayView< local_ordinal_type>& LIDList) const;
685 
710  getRemoteIndexList (const Teuchos::ArrayView<const global_ordinal_type> & GIDList,
711  const Teuchos::ArrayView< int> & nodeIDList) const;
712 
713  private:
725  typedef Kokkos::View<const global_ordinal_type*,
726  Kokkos::LayoutLeft,
727  Kokkos::HostSpace> global_indices_array_type;
728 
729  typedef Kokkos::View<const global_ordinal_type*,
730  device_type> global_indices_array_device_type;
731 
732  public:
752  global_indices_array_type getMyGlobalIndices () const;
753 
755  global_indices_array_device_type getMyGlobalIndicesDevice () const;
756 
757 
768  Teuchos::ArrayView<const global_ordinal_type> getLocalElementList() const;
769 
771 
773 
780  bool isNodeLocalElement (local_ordinal_type localIndex) const;
781 
788  bool isNodeGlobalElement (global_ordinal_type globalIndex) const;
789 
796  bool isUniform () const;
797 
809  bool isContiguous () const;
810 
831  bool isDistributed () const;
832 
858 
890 
896 
913 
915 
917 
919  Teuchos::RCP<const Teuchos::Comm<int> > getComm () const;
920 
921 
923 
925 
927  std::string description () const;
928 
950  void
951  describe (Teuchos::FancyOStream &out,
952  const Teuchos::EVerbosityLevel verbLevel =
953  Teuchos::Describable::verbLevel_default) const;
955 
957 
1006  Teuchos::RCP<const Map<local_ordinal_type, global_ordinal_type, Node> >
1007  removeEmptyProcesses () const;
1008 
1036  Teuchos::RCP<const Map<local_ordinal_type, global_ordinal_type, Node> >
1037  replaceCommWithSubset (const Teuchos::RCP<const Teuchos::Comm<int> >& newComm) const;
1039 
1040  private:
1045  std::string
1046  localDescribeToString (const Teuchos::EVerbosityLevel vl) const;
1047 
1055  void setupDirectory () const;
1056 
1071  bool checkIsDist() const;
1072 
1081  initialNonuniformDebugCheck(
1082  const char errorMessagePrefix[],
1083  const global_size_t numGlobalElements,
1084  const size_t numLocalElements,
1085  const global_ordinal_type indexBase,
1086  const Teuchos::RCP<const Teuchos::Comm<int>>& comm) const;
1087 
1088  void
1089  initWithNonownedHostIndexList(
1090  const char errorMessagePrefix[],
1091  const global_size_t numGlobalElements,
1092  const Kokkos::View<const global_ordinal_type*,
1093  Kokkos::LayoutLeft,
1094  Kokkos::HostSpace,
1095  Kokkos::MemoryUnmanaged>& entryList,
1096  const global_ordinal_type indexBase,
1097  const Teuchos::RCP<const Teuchos::Comm<int>>& comm);
1098 
1100  void lazyPushToHost() const;
1101 
1103  Teuchos::RCP<const Teuchos::Comm<int> > comm_;
1104 
1106  global_ordinal_type indexBase_;
1107 
1110  global_size_t numGlobalElements_;
1111 
1113  size_t numLocalElements_;
1114 
1116  global_ordinal_type minMyGID_;
1117 
1119  global_ordinal_type maxMyGID_;
1120 
1123  global_ordinal_type minAllGID_;
1124 
1127  global_ordinal_type maxAllGID_;
1128 
1135  global_ordinal_type firstContiguousGID_;
1136 
1150  global_ordinal_type lastContiguousGID_;
1151 
1157  bool uniform_;
1158 
1160  bool contiguous_;
1161 
1170  bool distributed_;
1171 
1201  mutable Kokkos::View<const global_ordinal_type*,
1202  Kokkos::LayoutLeft,
1203  device_type> lgMap_;
1204 
1212 #ifndef SWIG
1213  mutable Kokkos::View<const global_ordinal_type*,
1214  Kokkos::LayoutLeft,
1215  Kokkos::HostSpace> lgMapHost_;
1216 #endif
1217 
1219  typedef ::Tpetra::Details::FixedHashTable<global_ordinal_type,
1220  local_ordinal_type, device_type> global_to_local_table_type;
1221 
1234  global_to_local_table_type glMap_;
1235 
1237  typedef ::Tpetra::Details::FixedHashTable<
1238  global_ordinal_type, local_ordinal_type, Kokkos::HostSpace::device_type>
1239  global_to_local_table_host_type;
1240 
1246  mutable global_to_local_table_host_type glMapHost_;
1247 
1284  mutable Teuchos::RCP<
1285  Directory<
1287  >
1288  > directory_;
1289  }; // Map class
1290 
1304  template <class LocalOrdinal, class GlobalOrdinal>
1305  Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal> >
1306  createLocalMap (const size_t numElements,
1307  const Teuchos::RCP<const Teuchos::Comm<int> >& comm);
1308 
1323  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1324  Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >
1325  createLocalMapWithNode (const size_t numElements,
1326  const Teuchos::RCP<const Teuchos::Comm<int> >& comm);
1327 
1328 
1336  template <class LocalOrdinal, class GlobalOrdinal>
1337  Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal> >
1338  createUniformContigMap (const global_size_t numElements,
1339  const Teuchos::RCP<const Teuchos::Comm<int> >& comm);
1340 
1347  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1348  Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >
1349  createUniformContigMapWithNode (const global_size_t numElements,
1350  const Teuchos::RCP<const Teuchos::Comm<int> >& comm);
1351 
1352 
1359  template <class LocalOrdinal, class GlobalOrdinal>
1360  Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal> >
1361  createContigMap (const global_size_t numElements,
1362  const size_t localNumElements,
1363  const Teuchos::RCP<const Teuchos::Comm<int> >& comm);
1364 
1373  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1374  Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >
1375  createContigMapWithNode (const global_size_t numElements,
1376  const size_t localNumElements,
1377  const Teuchos::RCP<const Teuchos::Comm<int> >& comm);
1378 
1379 
1386  template <class LocalOrdinal, class GlobalOrdinal>
1387  Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal> >
1388  createNonContigMap (const Teuchos::ArrayView<const GlobalOrdinal>& elementList,
1389  const Teuchos::RCP<const Teuchos::Comm<int> >& comm);
1390 
1398  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1399  Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal,Node> >
1400  createNonContigMapWithNode (const Teuchos::ArrayView<const GlobalOrdinal> &elementList,
1401  const Teuchos::RCP<const Teuchos::Comm<int> > &comm);
1402 
1410 
1415  template<class LocalOrdinal, class GlobalOrdinal, class Node>
1416  Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal,Node> >
1417  createOneToOne (const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& M);
1418 
1424  template<class LocalOrdinal, class GlobalOrdinal, class Node>
1425  Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal,Node> >
1426  createOneToOne(const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > &M,
1427  const ::Tpetra::Details::TieBreak<LocalOrdinal,GlobalOrdinal> & tie_break);
1428 
1429 } // namespace Tpetra
1430 
1431 #include "Tpetra_Directory_decl.hpp"
1432 
1435 template <class LocalOrdinal, class GlobalOrdinal, class Node>
1438 { return map1.isSameAs (map2); }
1439 
1442 template <class LocalOrdinal, class GlobalOrdinal, class Node>
1445 { return ! map1.isSameAs (map2); }
1446 
1447 
1448 #endif // TPETRA_MAP_DECL_HPP
global_ordinal_type getMaxGlobalIndex() const
The maximum global index owned by the calling process.
LO local_ordinal_type
The type of local indices.
virtual ~Map()
Destructor (virtual for memory safety of derived classes).
GO global_ordinal_type
The type of global indices.
bool isSameAs(const Map< local_ordinal_type, global_ordinal_type, Node > &map) const
True if and only if map is identical to this Map.
global_indices_array_type getMyGlobalIndices() const
Return a view of the global indices owned by this process.
bool isContiguous() const
True if this Map is distributed contiguously, else false.
LookupStatus
Return status of Map remote index lookup (getRemoteIndexList()).
size_t getLocalNumElements() const
The number of elements belonging to the calling process.
bool isLocallyFitted(const Map< local_ordinal_type, global_ordinal_type, Node > &map) const
True if and only if map is locally fitted to this Map.
bool isUniform() const
Whether the range of global indices is uniform.
Declaration and definition of the Tpetra::Map class, an implementation detail of Tpetra::Map.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createContigMapWithNode(const global_size_t numElements, const size_t localNumElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Nonmember constructor for a (potentially) nonuniformly distributed, contiguous Map for a user-specifi...
global_indices_array_device_type getMyGlobalIndicesDevice() const
Return a view of the global indices owned by this process on the Map&#39;s device.
bool isDistributed() const
Whether this Map is globally distributed or locally replicated.
Teuchos::RCP< const Map< local_ordinal_type, global_ordinal_type, Node > > removeEmptyProcesses() const
Advanced methods.
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Accessors for the Teuchos::Comm and Kokkos Node objects.
bool operator!=(const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map1, const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map2)
True if map1 is not the same as (in the sense of isSameAs()) map2, else false.
typename device_type::memory_space memory_space
The Kokkos memory space.
local_ordinal_type getMinLocalIndex() const
The minimum local index.
&quot;Local&quot; part of Map suitable for Kokkos kernels.
Teuchos::ArrayView< const global_ordinal_type > getLocalElementList() const
Return a NONOWNING view of the global indices owned by this process.
size_t global_size_t
Global size_t object.
bool isCompatible(const Map< local_ordinal_type, global_ordinal_type, Node > &map) const
True if and only if map is compatible with this Map.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createOneToOne(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &M)
Nonmember constructor for a contiguous Map with user-defined weights and a user-specified, possibly nondefault Kokkos Node type.
global_ordinal_type getIndexBase() const
The index base for this Map.
global_ordinal_type getGlobalElement(local_ordinal_type localIndex) const
The global index corresponding to the given local index.
Forward declaration of Tpetra::Directory.
local_ordinal_type getMaxLocalIndex() const
The maximum local index on the calling process.
bool isNodeLocalElement(local_ordinal_type localIndex) const
Whether the given local index is valid for this Map on the calling process.
Teuchos::RCP< const Map< local_ordinal_type, global_ordinal_type, Node > > replaceCommWithSubset(const Teuchos::RCP< const Teuchos::Comm< int > > &newComm) const
Replace this Map&#39;s communicator with a subset communicator.
bool isNodeGlobalElement(global_ordinal_type globalIndex) const
Whether the given global index is owned by this Map on the calling process.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createUniformContigMap(const global_size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Non-member constructor for a uniformly distributed, contiguous Map with the default Kokkos Node...
bool operator==(const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map1, const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map2)
True if map1 is the same as (in the sense of isSameAs()) map2, else false.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createContigMap(const global_size_t numElements, const size_t localNumElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Non-member constructor for a (potentially) non-uniformly distributed, contiguous Map using the defaul...
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createLocalMap(const size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Nonmember constructor for a locally replicated Map with the default Kokkos Node.
local_map_type getLocalMap() const
Get the LocalMap for Kokkos-Kernels.
global_ordinal_type getMinGlobalIndex() const
The minimum global index owned by the calling process.
typename device_type::execution_space execution_space
The Kokkos execution space.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createUniformContigMapWithNode(const global_size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Non-member constructor for a uniformly distributed, contiguous Map with a user-specified Kokkos Node...
A parallel distribution of indices over processes.
global_ordinal_type getMaxAllGlobalIndex() const
The maximum global index over all processes in the communicator.
local_ordinal_type getLocalElement(global_ordinal_type globalIndex) const
The local index corresponding to the given global index.
Implement mapping from global ID to process ID and local ID.
bool isOneToOne() const
Whether the Map is one to one.
bool locallySameAs(const Map< local_ordinal_type, global_ordinal_type, node_type > &map) const
Is this Map locally the same as the input Map?
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createLocalMapWithNode(const size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Nonmember constructor for a locally replicated Map with a specified Kokkos Node.
Map & operator=(const Map< local_ordinal_type, global_ordinal_type, node_type > &)=default
Copy assigment (shallow copy).
typename node_type::device_type device_type
This class&#39; Kokkos::Device specialization.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createNonContigMap(const Teuchos::ArrayView< const GlobalOrdinal > &elementList, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Nonmember constructor for a non-contiguous Map using the default Kokkos::Device type.
LocalGlobal
Enum for local versus global allocation of Map entries.
::Tpetra::Details::LocalMap< local_ordinal_type, global_ordinal_type, device_type > local_map_type
Type of the &quot;local&quot; Map.
LookupStatus getRemoteIndexList(const Teuchos::ArrayView< const global_ordinal_type > &GIDList, const Teuchos::ArrayView< int > &nodeIDList, const Teuchos::ArrayView< local_ordinal_type > &LIDList) const
Return the process ranks and corresponding local indices for the given global indices.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Describe this object in a human-readable way to the given output stream.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createNonContigMapWithNode(const Teuchos::ArrayView< const GlobalOrdinal > &elementList, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Nonmember constructor for a noncontiguous Map with a user-specified, possibly nondefault Kokkos Node ...
node_type node_type
Legacy typedef that will go away at some point.
global_ordinal_type getMinAllGlobalIndex() const
The minimum global index over all processes in the communicator.
Forward declaration of Tpetra::Map.
global_size_t getGlobalNumElements() const
The number of elements in this Map.
Forward declaration for Tpetra::TieBreak.
std::string description() const
Implementation of Teuchos::Describable.
Map()
Default constructor (that does nothing).