All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Xpetra_ReorderedBlockedCrsMatrix.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Xpetra: A linear algebra interface package
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef XPETRA_REORDEREDBLOCKEDCRSMATRIX_HPP
47 #define XPETRA_REORDEREDBLOCKEDCRSMATRIX_HPP
48 
49 #include <Kokkos_DefaultNode.hpp>
50 
51 #include "Xpetra_ConfigDefs.hpp"
52 #include "Xpetra_Exceptions.hpp"
53 
54 #include "Xpetra_MapUtils.hpp"
55 
57 #include "Xpetra_CrsMatrixWrap.hpp"
59 
60 
65 namespace Xpetra {
66 
67  typedef std::string viewLabel_t;
68 
69  template <class Scalar,
70  class LocalOrdinal,
71  class GlobalOrdinal,
72  class Node>
74  public BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
75  public:
76  typedef Scalar scalar_type;
77  typedef LocalOrdinal local_ordinal_type;
78  typedef GlobalOrdinal global_ordinal_type;
79  typedef Node node_type;
80 
81  private:
82 #undef XPETRA_REORDEREDBLOCKEDCRSMATRIX_SHORT
83 #include "Xpetra_UseShortNames.hpp"
84 
85  public:
86 
88 
89 
91 
102  size_t npr,
106  : Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>(rangeMaps, domainMaps, npr,pftype) {
107  brm_ = brm;
108  fullOp_ = bmat;
109  }
110 
111  //protected:
112 
115 
117 
118  private:
120  RCP<const MapExtractor> fullRangeMapExtractor = fullOp_->getRangeMapExtractor();
121 
122  // number of sub blocks
123  size_t numBlocks = brm->GetNumBlocks();
124 
125  Teuchos::RCP<const Map> map = Teuchos::null;
126 
127  if(numBlocks == 0) {
128  // it is a leaf node
129  Teuchos::RCP<const Xpetra::BlockReorderLeaf> leaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(brm);
130 
131  // never extract Thyra style maps (since we have to merge them)
132  map = fullRangeMapExtractor->getMap(Teuchos::as<size_t>(leaf->GetIndex()), false);
133  } else {
134  // initialize vector for sub maps
135  std::vector<Teuchos::RCP<const Map> > subMaps (numBlocks, Teuchos::null);
136 
137  for(size_t i = 0; i < numBlocks; i++) {
138  Teuchos::RCP<const Xpetra::BlockReorderManager> blkMgr = brm->GetBlock(Teuchos::as<int>(i));
139  subMaps[i] = mergeSubBlockMaps(blkMgr);
140  TEUCHOS_ASSERT(subMaps[i].is_null()==false);
141  }
142 
143  map = MapUtils::concatenateMaps(subMaps);
144  }
145  TEUCHOS_ASSERT(map.is_null()==false);
146  return map;
147  }
148 
149  public:
151 
152 
154 
157  virtual void apply(const MultiVector& X, MultiVector& Y,
159  Scalar alpha = ScalarTraits<Scalar>::one(),
160  Scalar beta = ScalarTraits<Scalar>::zero()) const
161  {
162 
163  // Nested sub-operators should just use the provided X and B vectors
164  if(fullOp_->getNodeNumRows() != this->getNodeNumRows()) {
166  return;
167  }
168 
169  // Special handling for the top level of the nested operator
170 
171  // check whether input parameters are blocked or not
172  RCP<const MultiVector> refX = rcpFromRef(X);
173  RCP<const BlockedMultiVector> refbX = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(refX);
174  RCP<MultiVector> tmpY = rcpFromRef(Y);
175  RCP<BlockedMultiVector> tmpbY = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(tmpY);
176 
177 
178  bool bCopyResultX = false;
179  bool bCopyResultY = false;
180 
181  // TODO create a nested ReorderedBlockedMultiVector out of a nested BlockedMultiVector!
182  // probably not necessary, is it?
183  // check whether X and B are blocked but not ReorderedBlocked operators
184  /*if (refbX != Teuchos::null && fullOp_->getNodeNumRows() == this->getNodeNumRows()) {
185  RCP<const ReorderedBlockedMultiVector> rbCheck = Teuchos::rcp_dynamic_cast<const ReorderedBlockedMultiVector>(refbX);
186  if(rbCheck == Teuchos::null) {
187  RCP<const BlockedMultiVector> bX =
188  Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(Xpetra::buildReorderedBlockedMultiVector(brm_, refbX));
189  TEUCHOS_ASSERT(bX.is_null()==false);
190  refbX.swap(bX);
191  }
192  }
193  if (tmpbY != Teuchos::null && fullOp_->getNodeNumRows() == this->getNodeNumRows()) {
194  RCP<ReorderedBlockedMultiVector> rbCheck = Teuchos::rcp_dynamic_cast<ReorderedBlockedMultiVector>(tmpbY);
195  if(rbCheck == Teuchos::null) {
196  RCP<BlockedMultiVector> bY =
197  Teuchos::rcp_dynamic_cast<BlockedMultiVector>(Xpetra::buildReorderedBlockedMultiVector(brm_, tmpbY));
198  TEUCHOS_ASSERT(bY.is_null()==false);
199  tmpbY.swap(bY);
200  }
201  }*/
202 
203  // if X and B are not blocked, create a blocked version here and use the blocked vectors
204  // for the internal (nested) apply call.
205 
206  // Check whether "this" operator is the reordered variant of the underlying fullOp_.
207  // Note, that nested ReorderedBlockedCrsMatrices always have the same full operator "fullOp_"
208  // stored underneath for being able to "translate" the block ids.
209  if (refbX == Teuchos::null && fullOp_->getNodeNumRows() == this->getNodeNumRows()) {
210  // create a new (non-nested) blocked multi vector (using the blocked range map of fullOp_)
211  RCP<const BlockedMap> blkRgMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(fullOp_->getRangeMap());
212  TEUCHOS_ASSERT(blkRgMap.is_null()==false);
213  RCP<const BlockedMultiVector> bXtemp = Teuchos::rcp(new BlockedMultiVector(blkRgMap, refX));
214  TEUCHOS_ASSERT(bXtemp.is_null()==false);
216  Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(Xpetra::buildReorderedBlockedMultiVector(brm_, bXtemp));
217  TEUCHOS_ASSERT(bX.is_null()==false);
218  refbX.swap(bX);
219  bCopyResultX = true;
220  }
221 
222  if (tmpbY == Teuchos::null && fullOp_->getNodeNumRows() == this->getNodeNumRows()) {
223  // create a new (non-nested) blocked multi vector (using the blocked range map of fullOp_)
224  RCP<const BlockedMap> blkRgMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(fullOp_->getRangeMap());
225  TEUCHOS_ASSERT(blkRgMap.is_null()==false);
226  RCP<BlockedMultiVector> tmpbYtemp = Teuchos::rcp(new BlockedMultiVector(blkRgMap, tmpY));
227  TEUCHOS_ASSERT(tmpbYtemp.is_null()==false);
229  Teuchos::rcp_dynamic_cast<BlockedMultiVector>(Xpetra::buildReorderedBlockedMultiVector(brm_, tmpbYtemp));
230  TEUCHOS_ASSERT(bY.is_null()==false);
231  tmpbY.swap(bY);
232  bCopyResultY = true;
233  }
234 
235  TEUCHOS_ASSERT(refbX.is_null()==false);
236  TEUCHOS_ASSERT(tmpbY.is_null()==false);
237 
239 
240  if (bCopyResultX == true) {
241  RCP<const MultiVector> Xmerged = refbX->Merge();
242  RCP<MultiVector> nonconstX = Teuchos::rcp_const_cast<MultiVector>(refX);
244  }
245  if (bCopyResultY == true) {
246  RCP< MultiVector> Ymerged = tmpbY->Merge();
248  }
249 
250  }
251 
252  // @}
253 
254 
255 
257 
258 
261 
262  // @}
263 
265 
266 
268  std::string description() const { return "ReorderedBlockedCrsMatrix"; }
269 
273 
274  out << "Xpetra::ReorderedBlockedCrsMatrix: " << BlockedCrsMatrix::Rows() << " x " << BlockedCrsMatrix::Cols() << std::endl;
275 
277  out << "ReorderedBlockMatrix is fillComplete" << std::endl;
278 
279  out << "fullRowMap" << std::endl;
280  BlockedCrsMatrix::getRangeMap(0,false)->describe(out,verbLevel);
281 
282  //out << "fullColMap" << std::endl;
283  //fullcolmap_->describe(out,verbLevel);
284 
285  } else {
286  out << "Xpetra::ReorderedBlockedCrsMatrix is NOT fillComplete" << std::endl;
287  }
288 
289  for (size_t r = 0; r < BlockedCrsMatrix::Rows(); ++r)
290  for (size_t c = 0; c < BlockedCrsMatrix::Cols(); ++c) {
291  out << "Block(" << r << "," << c << ")" << std::endl;
292  BlockedCrsMatrix::getMatrix(r,c)->describe(out,verbLevel);
293  }
294  }
295 
297 
298  private:
301 
302 
303 };
304 
305 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
308 
309  // TODO distinguish between range and domain map extractor! provide MapExtractor as parameter!
310  RCP<const Xpetra::MapExtractor<Scalar,LocalOrdinal,GlobalOrdinal,Node> > fullRangeMapExtractor = bmat->getRangeMapExtractor();
311 
312  // number of sub blocks
313  size_t numBlocks = brm->GetNumBlocks();
314 
316 
317  if(numBlocks == 0) {
318  // it is a leaf node
319  Teuchos::RCP<const Xpetra::BlockReorderLeaf> leaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(brm);
320 
321  map = fullRangeMapExtractor->getMap(Teuchos::as<size_t>(leaf->GetIndex()), bThyraMode);
322  } else {
323  // initialize vector for sub maps
324  std::vector<Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > > subMaps (numBlocks, Teuchos::null);
325 
326  for(size_t i = 0; i < numBlocks; i++) {
327  Teuchos::RCP<const Xpetra::BlockReorderManager> blkMgr = brm->GetBlock(Teuchos::as<int>(i));
328  subMaps[i] = mergeSubBlockMaps(blkMgr,bmat,bThyraMode);
329  TEUCHOS_ASSERT(subMaps[i].is_null()==false);
330  }
331 
332 #if 1
333  // concatenate submaps
334  // for Thyra mode this map isn't important
336 
337  // create new BlockedMap (either in Thyra Mode or Xpetra mode)
338  map = Teuchos::rcp(new Xpetra::BlockedMap<LocalOrdinal,GlobalOrdinal,Node>(fullMap, subMaps, bThyraMode));
339 #else
340  // TAW: 11/27/16 we just concatenate the submaps to one monolithic Map object.
341  // Alternatively, we could create a new BlockedMap using the concatenated map and the submaps
342  // However, the block smoothers only need the concatenated map for creating MultiVectors...
343  // But for the Thyra mode version concatenating would not be ok for the whole map
344  map = MapUtils::concatenateMaps(subMaps);
345 #endif
346  }
347  TEUCHOS_ASSERT(map.is_null()==false);
348  return map;
349 }
350 
351 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
353 
360 
361  // number of sub blocks
362  size_t rowSz = rowMgr->GetNumBlocks();
363  size_t colSz = colMgr->GetNumBlocks();
364 
365  Teuchos::RCP<BlockedCrsMatrix> rbmat = Teuchos::null;
366 
367  if(rowSz == 0 && colSz == 0) {
368  // it is a leaf node
369  Teuchos::RCP<const Xpetra::BlockReorderLeaf> rowleaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(rowMgr);
370  Teuchos::RCP<const Xpetra::BlockReorderLeaf> colleaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(colMgr);
371 
372  // extract leaf node
373  Teuchos::RCP<Matrix> mat = bmat->getMatrix(rowleaf->GetIndex(), colleaf->GetIndex());
374 
375  if (mat == Teuchos::null) return Teuchos::null;
376 
377  // check, whether leaf node is of type Xpetra::CrsMatrixWrap
379  if(matwrap != Teuchos::null) {
380  // If the leaf node is of type Xpetra::CrsMatrixWrap wrap it into a 1x1 ReorderedBlockMatrix
381  // with the corresponding MapExtractors for translating Thyra to Xpetra GIDs if necessary
382  RCP<const MapExtractor> fullRangeMapExtractor = bmat->getRangeMapExtractor();
383  Teuchos::RCP<const Map> submap = fullRangeMapExtractor->getMap(rowleaf->GetIndex(),false);
384  std::vector<Teuchos::RCP<const Map> > rowSubMaps (1, submap);
385  Teuchos::RCP<const MapExtractor> rgMapExtractor = Teuchos::rcp(new MapExtractor(submap, rowSubMaps, false));
386 
387  RCP<const MapExtractor> fullDomainMapExtractor = bmat->getDomainMapExtractor();
388  Teuchos::RCP<const Map> submap2 = fullDomainMapExtractor->getMap(colleaf->GetIndex(),false);
389  std::vector<Teuchos::RCP<const Map> > colSubMaps (1, submap2);
390  Teuchos::RCP<const MapExtractor> doMapExtractor = Teuchos::rcp(new MapExtractor(submap2, colSubMaps, false));
391 
392  rbmat = Teuchos::rcp(new ReorderedBlockedCrsMatrix(rgMapExtractor,doMapExtractor, 33, rowMgr,bmat));
393  rbmat->setMatrix(0,0,mat);
394  } else {
395  // If leaf node is already wrapped into a blocked matrix do not wrap it again.
396  rbmat = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(mat);
397  TEUCHOS_ASSERT(rbmat != Teuchos::null);
398  }
399  TEUCHOS_ASSERT(mat->getNodeNumEntries() == rbmat->getNodeNumEntries());
400  } else {
401  // create the map extractors
402  // we cannot create block matrix in thyra mode since merged maps might not start with 0 GID
403  Teuchos::RCP<const MapExtractor> rgMapExtractor = Teuchos::null;
404  if(rowSz > 0) {
405  std::vector<Teuchos::RCP<const Map> > rowSubMaps (rowSz, Teuchos::null);
406  for(size_t i = 0; i < rowSz; i++) {
407  Teuchos::RCP<const Xpetra::BlockReorderManager> rowSubMgr = rowMgr->GetBlock(Teuchos::as<int>(i));
408  rowSubMaps[i] = mergeSubBlockMaps(rowSubMgr,bmat,false /*xpetra*/);
409  TEUCHOS_ASSERT(rowSubMaps[i].is_null()==false);
410  }
411  Teuchos::RCP<const Map> rgMergedSubMaps = MapUtils::concatenateMaps(rowSubMaps);
412  rgMapExtractor = Teuchos::rcp(new MapExtractor(rgMergedSubMaps, rowSubMaps, false));
413  } else {
414  Teuchos::RCP<const Xpetra::BlockReorderLeaf> rowleaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(rowMgr);
415  RCP<const MapExtractor> fullRangeMapExtractor = bmat->getRangeMapExtractor();
416  // TODO think about Thyra style maps: we cannot use thyra style maps when recombining several blocks!!!
417  // The GIDs might not start with 0 and may not be consecutive!
418  Teuchos::RCP<const Map> submap = fullRangeMapExtractor->getMap(rowleaf->GetIndex(),false);
419  std::vector<Teuchos::RCP<const Map> > rowSubMaps (1, submap);
420  rgMapExtractor = Teuchos::rcp(new MapExtractor(submap, rowSubMaps, false));
421  }
422 
423  Teuchos::RCP<const MapExtractor> doMapExtractor = Teuchos::null;
424  if(colSz > 0) {
425  std::vector<Teuchos::RCP<const Map> > colSubMaps (colSz, Teuchos::null);
426  for(size_t j = 0; j < colSz; j++) {
427  Teuchos::RCP<const Xpetra::BlockReorderManager> colSubMgr = colMgr->GetBlock(Teuchos::as<int>(j));
428  colSubMaps[j] = mergeSubBlockMaps(colSubMgr,bmat,false/*xpetra*/);
429  TEUCHOS_ASSERT(colSubMaps[j].is_null()==false);
430  }
431  Teuchos::RCP<const Map> doMergedSubMaps = MapUtils::concatenateMaps(colSubMaps);
432  doMapExtractor = Teuchos::rcp(new MapExtractor(doMergedSubMaps, colSubMaps, false));
433  } else {
434  Teuchos::RCP<const Xpetra::BlockReorderLeaf> colleaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(colMgr);
435  RCP<const MapExtractor> fullDomainMapExtractor = bmat->getDomainMapExtractor();
436  // TODO think about Thyra style maps: we cannot use thyra style maps when recombining several blocks!!!
437  // The GIDs might not start with 0 and may not be consecutive!
438  Teuchos::RCP<const Map> submap = fullDomainMapExtractor->getMap(colleaf->GetIndex(),false);
439  std::vector<Teuchos::RCP<const Map> > colSubMaps (1, submap);
440  doMapExtractor = Teuchos::rcp(new MapExtractor(submap, colSubMaps, false));
441  }
442 
443  rbmat = Teuchos::rcp(new ReorderedBlockedCrsMatrix(rgMapExtractor,doMapExtractor, 33, rowMgr,bmat));
444 
445  size_t cntNNZ = 0;
446 
447  if (rowSz == 0 && colSz > 0) {
448  for(size_t j = 0; j < colSz; j++) {
449  Teuchos::RCP<const Xpetra::BlockReorderManager> colSubMgr = colMgr->GetBlock(Teuchos::as<int>(j));
450  Teuchos::RCP<const Matrix> submat = mergeSubBlocks(rowMgr, colSubMgr, bmat);
451  rbmat->setMatrix(0,j,Teuchos::rcp_const_cast<Matrix>(submat));
452  if(submat != Teuchos::null) cntNNZ += submat->getNodeNumEntries();
453  }
454  } else if (rowSz > 0 && colSz == 0) {
455  for(size_t i = 0; i < rowSz; i++) {
456  Teuchos::RCP<const Xpetra::BlockReorderManager> rowSubMgr = rowMgr->GetBlock(Teuchos::as<int>(i));
457  Teuchos::RCP<const Matrix> submat = mergeSubBlocks(rowSubMgr, colMgr, bmat);
458  rbmat->setMatrix(i,0,Teuchos::rcp_const_cast<Matrix>(submat));
459  if(submat != Teuchos::null) cntNNZ += submat->getNodeNumEntries();
460  }
461  } else {
462  for(size_t i = 0; i < rowSz; i++) {
463  Teuchos::RCP<const Xpetra::BlockReorderManager> rowSubMgr = rowMgr->GetBlock(Teuchos::as<int>(i));
464  for(size_t j = 0; j < colSz; j++) {
465  Teuchos::RCP<const Xpetra::BlockReorderManager> colSubMgr = colMgr->GetBlock(Teuchos::as<int>(j));
466  Teuchos::RCP<const Matrix> submat = mergeSubBlocks(rowSubMgr, colSubMgr, bmat);
467  rbmat->setMatrix(i,j,Teuchos::rcp_const_cast<Matrix>(submat));
468  if(submat != Teuchos::null) cntNNZ += submat->getNodeNumEntries();
469  }
470  }
471  }
472  TEUCHOS_ASSERT(rbmat->getNodeNumEntries() == cntNNZ);
473  }
474  rbmat->fillComplete();
475  return rbmat;
476 }
477 
478  //MapExtractor(const std::vector<RCP<const Map> >& maps, const std::vector<RCP<const Map> >& thyramaps);
479 
480 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
482 
488 
489  TEUCHOS_ASSERT(bmat->getRangeMapExtractor()->getThyraMode() == true);
490  TEUCHOS_ASSERT(bmat->getDomainMapExtractor()->getThyraMode() == true);
491 
492  // number of sub blocks
493  size_t rowSz = rowMgr->GetNumBlocks();
494  size_t colSz = colMgr->GetNumBlocks();
495 
496  Teuchos::RCP<BlockedCrsMatrix> rbmat = Teuchos::null;
497 
498  if(rowSz == 0 && colSz == 0) {
499  // it is a leaf node
500  Teuchos::RCP<const Xpetra::BlockReorderLeaf> rowleaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(rowMgr);
501  Teuchos::RCP<const Xpetra::BlockReorderLeaf> colleaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(colMgr);
502 
503  // this matrix uses Thyra style GIDs as global row, range, domain and column indices
504  Teuchos::RCP<Matrix> mat = bmat->getMatrix(rowleaf->GetIndex(), colleaf->GetIndex());
505 
506  if(mat == Teuchos::null) return Teuchos::null; //std::cout << "Block " << rowleaf->GetIndex() << "," << colleaf->GetIndex() << " is zero" << std::endl;
507 
508  // check, whether leaf node is of type Xpetra::CrsMatrixWrap
510  if(matwrap != Teuchos::null) {
512  // build map extractors
513  RCP<const MapExtractor> fullRangeMapExtractor = bmat->getRangeMapExtractor();
514  // extract Xpetra and Thyra based GIDs
515  Teuchos::RCP<const Map> xpsubmap = fullRangeMapExtractor->getMap(rowleaf->GetIndex(),false);
516  Teuchos::RCP<const Map> thysubmap = fullRangeMapExtractor->getMap(rowleaf->GetIndex(),true);
517  std::vector<Teuchos::RCP<const Map> > rowXpSubMaps (1, xpsubmap);
518  std::vector<Teuchos::RCP<const Map> > rowTySubMaps (1, thysubmap);
519  // use expert constructor
520  Teuchos::RCP<const MapExtractor> rgMapExtractor = Teuchos::rcp(new MapExtractor(rowXpSubMaps, rowTySubMaps));
521 
522  RCP<const MapExtractor> fullDomainMapExtractor = bmat->getDomainMapExtractor();
523  // extract Xpetra and Thyra based GIDs
524  Teuchos::RCP<const Map> xpsubmap2 = fullDomainMapExtractor->getMap(colleaf->GetIndex(),false);
525  Teuchos::RCP<const Map> tysubmap2 = fullDomainMapExtractor->getMap(colleaf->GetIndex(),true);
526  std::vector<Teuchos::RCP<const Map> > colXpSubMaps (1, xpsubmap2);
527  std::vector<Teuchos::RCP<const Map> > colTySubMaps (1, tysubmap2);
528  // use expert constructor
529  Teuchos::RCP<const MapExtractor> doMapExtractor = Teuchos::rcp(new MapExtractor(colXpSubMaps, colTySubMaps));
530 
532  // build reordered block operator
533  rbmat = Teuchos::rcp(new ReorderedBlockedCrsMatrix(rgMapExtractor,doMapExtractor, 33, rowMgr,bmat));
534  rbmat->setMatrix(0,0,mat);
535  } else {
536  // If leaf node is already wrapped into a blocked matrix do not wrap it again.
537  rbmat = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(mat);
538  TEUCHOS_ASSERT(rbmat != Teuchos::null);
539  }
540  TEUCHOS_ASSERT(mat->getNodeNumEntries() == rbmat->getNodeNumEntries());
541  } else {
542  // create the map extractors
543  // we cannot create block matrix in thyra mode since merged maps might not start with 0 GID
544  Teuchos::RCP<const MapExtractor> rgMapExtractor = Teuchos::null;
545  if(rowSz > 0) {
546  std::vector<Teuchos::RCP<const Map> > rowXpSubMaps (rowSz, Teuchos::null);
547  std::vector<Teuchos::RCP<const Map> > rowTySubMaps (rowSz, Teuchos::null);
548  for(size_t i = 0; i < rowSz; i++) {
549  Teuchos::RCP<const Xpetra::BlockReorderManager> rowSubMgr = rowMgr->GetBlock(Teuchos::as<int>(i));
550  // extract Xpetra and Thyra based merged GIDs
551  rowXpSubMaps[i] = mergeSubBlockMaps(rowSubMgr,bmat,false);
552  rowTySubMaps[i] = mergeSubBlockMaps(rowSubMgr,bmat,true);
553  TEUCHOS_ASSERT(rowXpSubMaps[i].is_null()==false);
554  TEUCHOS_ASSERT(rowTySubMaps[i].is_null()==false);
555  }
556  // use expert constructor
557  rgMapExtractor = Teuchos::rcp(new MapExtractor(rowXpSubMaps, rowTySubMaps));
558  } else {
559  Teuchos::RCP<const Xpetra::BlockReorderLeaf> rowleaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(rowMgr);
560  RCP<const MapExtractor> fullRangeMapExtractor = bmat->getRangeMapExtractor();
561  // extract Xpetra and Thyra based GIDs
562  Teuchos::RCP<const Map> xpsubmap = fullRangeMapExtractor->getMap(rowleaf->GetIndex(),false);
563  Teuchos::RCP<const Map> thysubmap = fullRangeMapExtractor->getMap(rowleaf->GetIndex(),true);
564  std::vector<Teuchos::RCP<const Map> > rowXpSubMaps (1, xpsubmap);
565  std::vector<Teuchos::RCP<const Map> > rowTySubMaps (1, thysubmap);
566  // use expert constructor
567  rgMapExtractor = Teuchos::rcp(new MapExtractor(rowXpSubMaps, rowTySubMaps));
568  }
569 
570  Teuchos::RCP<const MapExtractor> doMapExtractor = Teuchos::null;
571  if(colSz > 0) {
572  std::vector<Teuchos::RCP<const Map> > colXpSubMaps (colSz, Teuchos::null);
573  std::vector<Teuchos::RCP<const Map> > colTySubMaps (colSz, Teuchos::null);
574  for(size_t j = 0; j < colSz; j++) {
575  Teuchos::RCP<const Xpetra::BlockReorderManager> colSubMgr = colMgr->GetBlock(Teuchos::as<int>(j));
576  // extract Xpetra and Thyra based merged GIDs
577  colXpSubMaps[j] = mergeSubBlockMaps(colSubMgr,bmat,false);
578  colTySubMaps[j] = mergeSubBlockMaps(colSubMgr,bmat,true);
579  TEUCHOS_ASSERT(colXpSubMaps[j].is_null()==false);
580  TEUCHOS_ASSERT(colTySubMaps[j].is_null()==false);
581  }
582  // use expert constructor
583  doMapExtractor = Teuchos::rcp(new MapExtractor(colXpSubMaps,colTySubMaps));
584  } else {
585  Teuchos::RCP<const Xpetra::BlockReorderLeaf> colleaf = Teuchos::rcp_dynamic_cast<const Xpetra::BlockReorderLeaf>(colMgr);
586  RCP<const MapExtractor> fullDomainMapExtractor = bmat->getDomainMapExtractor();
587  // extract Xpetra and Thyra based GIDs
588  Teuchos::RCP<const Map> xpsubmap = fullDomainMapExtractor->getMap(colleaf->GetIndex(),false);
589  Teuchos::RCP<const Map> tysubmap = fullDomainMapExtractor->getMap(colleaf->GetIndex(),true);
590  std::vector<Teuchos::RCP<const Map> > colXpSubMaps (1, xpsubmap);
591  std::vector<Teuchos::RCP<const Map> > colTySubMaps (1, tysubmap);
592  // use expert constructor
593  doMapExtractor = Teuchos::rcp(new MapExtractor(colXpSubMaps, colTySubMaps));
594  }
595 
596  // TODO matrix should have both rowMgr and colMgr??
597  rbmat = Teuchos::rcp(new ReorderedBlockedCrsMatrix(rgMapExtractor,doMapExtractor, 33, rowMgr,bmat));
598 
599  size_t cntNNZ = 0;
600 
601  if (rowSz == 0 && colSz > 0) {
602  for(size_t j = 0; j < colSz; j++) {
603  Teuchos::RCP<const Xpetra::BlockReorderManager> colSubMgr = colMgr->GetBlock(Teuchos::as<int>(j));
604  Teuchos::RCP<const Matrix> submat = mergeSubBlocksThyra(rowMgr, colSubMgr, bmat);
605  rbmat->setMatrix(0,j,Teuchos::rcp_const_cast<Matrix>(submat));
606  if(submat != Teuchos::null) cntNNZ += submat->getNodeNumEntries();
607  }
608  } else if (rowSz > 0 && colSz == 0) {
609  for(size_t i = 0; i < rowSz; i++) {
610  Teuchos::RCP<const Xpetra::BlockReorderManager> rowSubMgr = rowMgr->GetBlock(Teuchos::as<int>(i));
611  Teuchos::RCP<const Matrix> submat = mergeSubBlocksThyra(rowSubMgr, colMgr, bmat);
612  rbmat->setMatrix(i,0,Teuchos::rcp_const_cast<Matrix>(submat));
613  if(submat != Teuchos::null) cntNNZ += submat->getNodeNumEntries();
614  }
615  } else {
616  for(size_t i = 0; i < rowSz; i++) {
617  Teuchos::RCP<const Xpetra::BlockReorderManager> rowSubMgr = rowMgr->GetBlock(Teuchos::as<int>(i));
618  for(size_t j = 0; j < colSz; j++) {
619  Teuchos::RCP<const Xpetra::BlockReorderManager> colSubMgr = colMgr->GetBlock(Teuchos::as<int>(j));
620  Teuchos::RCP<const Matrix> submat = mergeSubBlocksThyra(rowSubMgr, colSubMgr, bmat);
621  rbmat->setMatrix(i,j,Teuchos::rcp_const_cast<Matrix>(submat));
622  if(submat != Teuchos::null) cntNNZ += submat->getNodeNumEntries();
623  }
624  }
625  }
626  TEUCHOS_ASSERT(rbmat->getNodeNumEntries() == cntNNZ);
627  }
628 
629  rbmat->fillComplete();
630  return rbmat;
631 }
632 
633 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
635  TEUCHOS_ASSERT(bmat->getRangeMapExtractor()->getThyraMode() == bmat->getDomainMapExtractor()->getThyraMode());
637  if(bmat->getRangeMapExtractor()->getThyraMode() == false) {
638  rbmat = mergeSubBlocks(brm, brm, bmat);
639  } else {
640  rbmat = mergeSubBlocksThyra(brm, brm, bmat);
641  }
642 
643  // TAW, 6/7/2016: rbmat might be Teuchos::null for empty blocks!
644  return rbmat;
645 }
646 
647 } //namespace Xpetra
648 
649 #define XPETRA_REORDEREDBLOCKEDCRSMATRIX_SHORT
650 #endif /* XPETRA_REORDEREDBLOCKEDCRSMATRIX_HPP */
virtual size_t getNodeNumEntries() const =0
Returns the local number of entries in this matrix.
bool is_null(const boost::shared_ptr< T > &p)
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)
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode. ...
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete.
RCP< const Map > getRangeMap() const
Returns the Map associated with the range of this operator.
Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > fullOp_
void swap(RCP< T > &r_ptr)
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
Xpetra utility class for common map-related routines.
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.
Teuchos::RCP< const Xpetra::BlockReorderManager > getBlockReorderManager()
Returns internal BlockReorderManager object.
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.
BlockedCrsMatrix(const Teuchos::RCP< const BlockedMap > &rangeMaps, const Teuchos::RCP< const BlockedMap > &domainMaps, size_t npr, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor.
int GetIndex() const
Get the index that is stored in this block/leaf.
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)
virtual void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)=0
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
virtual size_t GetNumBlocks() const
Returns number of subblocks.
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_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< Matrix > getMatrix(size_t r, size_t c) const
return block (r,c)
virtual const Teuchos::RCP< BlockReorderManager > GetBlock(int blockIndex)
Get a particular block. If there is no block at this index location return a new one.
virtual size_t Cols() const
number of column blocks
std::string viewLabel_t
static const EVerbosityLevel verbLevel_default
const RCP< const Map > getMap(size_t i, bool bThyraMode=false) const
std::string description() const
Return a simple one-line description of this object.
Teuchos::RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > buildReorderedBlockedCrsMatrix(Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > bmat)
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > mergeSubBlockMaps(Teuchos::RCP< const Xpetra::BlockReorderManager > brm)
void setMatrix(size_t r, size_t c, Teuchos::RCP< Matrix > mat)
set matrix block
virtual void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
Concrete implementation of Xpetra::Matrix.
ReorderedBlockedCrsMatrix(Teuchos::RCP< const MapExtractor > &rangeMaps, Teuchos::RCP< const MapExtractor > &domainMaps, size_t npr, Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > bmat, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor.
virtual size_t Rows() const
number of row blocks
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)
#define TEUCHOS_ASSERT(assertion_test)
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
Xpetra-specific matrix class.
virtual void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
Teuchos::RCP< const Xpetra::BlockReorderManager > brm_
bool is_null() const