Xpetra_BlockedCrsMatrix.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_BLOCKEDCRSMATRIX_HPP
47 #define XPETRA_BLOCKEDCRSMATRIX_HPP
48 
49 #include <Kokkos_DefaultNode.hpp>
50 
52 #include <Teuchos_Hashtable.hpp>
53 
54 #include "Xpetra_ConfigDefs.hpp"
55 #include "Xpetra_Exceptions.hpp"
56 
57 #include "Xpetra_MapFactory.hpp"
58 #include "Xpetra_MultiVector.hpp"
61 #include "Xpetra_BlockedVector.hpp"
62 #include "Xpetra_CrsGraph.hpp"
63 #include "Xpetra_CrsMatrix.hpp"
65 
66 #include "Xpetra_MapExtractor.hpp"
67 
68 #include "Xpetra_Matrix.hpp"
69 #include "Xpetra_MatrixFactory.hpp"
70 #include "Xpetra_CrsMatrixWrap.hpp"
71 
72 #ifdef HAVE_XPETRA_THYRA
73 #include <Thyra_ProductVectorSpaceBase.hpp>
74 #include <Thyra_VectorSpaceBase.hpp>
75 #include <Thyra_LinearOpBase.hpp>
76 #include <Thyra_BlockedLinearOpBase.hpp>
77 #include <Thyra_PhysicallyBlockedLinearOpBase.hpp>
78 #include "Xpetra_ThyraUtils.hpp"
79 #endif
80 
81 
86 namespace Xpetra {
87 
88  typedef std::string viewLabel_t;
89 
90  template <class Scalar,
91  class LocalOrdinal,
92  class GlobalOrdinal,
93  class Node>
95  public Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
96  public:
97  typedef Scalar scalar_type;
98  typedef LocalOrdinal local_ordinal_type;
99  typedef GlobalOrdinal global_ordinal_type;
100  typedef Node node_type;
101 
102  private:
103 #undef XPETRA_BLOCKEDCRSMATRIX_SHORT
104 #include "Xpetra_UseShortNames.hpp"
105 
106  public:
107 
109 
110 
112 
119  const Teuchos::RCP<const BlockedMap>& domainMaps,
120  size_t npr, Xpetra::ProfileType pftype = Xpetra::DynamicProfile) {
121  domainmaps_ = Teuchos::rcp(new MapExtractor(domainMaps));
122  rangemaps_ = Teuchos::rcp(new MapExtractor(rangeMaps));
123  bRangeThyraMode_ = rangeMaps->getThyraMode();
124  bDomainThyraMode_ = domainMaps->getThyraMode();
125 
126  blocks_.reserve(Rows()*Cols());
127 
128  // add CrsMatrix objects in row,column order
129  for (size_t r = 0; r < Rows(); ++r)
130  for (size_t c = 0; c < Cols(); ++c)
131  blocks_.push_back(MatrixFactory::Build(getRangeMap(r,bRangeThyraMode_), npr, pftype));
132 
133  // Default view
135  }
136 
138 
148  size_t npr, Xpetra::ProfileType pftype = Xpetra::DynamicProfile)
149  : domainmaps_(domainMaps), rangemaps_(rangeMaps)
150  {
151  bRangeThyraMode_ = rangeMaps->getThyraMode();
152  bDomainThyraMode_ = domainMaps->getThyraMode();
153 
154  blocks_.reserve(Rows()*Cols());
155 
156  // add CrsMatrix objects in row,column order
157  for (size_t r = 0; r < Rows(); ++r)
158  for (size_t c = 0; c < Cols(); ++c)
159  blocks_.push_back(MatrixFactory::Build(getRangeMap(r,bRangeThyraMode_), npr, pftype));
160 
161  // Default view
163  }
164 
165 #ifdef HAVE_XPETRA_THYRA
166 
173  BlockedCrsMatrix(const Teuchos::RCP<const Thyra::BlockedLinearOpBase<Scalar> >& thyraOp, const Teuchos::RCP<const Teuchos::Comm<int> >& comm)
174  : thyraOp_(thyraOp)
175  {
176  // extract information from Thyra blocked operator and rebuilt information
177  const Teuchos::RCP<const Thyra::ProductVectorSpaceBase<Scalar> > productRangeSpace = thyraOp->productRange();
178  const Teuchos::RCP<const Thyra::ProductVectorSpaceBase<Scalar> > productDomainSpace = thyraOp->productDomain();
179 
180  int numRangeBlocks = productRangeSpace->numBlocks();
181  int numDomainBlocks = productDomainSpace->numBlocks();
182 
183  // build range map extractor from Thyra::BlockedLinearOpBase object
184  std::vector<Teuchos::RCP<const Map> > subRangeMaps(numRangeBlocks);
185  for (size_t r=0; r<Teuchos::as<size_t>(numRangeBlocks); ++r) {
186  for (size_t c=0; c<Teuchos::as<size_t>(numDomainBlocks); ++c) {
187  if (thyraOp->blockExists(r,c)) {
188  // we only need at least one block in each block row to extract the range map
189  Teuchos::RCP<const Thyra::LinearOpBase<Scalar> > const_op = thyraOp->getBlock(r,c); // nonConst access is not allowed.
192  subRangeMaps[r] = xop->getRangeMap();
193  break;
194  }
195  }
196  }
197  Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > fullRangeMap = mergeMaps(subRangeMaps);
198 
199  // check whether the underlying Thyra operator uses Thyra-style numbering (default for most applications) or
200  // Xpetra style numbering
201  bool bRangeUseThyraStyleNumbering = false;
202  size_t numAllElements = 0;
203  for(size_t v = 0; v < subRangeMaps.size(); ++v) {
204  numAllElements += subRangeMaps[v]->getGlobalNumElements();
205  }
206  if ( fullRangeMap->getGlobalNumElements() != numAllElements) bRangeUseThyraStyleNumbering = true;
207  rangemaps_ = Teuchos::rcp(new Xpetra::MapExtractor<Scalar,LocalOrdinal,GlobalOrdinal,Node>(fullRangeMap, subRangeMaps, bRangeUseThyraStyleNumbering));
208 
209  // build domain map extractor from Thyra::BlockedLinearOpBase object
210  std::vector<Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > > subDomainMaps(numDomainBlocks);
211  for (size_t c=0; c<Teuchos::as<size_t>(numDomainBlocks); ++c) {
212  for (size_t r=0; r<Teuchos::as<size_t>(numRangeBlocks); ++r) {
213  if (thyraOp->blockExists(r,c)) {
214  // we only need at least one block in each block row to extract the range map
215  Teuchos::RCP<const Thyra::LinearOpBase<Scalar> > const_op = thyraOp->getBlock(r,c); // nonConst access is not allowed.
218  subDomainMaps[c] = xop->getDomainMap();
219  break;
220  }
221  }
222  }
223  Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > fullDomainMap = mergeMaps(subDomainMaps);
224  // plausibility check for numbering style (Xpetra or Thyra)
225  bool bDomainUseThyraStyleNumbering = false;
226  numAllElements = 0;
227  for(size_t v = 0; v < subDomainMaps.size(); ++v) {
228  numAllElements += subDomainMaps[v]->getGlobalNumElements();
229  }
230  if (fullDomainMap->getGlobalNumElements() != numAllElements) bDomainUseThyraStyleNumbering = true;
231  domainmaps_ = Teuchos::rcp(new Xpetra::MapExtractor<Scalar,LocalOrdinal,GlobalOrdinal,Node>(fullDomainMap, subDomainMaps, bDomainUseThyraStyleNumbering));
232 
233  // store numbering mode
234  bRangeThyraMode_ = bRangeUseThyraStyleNumbering;
235  bDomainThyraMode_ = bDomainUseThyraStyleNumbering;
236 
237  // add CrsMatrix objects in row,column order
238  blocks_.reserve(Rows()*Cols());
239  for (size_t r = 0; r < Rows(); ++r) {
240  for (size_t c = 0; c < Cols(); ++c) {
241  if(thyraOp->blockExists(r,c)) {
242  // TODO we do not support nested Thyra operators here!
243  Teuchos::RCP<const Thyra::LinearOpBase<Scalar> > const_op = thyraOp->getBlock(r,c); // nonConst access is not allowed.
244  Teuchos::RCP<Thyra::LinearOpBase<Scalar> > op = Teuchos::rcp_const_cast<Thyra::LinearOpBase<Scalar> >(const_op); // cast away const
249  blocks_.push_back(xwrap);
250  } else {
251  // add empty block
253  }
254  }
255  }
256  // Default view
258  }
259 
260  private:
262 
269  Teuchos::RCP<const Map> mergeMaps(std::vector<Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > > & subMaps) {
270  // TODO merging for Thyra mode is missing (similar to what we do in constructor of MapExtractor
271 
272  // merge submaps to global map
273  std::vector<GlobalOrdinal> gids;
274  for(size_t tt = 0; tt<subMaps.size(); ++tt) {
276 #if 1
277  Teuchos::ArrayView< const GlobalOrdinal > subMapGids = subMap->getNodeElementList();
278  gids.insert(gids.end(), subMapGids.begin(), subMapGids.end());
279 #else
280  size_t myNumElements = subMap->getNodeNumElements();
281  for(LocalOrdinal l = 0; l < Teuchos::as<LocalOrdinal>(myNumElements); ++l) {
282  GlobalOrdinal gid = subMap->getGlobalElement(l);
283  gids.push_back(gid);
284  }
285 #endif
286  }
287 
288  // we have to sort the matrix entries and get rid of the double entries
289  // since we use this to detect Thyra-style numbering or Xpetra-style
290  // numbering. In Thyra-style numbering mode, the Xpetra::MapExtractor builds
291  // the correct row maps.
293  std::sort(gids.begin(), gids.end());
294  gids.erase(std::unique(gids.begin(), gids.end()), gids.end());
295  Teuchos::ArrayView<GO> gidsView(&gids[0], gids.size());
296  Teuchos::RCP<Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > fullMap = Xpetra::MapFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(subMaps[0]->lib(), INVALID, gidsView, subMaps[0]->getIndexBase(), subMaps[0]->getComm());
297  return fullMap;
298  }
299 
300  public:
301 #endif
302 
304  virtual ~BlockedCrsMatrix() {}
305 
307 
308 
310 
311 
313 
338  void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal>& cols, const ArrayView<const Scalar>& vals) {
339  XPETRA_MONITOR("XpetraBlockedCrsMatrix::insertGlobalValues");
340  if (Rows() == 1 && Cols () == 1) {
341  getMatrix(0,0)->insertGlobalValues(globalRow, cols, vals);
342  return;
343  }
344  throw Xpetra::Exceptions::RuntimeError("insertGlobalValues not supported by BlockedCrsMatrix");
345  }
346 
348 
358  void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal>& cols, const ArrayView<const Scalar>& vals) {
359  XPETRA_MONITOR("XpetraBlockedCrsMatrix::insertLocalValues");
360  if (Rows() == 1 && Cols () == 1) {
361  getMatrix(0,0)->insertLocalValues(localRow, cols, vals);
362  return;
363  }
364  throw Xpetra::Exceptions::RuntimeError("insertLocalValues not supported by BlockedCrsMatrix");
365  }
366 
368  XPETRA_MONITOR("XpetraBlockedCrsMatrix::removeEmptyProcessesInPlace");
369  if (Rows() == 1 && Cols () == 1) {
370  getMatrix(0,0)->removeEmptyProcessesInPlace(newMap);
371  return;
372  }
373  throw Xpetra::Exceptions::RuntimeError("removeEmptyProcesses not supported by BlockedCrsMatrix");
374  }
375 
377 
385  void replaceGlobalValues(GlobalOrdinal globalRow,
386  const ArrayView<const GlobalOrdinal> &cols,
387  const ArrayView<const Scalar> &vals) {
388  XPETRA_MONITOR("XpetraBlockedCrsMatrix::replaceGlobalValues");
389  if (Rows() == 1 && Cols () == 1) {
390  getMatrix(0,0)->replaceGlobalValues(globalRow,cols,vals);
391  return;
392  }
393  throw Xpetra::Exceptions::RuntimeError("replaceGlobalValues not supported by BlockedCrsMatrix");
394  }
395 
397 
401  void replaceLocalValues(LocalOrdinal localRow,
402  const ArrayView<const LocalOrdinal> &cols,
403  const ArrayView<const Scalar> &vals) {
404  XPETRA_MONITOR("XpetraBlockedCrsMatrix::replaceLocalValues");
405  if (Rows() == 1 && Cols () == 1) {
406  getMatrix(0,0)->replaceLocalValues(localRow,cols,vals);
407  return;
408  }
409  throw Xpetra::Exceptions::RuntimeError("replaceLocalValues not supported by BlockedCrsMatrix");
410  }
411 
413  virtual void setAllToScalar(const Scalar& alpha) {
414  XPETRA_MONITOR("XpetraBlockedCrsMatrix::setAllToScalar");
415  for (size_t row = 0; row < Rows(); row++) {
416  for (size_t col = 0; col < Cols(); col++) {
417  if (!getMatrix(row,col).is_null()) {
418  getMatrix(row,col)->setAllToScalar(alpha);
419  }
420  }
421  }
422  }
423 
425  void scale(const Scalar& alpha) {
426  XPETRA_MONITOR("XpetraBlockedCrsMatrix::scale");
427  for (size_t row = 0; row < Rows(); row++) {
428  for (size_t col = 0; col < Cols(); col++) {
429  if (!getMatrix(row,col).is_null()) {
430  getMatrix(row,col)->scale(alpha);
431  }
432  }
433  }
434  }
435 
437 
439 
440 
448  void resumeFill(const RCP< ParameterList >& params = null) {
449  XPETRA_MONITOR("XpetraBlockedCrsMatrix::resumeFill");
450  for (size_t row = 0; row < Rows(); row++) {
451  for (size_t col = 0; col < Cols(); col++) {
452  if (!getMatrix(row,col).is_null()) {
453  getMatrix(row,col)->resumeFill(params);
454  }
455  }
456  }
457  }
458 
464  void fillComplete(const RCP<const Map>& domainMap, const RCP<const Map>& rangeMap, const RCP<ParameterList>& params = null) {
465  XPETRA_MONITOR("XpetraBlockedCrsMatrix::fillComplete");
466  if (Rows() == 1 && Cols () == 1) {
467  getMatrix(0,0)->fillComplete(domainMap, rangeMap, params);
468  return;
469  }
470  fillComplete(params);
471  }
472 
487  void fillComplete(const RCP<ParameterList>& params = null) {
488  XPETRA_MONITOR("XpetraBlockedCrsMatrix::fillComplete");
489  TEUCHOS_TEST_FOR_EXCEPTION(rangemaps_==Teuchos::null, Xpetra::Exceptions::RuntimeError,"BlockedCrsMatrix::fillComplete: rangemaps_ is not set. Error.");
490 
491  for (size_t r = 0; r < Rows(); ++r)
492  for (size_t c = 0; c < Cols(); ++c) {
493  if(getMatrix(r,c) != Teuchos::null) {
494  Teuchos::RCP<Matrix> Ablock = getMatrix(r,c);
495 
496  if (Ablock != Teuchos::null && !Ablock->isFillComplete())
497  Ablock->fillComplete(getDomainMap(c, bDomainThyraMode_), getRangeMap(r, bRangeThyraMode_), params);
498  }
499  }
500 
501 #if 0
502  // get full row map
503  RCP<const Map> rangeMap = rangemaps_->getFullMap();
504  fullrowmap_ = MapFactory::Build(rangeMap()->lib(), rangeMap()->getGlobalNumElements(), rangeMap()->getNodeElementList(), rangeMap()->getIndexBase(), rangeMap()->getComm());
505 
506  // build full col map
507  fullcolmap_ = Teuchos::null; // delete old full column map
508 
509  std::vector<GO> colmapentries;
510  for (size_t c = 0; c < Cols(); ++c) {
511  // copy all local column lids of all block rows to colset
512  std::set<GO> colset;
513  for (size_t r = 0; r < Rows(); ++r) {
514  Teuchos::RCP<CrsMatrix> Ablock = getMatrix(r,c);
515 
516  if (Ablock != Teuchos::null) {
517  Teuchos::ArrayView<const GO> colElements = Ablock->getColMap()->getNodeElementList();
518  Teuchos::RCP<const Map> colmap = Ablock->getColMap();
519  copy(colElements.getRawPtr(), colElements.getRawPtr() + colElements.size(), inserter(colset, colset.begin()));
520  }
521  }
522 
523  // remove duplicates (entries which are in column maps of more than one block row)
524  colmapentries.reserve(colmapentries.size() + colset.size());
525  copy(colset.begin(), colset.end(), back_inserter(colmapentries));
526  sort(colmapentries.begin(), colmapentries.end());
527  typename std::vector<GO>::iterator gendLocation;
528  gendLocation = std::unique(colmapentries.begin(), colmapentries.end());
529  colmapentries.erase(gendLocation,colmapentries.end());
530  }
531 
532  // sum up number of local elements
533  size_t numGlobalElements = 0;
534  Teuchos::reduceAll(*(rangeMap->getComm()), Teuchos::REDUCE_SUM, colmapentries.size(), Teuchos::outArg(numGlobalElements));
535 
536  // store global full column map
537  const Teuchos::ArrayView<const GO> aView = Teuchos::ArrayView<const GO>(colmapentries);
538  fullcolmap_ = MapFactory::Build(rangeMap->lib(), numGlobalElements, aView, 0, rangeMap->getComm());
539 #endif
540  }
541 
543 
545 
548  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalNumRows");
549  global_size_t globalNumRows = 0;
550 
551  for (size_t row = 0; row < Rows(); row++)
552  for (size_t col = 0; col < Cols(); col++)
553  if (!getMatrix(row,col).is_null()) {
554  globalNumRows += getMatrix(row,col)->getGlobalNumRows();
555  break; // we need only one non-null matrix in a row
556  }
557 
558  return globalNumRows;
559  }
560 
562 
565  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalNumCols");
566  global_size_t globalNumCols = 0;
567 
568  for (size_t col = 0; col < Cols(); col++)
569  for (size_t row = 0; row < Rows(); row++)
570  if (!getMatrix(row,col).is_null()) {
571  globalNumCols += getMatrix(row,col)->getGlobalNumCols();
572  break; // we need only one non-null matrix in a col
573  }
574 
575  return globalNumCols;
576  }
577 
579  size_t getNodeNumRows() const {
580  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNodeNumRows");
581  global_size_t nodeNumRows = 0;
582 
583  for (size_t row = 0; row < Rows(); ++row)
584  for (size_t col = 0; col < Cols(); col++)
585  if (!getMatrix(row,col).is_null()) {
586  nodeNumRows += getMatrix(row,col)->getNodeNumRows();
587  break; // we need only one non-null matrix in a row
588  }
589 
590  return nodeNumRows;
591  }
592 
595  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalNumEntries");
596  global_size_t globalNumEntries = 0;
597 
598  for (size_t row = 0; row < Rows(); ++row)
599  for (size_t col = 0; col < Cols(); ++col)
600  if (!getMatrix(row,col).is_null())
601  globalNumEntries += getMatrix(row,col)->getGlobalNumEntries();
602 
603  return globalNumEntries;
604  }
605 
607  size_t getNodeNumEntries() const {
608  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNodeNumEntries");
609  global_size_t nodeNumEntries = 0;
610 
611  for (size_t row = 0; row < Rows(); ++row)
612  for (size_t col = 0; col < Cols(); ++col)
613  if (!getMatrix(row,col).is_null())
614  nodeNumEntries += getMatrix(row,col)->getNodeNumEntries();
615 
616  return nodeNumEntries;
617  }
618 
620 
621  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const {
622  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNumEntriesInLocalRow");
623  if (Rows() == 1 && Cols () == 1) {
624  return getMatrix(0,0)->getNumEntriesInLocalRow(localRow);
625  }
626  throw Xpetra::Exceptions::RuntimeError("getNumEntriesInLocalRow not supported by BlockedCrsMatrix");
627  }
628 
630 
633  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalNumDiags");
634  if (Rows() == 1 && Cols () == 1) {
635  return getMatrix(0,0)->getGlobalNumDiags();
636  }
637  throw Xpetra::Exceptions::RuntimeError("getGlobalNumDiags() not supported by BlockedCrsMatrix");
638  }
639 
641 
643  size_t getNodeNumDiags() const {
644  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNodeNumDiags");
645  if (Rows() == 1 && Cols () == 1) {
646  return getMatrix(0,0)->getNodeNumDiags();
647  }
648  throw Xpetra::Exceptions::RuntimeError("getNodeNumDiags() not supported by BlockedCrsMatrix");
649  }
650 
652 
654  size_t getGlobalMaxNumRowEntries() const {
655  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalMaxNumRowEntries");
656  global_size_t globalMaxEntries = 0;
657 
658  for (size_t row = 0; row < Rows(); row++) {
659  global_size_t globalMaxEntriesBlockRows = 0;
660  for (size_t col = 0; col < Cols(); col++) {
661  if (!getMatrix(row,col).is_null()) {
662  globalMaxEntriesBlockRows += getMatrix(row,col)->getGlobalMaxNumRowEntries();
663  }
664  }
665  if(globalMaxEntriesBlockRows > globalMaxEntries)
666  globalMaxEntries = globalMaxEntriesBlockRows;
667  }
668  return globalMaxEntries;
669  }
670 
672 
674  size_t getNodeMaxNumRowEntries() const {
675  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNodeMaxNumRowEntries");
676  size_t localMaxEntries = 0;
677 
678  for (size_t row = 0; row < Rows(); row++) {
679  size_t localMaxEntriesBlockRows = 0;
680  for (size_t col = 0; col < Cols(); col++) {
681  if (!getMatrix(row,col).is_null()) {
682  localMaxEntriesBlockRows += getMatrix(row,col)->getNodeMaxNumRowEntries();
683  }
684  }
685  if(localMaxEntriesBlockRows > localMaxEntries)
686  localMaxEntries = localMaxEntriesBlockRows;
687  }
688  return localMaxEntries;
689  }
690 
692 
695  bool isLocallyIndexed() const {
696  XPETRA_MONITOR("XpetraBlockedCrsMatrix::isLocallyIndexed");
697  for (size_t i = 0; i < blocks_.size(); ++i)
698  if (blocks_[i] != Teuchos::null && !blocks_[i]->isLocallyIndexed())
699  return false;
700  return true;
701  }
702 
704 
707  bool isGloballyIndexed() const {
708  XPETRA_MONITOR("XpetraBlockedCrsMatrix::isGloballyIndexed");
709  for (size_t i = 0; i < blocks_.size(); i++)
710  if (blocks_[i] != Teuchos::null && !blocks_[i]->isGloballyIndexed())
711  return false;
712  return true;
713  }
714 
716  bool isFillComplete() const {
717  XPETRA_MONITOR("XpetraBlockedCrsMatrix::isFillComplete");
718  for (size_t i = 0; i < blocks_.size(); i++)
719  if (blocks_[i] != Teuchos::null && !blocks_[i]->isFillComplete())
720  return false;
721  return true;
722  }
723 
725 
739  virtual void getLocalRowCopy(LocalOrdinal LocalRow,
740  const ArrayView<LocalOrdinal>& Indices,
741  const ArrayView<Scalar>& Values,
742  size_t &NumEntries) const {
743  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalRowCopy");
744  if (Rows() == 1 && Cols () == 1) {
745  getMatrix(0,0)->getLocalRowCopy(LocalRow, Indices, Values, NumEntries);
746  return;
747  }
748  throw Xpetra::Exceptions::RuntimeError("getLocalRowCopy not supported by BlockedCrsMatrix" );
749  }
750 
752 
761  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal>& indices, ArrayView<const Scalar>& values) const {
762  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalRowView");
763  if (Rows() == 1 && Cols () == 1) {
764  getMatrix(0,0)->getGlobalRowView(GlobalRow, indices, values);
765  return;
766  }
767  throw Xpetra::Exceptions::RuntimeError("getGlobalRowView not supported by BlockedCrsMatrix");
768  }
769 
771 
780  void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal>& indices, ArrayView<const Scalar>& values) const {
781  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalRowView");
782  if (Rows() == 1 && Cols () == 1) {
783  getMatrix(0,0)->getLocalRowView(LocalRow, indices, values);
784  return;
785  }
786  throw Xpetra::Exceptions::RuntimeError("getLocalRowView not supported by BlockedCrsMatrix");
787  }
788 
790 
793  void getLocalDiagCopy(Vector& diag) const {
794  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalDiagCopy");
795 
796  //RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
797 
798  Teuchos::RCP<Vector> rcpdiag = Teuchos::rcpFromRef(diag);
799  Teuchos::RCP<BlockedVector> bdiag = Teuchos::rcp_dynamic_cast<BlockedVector>(rcpdiag);
800 
801  // special treatment for 1x1 block matrices
802  // ReorderedBlockedCrsMatrix object encapsulate single blocks in ReorderedBlocks
803  // BlockedVectors have Vector objects as Leaf objects.
804  if(Rows() == 1 && Cols() == 1 && bdiag.is_null() == true) {
806  rm->getLocalDiagCopy(diag);
807  return;
808  }
809 
810  TEUCHOS_TEST_FOR_EXCEPTION(bdiag.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::getLocalDiagCopy: diag must be a Blocked(Multi)Vector.");
811  TEUCHOS_TEST_FOR_EXCEPTION(bdiag->getNumVectors() != 1, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::getLocalDiagCopy: diag must be a Blocked(Multi)Vector with exactly one vector. However, the number of stored vectors is " << bdiag->getNumVectors());
812  TEUCHOS_TEST_FOR_EXCEPTION(bdiag->getBlockedMap()->getNumMaps() != this->Rows(), Xpetra::Exceptions::RuntimeError,
813  "BlockedCrsMatrix::getLocalDiagCopy(): the number of blocks in diag differ from the number of blocks in this operator." );
814  //XPETRA_TEST_FOR_EXCEPTION(bdiag->getMap()->isSameAs(*(getMap())) == false, Xpetra::Exceptions::RuntimeError,
815  // "BlockedCrsMatrix::getLocalDiagCopy(): the map of the vector diag is not compatible with the map of the blocked operator." );
816 
817  for (size_t row = 0; row < Rows(); row++) {
819  if (!rm.is_null()) {
820  Teuchos::RCP<Vector> rv = VectorFactory::Build(bdiag->getBlockedMap()->getMap(row,bdiag->getBlockedMap()->getThyraMode()));
821  rm->getLocalDiagCopy(*rv);
822  bdiag->setMultiVector(row,rv,bdiag->getBlockedMap()->getThyraMode());
823  }
824  }
825  }
826 
828  void leftScale (const Vector& x) {
829  XPETRA_MONITOR("XpetraBlockedCrsMatrix::leftScale");
830 
831  Teuchos::RCP<const Vector> rcpx = Teuchos::rcpFromRef(x);
832  Teuchos::RCP<const BlockedVector> bx = Teuchos::rcp_dynamic_cast<const BlockedVector>(rcpx);
833 
834  //RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
835 
836  // special treatment for 1xn block matrices
837  // ReorderedBlockedCrsMatrix object encapsulate single blocks in ReorderedBlocks
838  // BlockedVectors have Vector objects as Leaf objects.
839  if(Rows() == 1 && bx.is_null() == true) {
840  for (size_t col = 0; col < Cols(); ++col) {
841  Teuchos::RCP<Matrix> rm = getMatrix(0,col);
842  rm->leftScale(x);
843  }
844  return;
845  }
846 
847  TEUCHOS_TEST_FOR_EXCEPTION(bx.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::leftScale: x must be a Blocked(Multi)Vector.");
848  TEUCHOS_TEST_FOR_EXCEPTION(bx->getNumVectors() != 1, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::leftScale: x must be a Blocked(Multi)Vector with exactly one vector. However, the number of stored vectors is " << bx->getNumVectors());
849  TEUCHOS_TEST_FOR_EXCEPTION(bx->getBlockedMap()->getNumMaps() != this->Rows(), Xpetra::Exceptions::RuntimeError,
850  "BlockedCrsMatrix::leftScale(): the number of blocks in diag differ from the number of blocks in this operator." );
851 
852  for (size_t row = 0; row < Rows(); row++) {
853  Teuchos::RCP<const MultiVector> rmv = bx->getMultiVector(row);
854  Teuchos::RCP<const Vector> rscale = rmv->getVector(0);
855  XPETRA_TEST_FOR_EXCEPTION(rscale.is_null()==true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::leftScale: x must be a Vector.");
856  for (size_t col = 0; col < Cols(); ++col) {
857  Teuchos::RCP<Matrix> rm = getMatrix(row,col);
858  if (!rm.is_null()) {
859  rm->leftScale(*rscale);
860  }
861  }
862  }
863  }
864 
866  void rightScale (const Vector& x) {
867  XPETRA_MONITOR("XpetraBlockedCrsMatrix::rightScale");
868 
869  Teuchos::RCP<const Vector> rcpx = Teuchos::rcpFromRef(x);
870  Teuchos::RCP<const BlockedVector> bx = Teuchos::rcp_dynamic_cast<const BlockedVector>(rcpx);
871 
872  //RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
873 
874  // special treatment for nx1 block matrices
875  // ReorderedBlockedCrsMatrix object encapsulate single blocks in ReorderedBlocks
876  // BlockedVectors have Vector objects as Leaf objects.
877  if(Cols() == 1 && bx.is_null() == true) {
878  for (size_t row = 0; row < Rows(); ++row) {
879  Teuchos::RCP<Matrix> rm = getMatrix(row,0);
880  rm->rightScale(x);
881  }
882  return;
883  }
884 
885  TEUCHOS_TEST_FOR_EXCEPTION(bx.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::rightScale: x must be a Blocked(Multi)Vector.");
886  TEUCHOS_TEST_FOR_EXCEPTION(bx->getNumVectors() != 1, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::rightScale: x must be a Blocked(Multi)Vector with exactly one vector. However, the number of stored vectors is " << bx->getNumVectors());
887  TEUCHOS_TEST_FOR_EXCEPTION(bx->getBlockedMap()->getNumMaps() != this->Cols(), Xpetra::Exceptions::RuntimeError,
888  "BlockedCrsMatrix::rightScale(): the number of blocks in diag differ from the number of blocks in this operator." );
889 
890  for (size_t col = 0; col < Cols(); ++col) {
891  Teuchos::RCP<const MultiVector> rmv = bx->getMultiVector(col);
892  Teuchos::RCP<const Vector> rscale = rmv->getVector(0);
893  XPETRA_TEST_FOR_EXCEPTION(rscale.is_null()==true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::leftScale: x must be a Vector.");
894  for (size_t row = 0; row < Rows(); row++) {
895  Teuchos::RCP<Matrix> rm = getMatrix(row,col);
896  if (!rm.is_null()) {
897  rm->rightScale(*rscale);
898  }
899  }
900  }
901  }
902 
903 
906  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getFrobeniusNorm");
908  for (size_t col = 0; col < Cols(); ++col) {
909  for (size_t row = 0; row < Rows(); ++row) {
910  if(getMatrix(row,col)!=Teuchos::null) {
911  typename ScalarTraits<Scalar>::magnitudeType n = getMatrix(row,col)->getFrobeniusNorm();
912  ret += n * n;
913  }
914  }
915  }
917  }
918 
919 
921  virtual bool haveGlobalConstants() const {return true;}
922 
923 
925 
927 
928 
930 
950 
952 
953 
955 
958  virtual void apply(const MultiVector& X, MultiVector& Y,
960  Scalar alpha = ScalarTraits<Scalar>::one(),
961  Scalar beta = ScalarTraits<Scalar>::zero()) const
962  {
963  XPETRA_MONITOR("XpetraBlockedCrsMatrix::apply");
964  //using Teuchos::RCP;
965 
967  "apply() only supports the following modes: NO_TRANS and TRANS." );
968 
969  // check whether input parameters are blocked or not
970  RCP<const MultiVector> refX = rcpFromRef(X);
971  RCP<const BlockedMultiVector> refbX = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(refX);
972  //RCP<MultiVector> tmpY = rcpFromRef(Y);
973  //RCP<BlockedMultiVector> tmpbY = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(tmpY);
974 
975  // TODO get rid of me: adapt MapExtractor
976  bool bBlockedX = (refbX != Teuchos::null) ? true : false;
977 
978  // create (temporary) vectors for output
979  // In the end we call Y.update(alpha, *tmpY, beta). Therefore we need a new vector storing the temporary results
981 
982  //RCP<Teuchos::FancyOStream> out = rcp(new Teuchos::FancyOStream(rcp(&std::cout,false)));
983 
984  SC one = ScalarTraits<SC>::one();
985 
986  if (mode == Teuchos::NO_TRANS) {
987 
988  for (size_t row = 0; row < Rows(); row++) {
989  RCP<MultiVector> Yblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, true);
990  for (size_t col = 0; col < Cols(); col++) {
991 
992  // extract matrix block
993  RCP<Matrix> Ablock = getMatrix(row, col);
994 
995  if (Ablock.is_null())
996  continue;
997 
998  // check whether Ablock is itself a blocked operator
999  // If it is a blocked operator we have to provide Xpetra style GIDs, i.e. we have to transform GIDs
1000  bool bBlockedSubMatrix = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Ablock) == Teuchos::null ? false : true;
1001 
1002  // input/output vectors for local block operation
1003  RCP<const MultiVector> Xblock = Teuchos::null; // subpart of X vector to be applied to subblock of A
1004 
1005  // extract sub part of X using Xpetra or Thyra GIDs
1006  // if submatrix is again blocked, we extract it using Xpetra style gids. If it is a single
1007  // block matrix we use the Thyra or Xpetra style GIDs that are used to store the matrix
1008  if(bBlockedX) Xblock = domainmaps_->ExtractVector(refbX, col, bDomainThyraMode_);
1009  else Xblock = domainmaps_->ExtractVector(refX, col, bBlockedSubMatrix == true ? false : bDomainThyraMode_);
1010 
1011  RCP<MultiVector> tmpYblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, false); // subpart of Y vector containing part of solution of Xblock applied to Ablock
1012  Ablock->apply(*Xblock, *tmpYblock);
1013 
1014  Yblock->update(one, *tmpYblock, one);
1015  }
1016  rangemaps_->InsertVector(Yblock, row, tmpY, bRangeThyraMode_);
1017  }
1018 
1019  } else if (mode == Teuchos::TRANS) {
1020  // TODO: test me!
1021  for (size_t col = 0; col < Cols(); col++) {
1022  RCP<MultiVector> Yblock = domainmaps_->getVector(col, Y.getNumVectors(), bDomainThyraMode_, true);
1023 
1024  for (size_t row = 0; row<Rows(); row++) {
1025  RCP<Matrix> Ablock = getMatrix(row, col);
1026 
1027  if (Ablock.is_null())
1028  continue;
1029 
1030  // check whether Ablock is itself a blocked operator
1031  // If it is a blocked operator we have to provide Xpetra style GIDs, i.e. we have to transform GIDs
1032  bool bBlockedSubMatrix = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Ablock) == Teuchos::null ? false : true;
1033 
1034  RCP<const MultiVector> Xblock = Teuchos::null;
1035 
1036  // extract sub part of X using Xpetra or Thyra GIDs
1037  if(bBlockedX) Xblock = rangemaps_->ExtractVector(refbX, row, bRangeThyraMode_);
1038  else Xblock = rangemaps_->ExtractVector(refX, row, bBlockedSubMatrix == true ? false : bRangeThyraMode_);
1039  RCP<MultiVector> tmpYblock = domainmaps_->getVector(col, Y.getNumVectors(), bDomainThyraMode_, false);
1040  Ablock->apply(*Xblock, *tmpYblock, Teuchos::TRANS);
1041 
1042  Yblock->update(one, *tmpYblock, one);
1043  }
1044  domainmaps_->InsertVector(Yblock, col, tmpY, bDomainThyraMode_);
1045  }
1046  }
1047  Y.update(alpha, *tmpY, beta);
1048  }
1049 
1051  RCP<const Map > getFullDomainMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getFullDomainMap()"); return domainmaps_->getFullMap(); }
1052 
1054  RCP<const BlockedMap > getBlockedDomainMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getBlockedDomainMap()"); return domainmaps_->getBlockedMap(); }
1055 
1057  RCP<const Map > getDomainMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMap()"); return domainmaps_->getMap(); /*domainmaps_->getFullMap();*/ }
1058 
1060  RCP<const Map > getDomainMap(size_t i) const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMap(size_t)"); return domainmaps_->getMap(i, bDomainThyraMode_); }
1061 
1063  RCP<const Map > getDomainMap(size_t i, bool bThyraMode) const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMap(size_t,bool)"); return domainmaps_->getMap(i, bThyraMode); }
1064 
1066  RCP<const Map > getFullRangeMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap()"); return rangemaps_->getFullMap(); }
1067 
1069  RCP<const BlockedMap > getBlockedRangeMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getBlockedRangeMap()"); return rangemaps_->getBlockedMap(); }
1070 
1072  RCP<const Map > getRangeMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap()"); return rangemaps_->getMap(); /*rangemaps_->getFullMap();*/ }
1073 
1075  RCP<const Map > getRangeMap(size_t i) const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap(size_t)"); return rangemaps_->getMap(i, bRangeThyraMode_); }
1076 
1078  RCP<const Map > getRangeMap(size_t i, bool bThyraMode) const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap(size_t,bool)"); return rangemaps_->getMap(i, bThyraMode); }
1079 
1081  RCP<const MapExtractor> getRangeMapExtractor() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMapExtractor()"); return rangemaps_; }
1082 
1084  RCP<const MapExtractor> getDomainMapExtractor() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMapExtractor()"); return domainmaps_; }
1085 
1087 
1089  //{@
1090 
1092 
1098  virtual void bgs_apply(const MultiVector& X, MultiVector& Y, size_t row,
1100  Scalar alpha = ScalarTraits<Scalar>::one(),
1101  Scalar beta = ScalarTraits<Scalar>::zero()) const
1102  {
1103  XPETRA_MONITOR("XpetraBlockedCrsMatrix::bgs_apply");
1104  //using Teuchos::RCP;
1105 
1107  "apply() only supports the following modes: NO_TRANS and TRANS." );
1108 
1109  // check whether input parameters are blocked or not
1110  RCP<const MultiVector> refX = rcpFromRef(X);
1111  RCP<const BlockedMultiVector> refbX = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(refX);
1112  //RCP<MultiVector> tmpY = rcpFromRef(Y);
1113  //RCP<BlockedMultiVector> tmpbY = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(tmpY);
1114 
1115  bool bBlockedX = (refbX != Teuchos::null) ? true : false;
1116 
1117  // create (temporary) vectors for output
1118  // In the end we call Y.update(alpha, *tmpY, beta). Therefore we need a new vector storing the temporary results
1120 
1121  SC one = ScalarTraits<SC>::one();
1122 
1123  if (mode == Teuchos::NO_TRANS) {
1124  RCP<MultiVector> Yblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, true);
1125  for (size_t col = 0; col < Cols(); col++) {
1126 
1127  // extract matrix block
1128  RCP<Matrix> Ablock = getMatrix(row, col);
1129 
1130  if (Ablock.is_null())
1131  continue;
1132 
1133  // check whether Ablock is itself a blocked operator
1134  // If it is a blocked operator we have to provide Xpetra style GIDs, i.e. we have to transform GIDs
1135  bool bBlockedSubMatrix = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Ablock) == Teuchos::null ? false : true;
1136 
1137  // input/output vectors for local block operation
1138  RCP<const MultiVector> Xblock = Teuchos::null; // subpart of X vector to be applied to subblock of A
1139 
1140  // extract sub part of X using Xpetra or Thyra GIDs
1141  // if submatrix is again blocked, we extract it using Xpetra style gids. If it is a single
1142  // block matrix we use the Thyra or Xpetra style GIDs that are used to store the matrix
1143  if(bBlockedX) Xblock = domainmaps_->ExtractVector(refbX, col, bDomainThyraMode_);
1144  else Xblock = domainmaps_->ExtractVector(refX, col, bBlockedSubMatrix == true ? false : bDomainThyraMode_);
1145 
1146  RCP<MultiVector> tmpYblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, false); // subpart of Y vector containing part of solution of Xblock applied to Ablock
1147  Ablock->apply(*Xblock, *tmpYblock);
1148 
1149  Yblock->update(one, *tmpYblock, one);
1150  }
1151  rangemaps_->InsertVector(Yblock, row, tmpY, bRangeThyraMode_);
1152  } else {
1153  TEUCHOS_TEST_FOR_EXCEPTION(true,Xpetra::Exceptions::NotImplemented,"Xpetar::BlockedCrsMatrix::bgs_apply: not implemented for transpose case.");
1154  }
1155  Y.update(alpha, *tmpY, beta);
1156  }
1157 
1158 
1160 
1161 
1163  //{@
1164 
1167  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getMap");
1168  if (Rows() == 1 && Cols () == 1) {
1169  return getMatrix(0,0)->getMap();
1170  }
1171  throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::getMap(): operation not supported.");
1172  }
1173 
1175  void doImport(const Matrix &source, const Import& importer, CombineMode CM) {
1176  XPETRA_MONITOR("XpetraBlockedCrsMatrix::doImport");
1177  if (Rows() == 1 && Cols () == 1) {
1178  getMatrix(0,0)->doImport(source, importer, CM);
1179  return;
1180  }
1181  throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doImport(): operation not supported.");
1182  }
1183 
1185  void doExport(const Matrix& dest, const Import& importer, CombineMode CM) {
1186  XPETRA_MONITOR("XpetraBlockedCrsMatrix::doExport");
1187  if (Rows() == 1 && Cols () == 1) {
1188  getMatrix(0,0)->doExport(dest, importer, CM);
1189  return;
1190  }
1191  throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doExport(): operation not supported.");
1192  }
1193 
1195  void doImport(const Matrix& source, const Export& exporter, CombineMode CM) {
1196  XPETRA_MONITOR("XpetraBlockedCrsMatrix::doImport");
1197  if (Rows() == 1 && Cols () == 1) {
1198  getMatrix(0,0)->doImport(source, exporter, CM);
1199  return;
1200  }
1201  throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doImport(): operation not supported.");
1202  }
1203 
1205  void doExport(const Matrix& dest, const Export& exporter, CombineMode CM) {
1206  XPETRA_MONITOR("XpetraBlockedCrsMatrix::doExport");
1207  if (Rows() == 1 && Cols () == 1) {
1208  getMatrix(0,0)->doExport(dest, exporter, CM);
1209  return;
1210  }
1211  throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doExport(): operation not supported.");
1212  }
1213 
1214  // @}
1215 
1217 
1218 
1220  std::string description() const { return "Xpetra_BlockedCrsMatrix.description()"; }
1221 
1224  out << "Xpetra::BlockedCrsMatrix: " << Rows() << " x " << Cols() << std::endl;
1225 
1226  if (isFillComplete()) {
1227  out << "BlockMatrix is fillComplete" << std::endl;
1228 
1229  /*if(fullrowmap_ != Teuchos::null) {
1230  out << "fullRowMap" << std::endl;
1231  fullrowmap_->describe(out,verbLevel);
1232  } else {
1233  out << "fullRowMap not set. Check whether block matrix is properly fillCompleted!" << std::endl;
1234  }*/
1235 
1236  //out << "fullColMap" << std::endl;
1237  //fullcolmap_->describe(out,verbLevel);
1238 
1239  } else {
1240  out << "BlockMatrix is NOT fillComplete" << std::endl;
1241  }
1242 
1243  for (size_t r = 0; r < Rows(); ++r)
1244  for (size_t c = 0; c < Cols(); ++c) {
1245  if(getMatrix(r,c)!=Teuchos::null) {
1246  out << "Block(" << r << "," << c << ")" << std::endl;
1247  getMatrix(r,c)->describe(out,verbLevel);
1248  } else out << "Block(" << r << "," << c << ") = null" << std::endl;
1249  }
1250  }
1251 
1254  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getCrsGraph");
1255  if (Rows() == 1 && Cols () == 1) {
1256  return getMatrix(0,0)->getCrsGraph();
1257  }
1258  throw Xpetra::Exceptions::RuntimeError("getCrsGraph() not supported by BlockedCrsMatrix");
1259  }
1260 
1262 
1264 
1265 
1267  virtual size_t Rows() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::Rows"); return rangemaps_->NumMaps(); }
1268 
1270  virtual size_t Cols() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::Cols"); return domainmaps_->NumMaps(); }
1271 
1274  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getCrsMatrix");
1275  TEUCHOS_TEST_FOR_EXCEPTION(Rows()!=1, std::out_of_range, "Can only unwrap a 1x1 blocked matrix. The matrix has " << Rows() << " block rows, though.");
1276  TEUCHOS_TEST_FOR_EXCEPTION(Cols()!=1, std::out_of_range, "Can only unwrap a 1x1 blocked matrix. The matrix has " << Cols() << " block columns, though.");
1277 
1278  RCP<Matrix> mat = getMatrix(0,0);
1279  RCP<BlockedCrsMatrix> bmat = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(mat);
1280  if (bmat == Teuchos::null) return mat;
1281  return bmat->getCrsMatrix();
1282  }
1283 
1286  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getInnermostCrsMatrix");
1287  size_t row = Rows()+1, col = Cols()+1;
1288  for (size_t r = 0; r < Rows(); ++r)
1289  for(size_t c = 0; c < Cols(); ++c)
1290  if (getMatrix(r,c) != Teuchos::null) {
1291  row = r;
1292  col = c;
1293  break;
1294  }
1295  TEUCHOS_TEST_FOR_EXCEPTION(row == Rows()+1 || col == Cols()+1, Xpetra::Exceptions::Incompatible, "Xpetra::BlockedCrsMatrix::getInnermostCrsMatrix: Could not find a non-zero sub-block in blocked operator.")
1296  RCP<Matrix> mm = getMatrix(row,col);
1297  RCP<BlockedCrsMatrix> bmat = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(mm);
1298  if (bmat == Teuchos::null) return mm;
1299  return bmat->getInnermostCrsMatrix();
1300  }
1301 
1303  Teuchos::RCP<Matrix> getMatrix(size_t r, size_t c) const {
1304  XPETRA_MONITOR("XpetraBlockedCrsMatrix::getMatrix");
1305  TEUCHOS_TEST_FOR_EXCEPTION(r > Rows(), std::out_of_range, "Error, r = " << Rows() << " is too big");
1306  TEUCHOS_TEST_FOR_EXCEPTION(c > Cols(), std::out_of_range, "Error, c = " << Cols() << " is too big");
1307 
1308  // transfer strided/blocked map information
1309  if (blocks_[r*Cols()+c] != Teuchos::null &&
1310  blocks_[r*Cols()+c]->IsView("stridedMaps") == false)
1311  blocks_[r*Cols()+c]->CreateView("stridedMaps", getRangeMap(r,bRangeThyraMode_), getDomainMap(c,bDomainThyraMode_));
1312  return blocks_[r*Cols()+c];
1313  }
1314 
1316  //void setMatrix(size_t r, size_t c, Teuchos::RCP<CrsMatrix> mat) {
1317  void setMatrix(size_t r, size_t c, Teuchos::RCP<Matrix> mat) {
1318  XPETRA_MONITOR("XpetraBlockedCrsMatrix::setMatrix");
1319  // TODO: if filled -> return error
1320 
1321  TEUCHOS_TEST_FOR_EXCEPTION(r > Rows(), std::out_of_range, "Error, r = " << Rows() << " is too big");
1322  TEUCHOS_TEST_FOR_EXCEPTION(c > Cols(), std::out_of_range, "Error, c = " << Cols() << " is too big");
1323 
1324  // set matrix
1325  blocks_[r*Cols() + c] = mat;
1326  }
1327 
1329  // NOTE: This is a rather expensive operation, since all blocks are copied
1330  // into a new big CrsMatrix
1332  XPETRA_MONITOR("XpetraBlockedCrsMatrix::Merge");
1333  using Teuchos::RCP;
1334  using Teuchos::rcp_dynamic_cast;
1335  Scalar one = ScalarTraits<SC>::one();
1336 
1338  "BlockedCrsMatrix::Merge: only implemented for Xpetra-style or Thyra-style numbering. No mixup allowed!" );
1339 
1341  "BlockedCrsMatrix::Merge: BlockMatrix must be fill-completed." );
1342 
1343  RCP<Matrix> sparse = MatrixFactory::Build(getRangeMapExtractor()->getFullMap(), 33);
1344 
1345  if(bRangeThyraMode_ == false) {
1346  // Xpetra mode
1347  for (size_t i = 0; i < Rows(); i++) {
1348  for (size_t j = 0; j < Cols(); j++) {
1349  if (getMatrix(i,j) != Teuchos::null) {
1350  RCP<const Matrix> mat = getMatrix(i,j);
1351 
1352  // recursively call Merge routine
1353  RCP<const BlockedCrsMatrix> bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1354  if (bMat != Teuchos::null) mat = bMat->Merge();
1355 
1356  bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1358  "BlockedCrsMatrix::Merge: Merging of blocked sub-operators failed?!" );
1359 
1360  // jump over empty blocks
1361  if(mat->getNodeNumEntries() == 0) continue;
1362 
1363  this->Add(*mat, one, *sparse, one);
1364  }
1365  }
1366  }
1367  } else {
1368  // Thyra mode
1369  for (size_t i = 0; i < Rows(); i++) {
1370  for (size_t j = 0; j < Cols(); j++) {
1371  if (getMatrix(i,j) != Teuchos::null) {
1372  RCP<const Matrix> mat = getMatrix(i,j);
1373  // recursively call Merge routine
1374  RCP<const BlockedCrsMatrix> bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1375  if (bMat != Teuchos::null) mat = bMat->Merge();
1376 
1377  bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1379  "BlockedCrsMatrix::Merge: Merging of blocked sub-operators failed?!" );
1380 
1381  // check whether we have a CrsMatrix block (no blocked operator)
1382  RCP<const CrsMatrixWrap> crsMat = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(mat);
1383  TEUCHOS_ASSERT(crsMat != Teuchos::null);
1384 
1385  // these are the thyra style maps of the matrix
1386  RCP<const Map> trowMap = mat->getRowMap();
1387  RCP<const Map> tcolMap = mat->getColMap();
1388  RCP<const Map> tdomMap = mat->getDomainMap();
1389 
1390  // get Xpetra maps
1391  RCP<const Map> xrowMap = getRangeMapExtractor()->getMap(i,false);
1392  RCP<const Map> xdomMap = getDomainMapExtractor()->getMap(j,false);
1393 
1394  // generate column map with Xpetra GIDs
1395  // We have to do this separately for each block since the column
1396  // map of each block might be different in the same block column
1398  *tcolMap,
1399  *tdomMap,
1400  *xdomMap);
1401 
1402  // jump over empty blocks
1403  if(mat->getNodeNumEntries() == 0) continue;
1404 
1405  size_t maxNumEntries = mat->getNodeMaxNumRowEntries();
1406 
1407  size_t numEntries;
1408  Array<GO> inds (maxNumEntries);
1409  Array<GO> inds2(maxNumEntries);
1410  Array<SC> vals (maxNumEntries);
1411 
1412  // loop over all rows and add entries
1413  for (size_t k = 0; k < mat->getNodeNumRows(); k++) {
1414  GlobalOrdinal rowTGID = trowMap->getGlobalElement(k);
1415  crsMat->getCrsMatrix()->getGlobalRowCopy(rowTGID, inds(), vals(), numEntries);
1416 
1417  // create new indices array
1418  for (size_t l = 0; l < numEntries; ++l) {
1419  LocalOrdinal lid = tcolMap->getLocalElement(inds[l]);
1420  inds2[l] = xcolMap->getGlobalElement(lid);
1421  }
1422 
1423  GlobalOrdinal rowXGID = xrowMap->getGlobalElement(k);
1424  sparse->insertGlobalValues(
1425  rowXGID, inds2(0, numEntries),
1426  vals(0, numEntries));
1427  }
1428  }
1429  }
1430  }
1431  }
1432 
1433  sparse->fillComplete(getDomainMap(), getRangeMap());
1434 
1436  "BlockedCrsMatrix::Merge: Local number of entries of merged matrix does not coincide with local number of entries of blocked operator." );
1437 
1439  "BlockedCrsMatrix::Merge: Global number of entries of merged matrix does not coincide with global number of entries of blocked operator." );
1440 
1441  return sparse;
1442  }
1444 
1445 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
1446  typedef typename CrsMatrix::local_matrix_type local_matrix_type;
1447 
1449  local_matrix_type getLocalMatrix () const {
1450  if (Rows() == 1 && Cols () == 1) {
1451  return getMatrix(0,0)->getLocalMatrix();
1452  }
1453  throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::getLocalMatrix(): operation not supported.");
1454  }
1455 #endif
1456 
1457 #ifdef HAVE_XPETRA_THYRA
1460  Xpetra::ThyraUtils<Scalar,LO,GO,Node>::toThyra(Teuchos::rcpFromRef(*this));
1462 
1464  Teuchos::rcp_dynamic_cast<Thyra::BlockedLinearOpBase<Scalar> >(thOp);
1466  return thbOp;
1467  }
1468 #endif
1469 
1470  private:
1471 
1474 
1476 
1486  void Add(const Matrix& A, const Scalar scalarA, Matrix& B, const Scalar scalarB) const {
1487  XPETRA_MONITOR("XpetraBlockedCrsMatrix::Add");
1489  "Matrix A is not completed");
1490  using Teuchos::Array;
1491  using Teuchos::ArrayView;
1492 
1493  B.scale(scalarB);
1494 
1495  Scalar one = ScalarTraits<SC>::one();
1496  Scalar zero = ScalarTraits<SC>::zero();
1497 
1498  if (scalarA == zero)
1499  return;
1500 
1501  Teuchos::RCP<const Matrix> rcpA = Teuchos::rcpFromRef(A);
1502  Teuchos::RCP<const CrsMatrixWrap> rcpAwrap = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(rcpA);
1503  TEUCHOS_TEST_FOR_EXCEPTION(rcpAwrap == Teuchos::null, Xpetra::Exceptions::BadCast,
1504  "BlockedCrsMatrix::Add: matrix A must be of type CrsMatrixWrap.");
1505  Teuchos::RCP<const CrsMatrix> crsA = rcpAwrap->getCrsMatrix();
1506 
1507  size_t maxNumEntries = crsA->getNodeMaxNumRowEntries();
1508 
1509  size_t numEntries;
1510  Array<GO> inds(maxNumEntries);
1511  Array<SC> vals(maxNumEntries);
1512 
1513  RCP<const Map> rowMap = crsA->getRowMap();
1514  RCP<const Map> colMap = crsA->getColMap();
1515 
1516  ArrayView<const GO> rowGIDs = crsA->getRowMap()->getNodeElementList();
1517  for (size_t i = 0; i < crsA->getNodeNumRows(); i++) {
1518  GO row = rowGIDs[i];
1519  crsA->getGlobalRowCopy(row, inds(), vals(), numEntries);
1520 
1521  if (scalarA != one)
1522  for (size_t j = 0; j < numEntries; ++j)
1523  vals[j] *= scalarA;
1524 
1525  B.insertGlobalValues(row, inds(0, numEntries), vals(0, numEntries)); // insert should be ok, since blocks in BlockedCrsOpeartor do not overlap!
1526  }
1527  }
1528 
1530 
1531  // Default view is created after fillComplete()
1532  // Because ColMap might not be available before fillComplete().
1534  XPETRA_MONITOR("XpetraBlockedCrsMatrix::CreateDefaultView");
1535 
1536  // Create default view
1537  this->defaultViewLabel_ = "point";
1539 
1540  // Set current view
1541  this->currentViewLabel_ = this->GetDefaultViewLabel();
1542  }
1543 
1544  private:
1545  Teuchos::RCP<const MapExtractor> domainmaps_; // full domain map together with all partial domain maps
1546  Teuchos::RCP<const MapExtractor> rangemaps_; // full range map together with all partial domain maps
1547 
1548  std::vector<Teuchos::RCP<Matrix> > blocks_; // row major matrix block storage
1549 #ifdef HAVE_XPETRA_THYRA
1551 #endif
1554 
1555 };
1556 
1557 } //namespace Xpetra
1558 
1559 #define XPETRA_BLOCKEDCRSMATRIX_SHORT
1560 #endif /* XPETRA_BLOCKEDCRSMATRIX_HPP */
Teuchos::RCP< Matrix > Merge() const
merge BlockedCrsMatrix blocks in a CrsMatrix
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
virtual size_t getNodeMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on this node.
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const =0
Get this Map&#39;s Comm object.
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
static Teuchos::RCP< Map< LocalOrdinal, GlobalOrdinal, Node > > Build(UnderlyingLib lib, global_size_t numGlobalElements, GlobalOrdinal indexBase, const Teuchos::RCP< const Teuchos::Comm< int > > &comm, LocalGlobal lg=Xpetra::GloballyDistributed, const Teuchos::RCP< Node > &node=defaultArgNode())
Map constructor with Xpetra-defined contiguous uniform distribution.
const Teuchos::RCP< const Map > getMap() const
Implements DistObject interface.
bool isLocallyIndexed() const
If matrix indices of all matrix blocks are in the local range, this function returns true...
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const =0
Returns the Map that describes the column distribution in this matrix.
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode. ...
RCP< const Map > getRangeMap(size_t i, bool bThyraMode) const
Returns the Map associated with the i&#39;th block range of this operator.
RCP< const BlockedMap > getBlockedDomainMap() const
Returns the BlockedMap associated with the domain of this operator.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete.
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
GlobalOrdinal GO
RCP< const Map > getDomainMap(size_t i) const
Returns the Map associated with the i&#39;th block domain of this operator.
global_size_t getGlobalNumDiags() const
Returns the number of global diagonal entries, based on global row/column index comparisons.
bool is_null(const std::shared_ptr< T > &p)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const MapExtractor > rangemaps_
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const =0
Returns the Map that describes the row distribution in this matrix.
Xpetra namespace
void doImport(const Matrix &source, const Import &importer, CombineMode CM)
Import.
Exception throws to report errors in the internal logical of the program.
size_type size() const
virtual ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
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.
virtual LocalOrdinal getLocalElement(GlobalOrdinal globalIndex) const =0
The local index corresponding to the given global index.
virtual ~BlockedCrsMatrix()
Destructor.
Teuchos::RCP< Matrix > getCrsMatrix() const
return unwrap 1x1 blocked operators
void doExport(const Matrix &dest, const Import &importer, CombineMode CM)
Export.
virtual bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
void resumeFill(const RCP< ParameterList > &params=null)
BlockedCrsMatrix(const Teuchos::RCP< const BlockedMap > &rangeMaps, const Teuchos::RCP< const BlockedMap > &domainMaps, size_t npr, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor.
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.
size_t getNodeNumDiags() const
Returns the number of local diagonal entries, based on global row/column index comparisons.
void rightScale(const Vector &x)
Right scale matrix using the given vector entries.
RCP< const BlockedMap > getBlockedRangeMap() const
Returns the BlockedMap associated with the range of this operator.
virtual void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const =0
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
bool is_null() const
RCP< const Map > getFullRangeMap() const
Returns the Map associated with the full range of this operator.
RCP< const MapExtractor > getRangeMapExtractor() const
Returns map extractor class for range map.
Exception indicating invalid cast attempted.
iterator begin() const
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
viewLabel_t defaultViewLabel_
virtual void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Insert matrix entries, using global IDs.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void doExport(const Matrix &dest, const Export &exporter, CombineMode CM)
Export (using an Importer).
void leftScale(const Vector &x)
Left scale matrix using the given vector entries.
Teuchos::RCP< Matrix > getInnermostCrsMatrix()
helper routine recursively returns the first inner-most non-null matrix block from a (nested) blocked...
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
RCP< const Map > getDomainMap(size_t i, bool bThyraMode) const
Returns the Map associated with the i&#39;th block domain of this operator.
virtual bool isFillComplete() const =0
Returns true if fillComplete() has been called and the matrix is in compute mode. ...
const viewLabel_t & GetDefaultViewLabel() const
void getLocalDiagCopy(Vector &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
void Add(const Matrix &A, const Scalar scalarA, Matrix &B, const Scalar scalarB) const
Add a Xpetra::CrsMatrix to another: B = B*scalarB + A*scalarA.
static Teuchos::RCP< Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > transformThyra2XpetraGIDs(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &input, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &nonOvlInput, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &nonOvlReferenceInput)
replace set of global ids by new global ids
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
void CreateView(viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap)
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const =0
The Map describing the parallel distribution of this object.
viewLabel_t currentViewLabel_
TEUCHOS_DEPRECATED void reduceAll(const Comm< Ordinal > &comm, const EReductionType reductType, const Packet &send, Packet *globalReduct)
void doImport(const Matrix &source, const Export &exporter, CombineMode CM)
Import (using an Exporter).
virtual GlobalOrdinal getGlobalElement(LocalOrdinal localIndex) const =0
The global index corresponding to the given local index.
Exception throws when you call an unimplemented method of Xpetra.
RCP< const Map > getRangeMap(size_t i) const
Returns the Map associated with the i&#39;th block range of this operator.
std::string viewLabel_t
size_t global_size_t
Global size_t object.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
static RCP< Vector > Build(const Teuchos::RCP< const Map > &map, bool zeroOut=true)
Constructor specifying the number of non-zeros for all rows.
static const EVerbosityLevel verbLevel_default
std::vector< Teuchos::RCP< Matrix > > blocks_
virtual void bgs_apply(const MultiVector &X, MultiVector &Y, size_t row, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Special multiplication routine (for BGS smoother)
static magnitudeType magnitude(T a)
virtual void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
T * getRawPtr() const
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
RCP< const Map > getFullDomainMap() const
Returns the Map associated with the full domain of this operator.
#define XPETRA_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void setMatrix(size_t r, size_t c, Teuchos::RCP< Matrix > mat)
set matrix block
bool bDomainThyraMode_
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub...
BlockedCrsMatrix(Teuchos::RCP< const MapExtractor > &rangeMaps, Teuchos::RCP< const MapExtractor > &domainMaps, size_t npr, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor.
Exception throws to report incompatible objects (like maps).
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
RCP< const Map > getRangeMap() const
Returns the Map associated with the range of this operator.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
Concrete implementation of Xpetra::Matrix.
Teuchos::RCP< Matrix > getMatrix(size_t r, size_t c) const
return block (r,c)
std::string description() const
Return a simple one-line description of this object.
static Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Constructor specifying the number of non-zeros for all rows.
bool IsView(const viewLabel_t viewLabel) const
CombineMode
Xpetra::Combine Mode enumerable type.
global_size_t getGlobalNumRows() const
Returns the number of global rows.
#define XPETRA_MONITOR(funcName)
RCP< const MapExtractor > getDomainMapExtractor() const
Returns map extractor for domain map.
static RCP< Matrix > Build(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor specifying the number of non-zeros for all rows.
virtual void scale(const Scalar &alpha)=0
Scale the current values of a matrix, this = alpha*this.
void fillComplete(const RCP< ParameterList > &params=null)
Signal that data entry is complete.
#define TEUCHOS_ASSERT(assertion_test)
virtual size_t Rows() const
number of row blocks
iterator end() const
virtual size_t Cols() const
number of column blocks
virtual size_t getNumVectors() const =0
Number of columns in the multivector.
RCP< const Map > getDomainMap() const
Returns the Map associated with the domain of this operator.
virtual size_t getNodeNumRows() const =0
Returns the number of matrix rows owned on the calling node.
virtual UnderlyingLib lib() const =0
Get the library used by this object (Tpetra or Epetra?)
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
virtual Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const =0
Return a Vector which is a const view of column j.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
Xpetra-specific matrix class.
Teuchos::RCP< const MapExtractor > domainmaps_
bool bRangeThyraMode_
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub...
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.