Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_StridedMapFactory_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 // WARNING: This code is experimental. Backwards compatibility should not be expected.
11 
12 #ifndef XPETRA_STRIDEDMAPFACTORY_DEF_HPP
13 #define XPETRA_STRIDEDMAPFACTORY_DEF_HPP
14 
16 
17 #include "Xpetra_Exceptions.hpp"
18 
19 namespace Xpetra {
20 
21 template <class LocalOrdinal, class GlobalOrdinal, class Node>
22 RCP<Xpetra::StridedMap<LocalOrdinal, GlobalOrdinal, Node>>
25  global_size_t numGlobalElements,
26  GlobalOrdinal indexBase,
27  std::vector<size_t>& stridingInfo,
28  const Teuchos::RCP<const Teuchos::Comm<int>>& comm,
29  LocalOrdinal stridedBlockId,
30  GlobalOrdinal offset,
31  LocalGlobal lg) {
32  return rcp(new Xpetra::StridedMap<LocalOrdinal, GlobalOrdinal, Node>(lib, numGlobalElements, indexBase, stridingInfo, comm, stridedBlockId, offset, lg));
33 }
34 
35 template <class LocalOrdinal, class GlobalOrdinal, class Node>
36 RCP<Xpetra::StridedMap<LocalOrdinal, GlobalOrdinal, Node>>
39  global_size_t numGlobalElements,
40  size_t numLocalElements,
41  GlobalOrdinal indexBase,
42  std::vector<size_t>& stridingInfo,
43  const Teuchos::RCP<const Teuchos::Comm<int>>& comm,
44  LocalOrdinal stridedBlockId,
45  GlobalOrdinal offset) {
46  return rcp(new StridedMap(lib, numGlobalElements, numLocalElements, indexBase, stridingInfo, comm, stridedBlockId, offset));
47 }
48 
49 template <class LocalOrdinal, class GlobalOrdinal, class Node>
50 RCP<Xpetra::StridedMap<LocalOrdinal, GlobalOrdinal, Node>>
52  Build(const RCP<const Map>& map,
53  std::vector<size_t>& stridingInfo,
54  LocalOrdinal stridedBlockId,
55  GlobalOrdinal offset) {
56  return rcp(new StridedMap(map, stridingInfo, map->getIndexBase(), stridedBlockId, offset));
57 }
58 
59 template <class LocalOrdinal, class GlobalOrdinal, class Node>
60 RCP<Xpetra::StridedMap<LocalOrdinal, GlobalOrdinal, Node>>
62  Build(const RCP<const StridedMap>& map, LocalOrdinal stridedBlockId) {
63  TEUCHOS_TEST_FOR_EXCEPTION(stridedBlockId < 0, Exceptions::RuntimeError,
64  "Xpetra::StridedMapFactory::Build: constructor expects stridedBlockId > -1.");
65 
66  TEUCHOS_TEST_FOR_EXCEPTION(map->getStridedBlockId() != -1, Exceptions::RuntimeError,
67  "Xpetra::StridedMapFactory::Build: constructor expects a full map (stridedBlockId == -1).");
68 
69  std::vector<size_t> stridingInfo = map->getStridingData();
70 
71  Teuchos::ArrayView<const GlobalOrdinal> dofGids = map->getLocalElementList();
72 
73  // determine nStridedOffset
74  size_t nStridedOffset = 0;
75  for (int j = 0; j < map->getStridedBlockId(); j++) {
76  nStridedOffset += stridingInfo[j];
77  }
78 
79  const size_t numMyBlockDofs = (stridingInfo[stridedBlockId] * map->getLocalNumElements()) / map->getFixedBlockSize();
80 
81  std::vector<GlobalOrdinal> subBlockDofGids(numMyBlockDofs);
82 
83  // TODO fill vector with dofs
84  LocalOrdinal ind = 0;
85  for (typename Teuchos::ArrayView<const GlobalOrdinal>::iterator it = dofGids.begin(); it != dofGids.end(); ++it) {
86  if (map->GID2StridingBlockId(*it) == Teuchos::as<size_t>(stridedBlockId)) {
87  subBlockDofGids[ind++] = *it;
88  }
89  }
90 
91  const Teuchos::ArrayView<const GlobalOrdinal> subBlockDofGids_view(&subBlockDofGids[0], subBlockDofGids.size());
92 
93  return rcp(new StridedMap(map->lib(),
94  Teuchos::OrdinalTraits<global_size_t>::invalid(),
95  subBlockDofGids_view,
96  map->getIndexBase(),
97  stridingInfo,
98  map->getComm(),
99  stridedBlockId));
100 }
101 
102 template <class LocalOrdinal, class GlobalOrdinal, class Node>
103 RCP<Xpetra::StridedMap<LocalOrdinal, GlobalOrdinal, Node>>
105  Build(const StridedMap& map) {
106  XPETRA_MONITOR("MapFactory::Build");
107 
108  LocalOrdinal N = map.getLocalNumElements();
109  Teuchos::ArrayView<const GlobalOrdinal> oldElements = map.getLocalElementList();
110  Teuchos::Array<GlobalOrdinal> newElements(map.getLocalNumElements());
111 
112  for (LocalOrdinal i = 0; i < N; i++)
113  newElements[i] = oldElements[i];
114 
115  std::vector<size_t> strData = map.getStridingData();
116  return rcp(new StridedMap(map.lib(), map.getGlobalNumElements(), newElements, map.getIndexBase(), strData, map.getComm(), map.getStridedBlockId()));
117 }
118 
119 template <class LocalOrdinal, class GlobalOrdinal, class Node>
120 RCP<Xpetra::StridedMap<LocalOrdinal, GlobalOrdinal, Node>>
123  global_size_t numGlobalElements,
124  const Teuchos::ArrayView<const GlobalOrdinal>& elementList,
125  GlobalOrdinal indexBase,
126  std::vector<size_t>& stridingInfo,
127  const Teuchos::RCP<const Teuchos::Comm<int>>& comm,
128  LocalOrdinal stridedBlockId, // FIXME (mfh 03 Sep 2014) This breaks if LocalOrdinal is unsigned
129  GlobalOrdinal /* offset */) {
130  return rcp(new StridedMap(lib, numGlobalElements, elementList, indexBase, stridingInfo, comm, stridedBlockId));
131 }
132 
133 } // namespace Xpetra
134 
135 #endif // XPETRA_STRIDEDMAPFACTORY_DEF_HPP
136 
137 // TODO: removed unused methods
Teuchos::ArrayView< const GlobalOrdinal > getLocalElementList() const
Return a list of the global indices owned by this node.
LocalOrdinal getStridedBlockId() const
size_t getLocalNumElements() const
Returns the number of elements belonging to the calling node.
Exception throws to report errors in the internal logical of the program.
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Get the Comm object for this Map.
size_t global_size_t
Global size_t object.
GlobalOrdinal getIndexBase() const
Returns the index base for this Map.
global_size_t getGlobalNumElements() const
Returns the number of elements in this Map.
std::vector< size_t > getStridingData() const
UnderlyingLib lib() const
Get the library used by this object (Tpetra or Epetra?)
#define XPETRA_MONITOR(funcName)
static RCP< Xpetra::StridedMap< LocalOrdinal, GlobalOrdinal, Node > > Build(UnderlyingLib lib, global_size_t numGlobalElements, GlobalOrdinal indexBase, std::vector< size_t > &stridingInfo, const Teuchos::RCP< const Teuchos::Comm< int >> &comm, LocalOrdinal stridedBlockId=-1, GlobalOrdinal offset=0, LocalGlobal lg=Xpetra::GloballyDistributed)
Map constructor with Xpetra-defined contiguous uniform distribution.
Class that stores a strided map.