Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_ReorderedBlockedMultiVector.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 XPETRA_REORDEREDBLOCKEDMULTIVECTOR_HPP
11 #define XPETRA_REORDEREDBLOCKEDMULTIVECTOR_HPP
12 
13 #include <Tpetra_KokkosCompat_DefaultNode.hpp>
14 
15 #include "Xpetra_ConfigDefs.hpp"
16 #include "Xpetra_Exceptions.hpp"
17 
18 #include "Xpetra_MapUtils.hpp"
19 
21 #include "Xpetra_BlockedMap.hpp"
22 #include "Xpetra_BlockedMultiVector.hpp"
23 
28 namespace Xpetra {
29 
30 typedef std::string viewLabel_t;
31 
32 template <class Scalar,
33  class LocalOrdinal,
34  class GlobalOrdinal,
35  class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
36 class ReorderedBlockedMultiVector : public BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
37  public:
38  typedef Scalar scalar_type;
39  typedef LocalOrdinal local_ordinal_type;
40  typedef GlobalOrdinal global_ordinal_type;
41  typedef Node node_type;
42 
43  private:
44 #undef XPETRA_REORDEREDBLOCKEDMULTIVECTOR_SHORT
45 #include "Xpetra_UseShortNames.hpp"
46 
47  public:
49 
50 
52 
59  ReorderedBlockedMultiVector(Teuchos::RCP<const BlockedMap>& rangeMap,
60  Teuchos::RCP<const Xpetra::BlockReorderManager> brm,
62  : Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(rangeMap, bvec->getNumVectors(), false) {
63  brm_ = brm;
64  fullVec_ = bvec;
65  }
66 
67  // protected:
68 
71  brm_ = Teuchos::null;
72  fullVec_ = Teuchos::null;
73  }
74 
76 
77  private:
78  Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> mergeSubBlockMaps(Teuchos::RCP<const Xpetra::BlockReorderManager> brm) {
79  RCP<const BlockedMap> bmap = fullVec_->getBlockedMap();
80 
81  // number of sub blocks
82  size_t numBlocks = brm->GetNumBlocks();
83 
84  Teuchos::RCP<const Map> map = Teuchos::null;
85 
86  if (numBlocks == 0) {
87  // it is a leaf node
88  Teuchos::RCP<const Xpetra::BlockReorderLeaf> leaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(brm);
89 
90  // never extract Thyra style maps (since we have to merge them)
91  map = bmap->getMap(Teuchos::as<size_t>(leaf->GetIndex()), false);
92  } else {
93  // initialize vector for sub maps
94  std::vector<Teuchos::RCP<const Map>> subMaps(numBlocks, Teuchos::null);
95 
96  for (size_t i = 0; i < numBlocks; i++) {
97  Teuchos::RCP<const Xpetra::BlockReorderManager> blkMgr = brm->GetBlock(Teuchos::as<int>(i));
98  subMaps[i] = mergeSubBlockMaps(blkMgr);
99  TEUCHOS_ASSERT(subMaps[i].is_null() == false);
100  }
101 
102  map = MapUtils::concatenateMaps(subMaps);
103  }
104  TEUCHOS_ASSERT(map.is_null() == false);
105  return map;
106  }
107 
108  public:
110 
111 
113  std::string description() const { return "ReorderedBlockedMultiVector"; }
114 
116  void describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default) const {
117  TEUCHOS_ASSERT(brm_ != Teuchos::null);
118  out << description() << ": " << brm_->toString() << std::endl;
119  fullVec_->describe(out, verbLevel);
120  }
121 
123 
124  private:
125  Teuchos::RCP<const Xpetra::BlockReorderManager> brm_;
126  Teuchos::RCP<const Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> fullVec_;
127 };
128 
129 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
130 Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> mergeSubBlockMaps(Teuchos::RCP<const Xpetra::BlockReorderManager> brm, Teuchos::RCP<const Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> bvec, bool bThyraMode) {
132 
133  // TODO distinguish between range and domain map extractor! provide MapExtractor as parameter!
134  RCP<const Xpetra::BlockedMap<LocalOrdinal, GlobalOrdinal, Node>> bmap = bvec->getBlockedMap();
135 
136  // number of sub blocks
137  size_t numBlocks = brm->GetNumBlocks();
138 
139  Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> map = Teuchos::null;
140 
141  if (numBlocks == 0) {
142  // it is a leaf node
143  Teuchos::RCP<const Xpetra::BlockReorderLeaf> leaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(brm);
144 
145  map = bmap->getMap(Teuchos::as<size_t>(leaf->GetIndex()), bThyraMode);
146  } else {
147  // initialize vector for sub maps
148  std::vector<Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>> subMaps(numBlocks, Teuchos::null);
149 
150  for (size_t i = 0; i < numBlocks; i++) {
151  Teuchos::RCP<const Xpetra::BlockReorderManager> blkMgr = brm->GetBlock(Teuchos::as<int>(i));
152  subMaps[i] = mergeSubBlockMaps(blkMgr, bvec, bThyraMode);
153  TEUCHOS_ASSERT(subMaps[i].is_null() == false);
154  }
155 #if 1
156  // concatenate submaps
157  // for Thyra mode this map isn't important
158  Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> fullMap = MapUtils::concatenateMaps(subMaps);
159 
160  // create new BlockedMap (either in Thyra Mode or Xpetra mode)
161  map = Teuchos::rcp(new Xpetra::BlockedMap<LocalOrdinal, GlobalOrdinal, Node>(fullMap, subMaps, bThyraMode));
162 #else
163  // TAW: 11/27/16 we just concatenate the submaps to one monolithic Map object.
164  // Alternatively, we could create a new BlockedMap using the concatenated map and the submaps
165  // However, the block smoothers only need the concatenated map for creating MultiVectors...
166  // But for the Thyra mode version concatenating would not be ok for the whole map
167  map = MapUtils::concatenateMaps(subMaps);
168 #endif
169  }
170  TEUCHOS_ASSERT(map.is_null() == false);
171  return map;
172 }
173 
174 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
175 Teuchos::RCP<const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> mergeSubBlocks(Teuchos::RCP<const Xpetra::BlockReorderManager> rowMgr, Teuchos::RCP<const Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> bvec) {
182 
183  // number of sub blocks
184  size_t rowSz = rowMgr->GetNumBlocks();
185 
186  Teuchos::RCP<BlockedMultiVector> rbvec = Teuchos::null;
187 
188  if (rowSz == 0) {
189  // it is a leaf node
190  Teuchos::RCP<const Xpetra::BlockReorderLeaf> rowleaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(rowMgr);
191 
192  // extract leaf node
193  Teuchos::RCP<MultiVector> vec = bvec->getMultiVector(rowleaf->GetIndex(), false);
194 
195  TEUCHOS_ASSERT(vec != Teuchos::null);
196 
197  // check, whether leaf node is of type Xpetra::CrsMatrixWrap
198  Teuchos::RCP<BlockedMultiVector> subBVec = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(vec);
199  if (subBVec == Teuchos::null) {
200  // DEBUG
201  /*{
202  std::cout << "MultiVector:" << std::endl;
203  Teuchos::ArrayRCP<const Scalar> vData = vec->getData(0);
204  for(size_t j=0; j< vec->getMap()->getLocalNumElements(); j++) {
205  std::cout << j << ": " << vec->getMap()->getGlobalElement(j) << ": " << vData[j] << std::endl;
206  }
207  }*/
208  // END DEBUG
209 
210  // If the leaf node is of type Xpetra::MultiVector. Wrap it into a ReorderedBlockMultiVector
211  // with the corresponding MapExtractors for translating Thyra to Xpetra GIDs if necessary
212  RCP<const BlockedMap> fullBlockedMap = bvec->getBlockedMap();
213  Teuchos::RCP<const Map> submap = fullBlockedMap->getMap(rowleaf->GetIndex(), false);
214  std::vector<Teuchos::RCP<const Map>> rowSubMaps(1, submap);
215  Teuchos::RCP<const BlockedMap> bbmap = Teuchos::rcp(new BlockedMap(submap, rowSubMaps, false));
216 
217  rbvec = Teuchos::rcp(new ReorderedBlockedMultiVector(bbmap, rowMgr, bvec));
218  rbvec->setMultiVector(0, Teuchos::rcp_const_cast<MultiVector>(vec), false);
219 
220  } else {
221  // If leaf node is already wrapped into a blocked matrix do not wrap it again.
222  rbvec = subBVec;
223  TEUCHOS_ASSERT(rbvec != Teuchos::null);
224  }
225  } else {
226  // create the map extractors
227  // we cannot create block matrix in thyra mode since merged maps might not start with 0 GID
228  Teuchos::RCP<const BlockedMap> rgBlockedMap = Teuchos::null;
229  std::vector<Teuchos::RCP<const Map>> rowSubMaps(rowSz, Teuchos::null);
230  for (size_t i = 0; i < rowSz; i++) {
231  Teuchos::RCP<const Xpetra::BlockReorderManager> rowSubMgr = rowMgr->GetBlock(Teuchos::as<int>(i));
232  rowSubMaps[i] = mergeSubBlockMaps(rowSubMgr, bvec, false /*xpetra*/);
233  TEUCHOS_ASSERT(rowSubMaps[i].is_null() == false);
234  }
235  Teuchos::RCP<const Map> rgMergedSubMaps = MapUtils::concatenateMaps(rowSubMaps);
236  rgBlockedMap = Teuchos::rcp(new BlockedMap(rgMergedSubMaps, rowSubMaps, false));
237  rbvec = Teuchos::rcp(new ReorderedBlockedMultiVector(rgBlockedMap, rowMgr, bvec));
238 
239  for (size_t i = 0; i < rowSz; i++) {
240  Teuchos::RCP<const Xpetra::BlockReorderManager> rowSubMgr = rowMgr->GetBlock(Teuchos::as<int>(i));
241  Teuchos::RCP<const MultiVector> subvec = mergeSubBlocks(rowSubMgr, bvec);
242  rbvec->setMultiVector(i, Teuchos::rcp_const_cast<MultiVector>(subvec), false);
243  }
244  }
245  return rbvec;
246 }
247 
248 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
249 Teuchos::RCP<const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> mergeSubBlocksThyra(Teuchos::RCP<const Xpetra::BlockReorderManager> rowMgr, Teuchos::RCP<const Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> bvec) {
255 
256  TEUCHOS_ASSERT(bvec->getBlockedMap()->getThyraMode() == true);
257 
258  // number of sub blocks
259  size_t rowSz = rowMgr->GetNumBlocks();
260 
261  Teuchos::RCP<BlockedMultiVector> rbvec = Teuchos::null;
262 
263  if (rowSz == 0) {
264  // it is a leaf node
265  Teuchos::RCP<const Xpetra::BlockReorderLeaf> rowleaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(rowMgr);
266 
267  // this MultiVector uses Thyra style GIDs as global row indices
268  Teuchos::RCP<MultiVector> vec = bvec->getMultiVector(rowleaf->GetIndex(), true);
269 
270  TEUCHOS_ASSERT(vec.is_null() == false);
271 
272  // check, whether leaf node is of type Xpetra::CrsMatrixWrap
273  Teuchos::RCP<BlockedMultiVector> bbvec = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(vec);
274  if (bbvec == Teuchos::null) {
276  // build blocked map
277  RCP<const BlockedMap> fullBlockedRangeMap = bvec->getBlockedMap();
278  // extract Xpetra and Thyra based GIDs
279  Teuchos::RCP<const Map> xpsubmap = fullBlockedRangeMap->getMap(rowleaf->GetIndex(), false);
280  Teuchos::RCP<const Map> thysubmap = fullBlockedRangeMap->getMap(rowleaf->GetIndex(), true);
281  std::vector<Teuchos::RCP<const Map>> rowXpSubMaps(1, xpsubmap);
282  std::vector<Teuchos::RCP<const Map>> rowTySubMaps(1, thysubmap);
283  // use expert constructor
284  Teuchos::RCP<const BlockedMap> rgBlockedMap = Teuchos::rcp(new BlockedMap(rowXpSubMaps, rowTySubMaps));
285 
287  // build reordered blocked multi vector
288  rbvec = Teuchos::rcp(new ReorderedBlockedMultiVector(rgBlockedMap, rowMgr, bvec));
289  rbvec->setMultiVector(0, vec, true);
290  } else {
291  // If leaf node is already wrapped into a blocked matrix do not wrap it again.
292  rbvec = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(vec);
293  }
294  } else {
295  // create the blocked map
296  // we cannot create block multivector in thyra mode since merged maps might not start with 0 GID
297 
298  std::vector<Teuchos::RCP<const Map>> rowXpSubMaps(rowSz, Teuchos::null);
299  std::vector<Teuchos::RCP<const Map>> rowTySubMaps(rowSz, Teuchos::null);
300  for (size_t i = 0; i < rowSz; i++) {
301  Teuchos::RCP<const Xpetra::BlockReorderManager> rowSubMgr = rowMgr->GetBlock(Teuchos::as<int>(i));
302  // extract Xpetra and Thyra based merged GIDs
303  rowXpSubMaps[i] = mergeSubBlockMaps(rowSubMgr, bvec, false);
304  rowTySubMaps[i] = mergeSubBlockMaps(rowSubMgr, bvec, true);
305  TEUCHOS_ASSERT(rowXpSubMaps[i].is_null() == false);
306  TEUCHOS_ASSERT(rowTySubMaps[i].is_null() == false);
307  }
308  // use expert constructor
309  Teuchos::RCP<const BlockedMap> rgBlockedMap = Teuchos::rcp(new BlockedMap(rowXpSubMaps, rowTySubMaps));
310 
311  rbvec = Teuchos::rcp(new ReorderedBlockedMultiVector(rgBlockedMap, rowMgr, bvec));
312 
313  for (size_t i = 0; i < rowSz; i++) {
314  Teuchos::RCP<const Xpetra::BlockReorderManager> rowSubMgr = rowMgr->GetBlock(Teuchos::as<int>(i));
315  Teuchos::RCP<const MultiVector> subvec = mergeSubBlocksThyra(rowSubMgr, bvec);
316  rbvec->setMultiVector(i, Teuchos::rcp_const_cast<MultiVector>(subvec), true);
317  }
318  }
319  return rbvec;
320 }
321 
322 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
323 Teuchos::RCP<const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> buildReorderedBlockedMultiVector(Teuchos::RCP<const Xpetra::BlockReorderManager> brm, Teuchos::RCP<const Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> bvec) {
324  Teuchos::RCP<const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> rbvec = Teuchos::null;
325  if (bvec->getBlockedMap()->getThyraMode() == false) {
326  rbvec = mergeSubBlocks(brm, bvec);
327  } else {
328  rbvec = mergeSubBlocksThyra(brm, bvec);
329  }
330  TEUCHOS_ASSERT(rbvec.is_null() == false);
331  return rbvec;
332 }
333 
334 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
335 Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> buildReorderedBlockedMultiVector(Teuchos::RCP<const Xpetra::BlockReorderManager> brm, Teuchos::RCP<Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> bvec) {
338  Teuchos::RCP<const MultiVector> rbvec = Teuchos::null;
339  if (bvec->getBlockedMap()->getThyraMode() == false) {
340  rbvec = mergeSubBlocks(brm, Teuchos::rcp_const_cast<const BlockedMultiVector>(bvec));
341  } else {
342  rbvec = mergeSubBlocksThyra(brm, Teuchos::rcp_const_cast<const BlockedMultiVector>(bvec));
343  }
344  TEUCHOS_ASSERT(rbvec.is_null() == false);
345  return Teuchos::rcp_const_cast<MultiVector>(rbvec);
346 }
347 
348 } // namespace Xpetra
349 
350 #define XPETRA_REORDEREDBLOCKEDMULTIVECTOR_SHORT
351 #endif /* XPETRA_REORDEREDBLOCKEDMULTIVECTOR_HPP */
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > mergeSubBlockMaps(Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> bmat, bool bThyraMode)
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > mergeSubBlockMaps(Teuchos::RCP< const Xpetra::BlockReorderManager > brm)
Xpetra utility class for common map-related routines.
Teuchos::RCP< const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > buildReorderedBlockedMultiVector(Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >> bvec)
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.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
ReorderedBlockedMultiVector(Teuchos::RCP< const BlockedMap > &rangeMap, Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >> bvec)
Constructor.
Teuchos::RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mergeSubBlocks(Teuchos::RCP< const Xpetra::BlockReorderManager > rowMgr, Teuchos::RCP< const Xpetra::BlockReorderManager > colMgr, Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> bmat)
Teuchos::RCP< const Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > fullVec_
virtual size_t getNumVectors() const
Number of columns in the multivector.
Teuchos::RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mergeSubBlocksThyra(Teuchos::RCP< const Xpetra::BlockReorderManager > rowMgr, Teuchos::RCP< const Xpetra::BlockReorderManager > colMgr, Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> bmat)
Teuchos::RCP< const Xpetra::BlockReorderManager > brm_
std::string description() const
Return a simple one-line description of this object.