Xpetra Version of the Day
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
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"
59#include "Xpetra_BlockedMultiVector.hpp"
60#include "Xpetra_MultiVectorFactory.hpp"
61#include "Xpetra_BlockedVector.hpp"
62#include "Xpetra_CrsGraph.hpp"
63#include "Xpetra_CrsMatrix.hpp"
65
66#include "Xpetra_MapExtractor.hpp"
68
69#include "Xpetra_Matrix.hpp"
71#include "Xpetra_CrsMatrixWrap.hpp"
72
73#ifdef HAVE_XPETRA_THYRA
74#include <Thyra_ProductVectorSpaceBase.hpp>
75#include <Thyra_VectorSpaceBase.hpp>
76#include <Thyra_LinearOpBase.hpp>
77#include <Thyra_BlockedLinearOpBase.hpp>
78#include <Thyra_PhysicallyBlockedLinearOpBase.hpp>
79#include "Xpetra_ThyraUtils.hpp"
80#endif
81
82#include "Xpetra_VectorFactory.hpp"
83
84
89namespace Xpetra {
90
91 typedef std::string viewLabel_t;
92
93 template <class Scalar,
94 class LocalOrdinal,
95 class GlobalOrdinal,
98 public Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
99 public:
100 typedef Scalar scalar_type;
101 typedef LocalOrdinal local_ordinal_type;
102 typedef GlobalOrdinal global_ordinal_type;
103 typedef Node node_type;
104
105 private:
106#undef XPETRA_BLOCKEDCRSMATRIX_SHORT
108
109 public:
110
112
113
115
121 const Teuchos::RCP<const BlockedMap>& domainMaps,
122 size_t numEntriesPerRow)
123 : is_diagonal_(true)
124 {
127 bRangeThyraMode_ = rangeMaps->getThyraMode();
128 bDomainThyraMode_ = domainMaps->getThyraMode();
129
130 blocks_.reserve(Rows()*Cols());
131
132 // add CrsMatrix objects in row,column order
133 for (size_t r = 0; r < Rows(); ++r)
134 for (size_t c = 0; c < Cols(); ++c)
135 blocks_.push_back(MatrixFactory::Build(getRangeMap(r, bRangeThyraMode_), numEntriesPerRow));
136
137 // Default view
139 }
140
142
150 Teuchos::RCP<const MapExtractor>& domainMapExtractor,
151 size_t numEntriesPerRow)
152 : is_diagonal_(true), domainmaps_(domainMapExtractor), rangemaps_(rangeMapExtractor)
153 {
154 bRangeThyraMode_ = rangeMapExtractor->getThyraMode();
155 bDomainThyraMode_ = domainMapExtractor->getThyraMode();
156
157 blocks_.reserve(Rows()*Cols());
158
159 // add CrsMatrix objects in row,column order
160 for (size_t r = 0; r < Rows(); ++r)
161 for (size_t c = 0; c < Cols(); ++c)
162 blocks_.push_back(MatrixFactory::Build(getRangeMap(r, bRangeThyraMode_), numEntriesPerRow));
163
164 // Default view
166 }
167
168#ifdef HAVE_XPETRA_THYRA
170
175 BlockedCrsMatrix(const Teuchos::RCP<const Thyra::BlockedLinearOpBase<Scalar> >& thyraOp, const Teuchos::RCP<const Teuchos::Comm<int> >& /* comm */)
176 : is_diagonal_(true), thyraOp_(thyraOp)
177 {
178 // extract information from Thyra blocked operator and rebuilt information
179 const Teuchos::RCP<const Thyra::ProductVectorSpaceBase<Scalar> > productRangeSpace = thyraOp->productRange();
180 const Teuchos::RCP<const Thyra::ProductVectorSpaceBase<Scalar> > productDomainSpace = thyraOp->productDomain();
181
182 int numRangeBlocks = productRangeSpace->numBlocks();
183 int numDomainBlocks = productDomainSpace->numBlocks();
184
185 // build range map extractor from Thyra::BlockedLinearOpBase object
186 std::vector<Teuchos::RCP<const Map> > subRangeMaps(numRangeBlocks);
187 for (size_t r=0; r<Teuchos::as<size_t>(numRangeBlocks); ++r) {
188 for (size_t c=0; c<Teuchos::as<size_t>(numDomainBlocks); ++c) {
189 if (thyraOp->blockExists(r,c)) {
190 // we only need at least one block in each block row to extract the range map
191 Teuchos::RCP<const Thyra::LinearOpBase<Scalar> > const_op = thyraOp->getBlock(r,c); // nonConst access is not allowed.
194 subRangeMaps[r] = xop->getRangeMap();
195 if(r!=c) is_diagonal_ = false;
196 break;
197 }
198 }
199 }
200 Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > fullRangeMap = mergeMaps(subRangeMaps);
201
202 // check whether the underlying Thyra operator uses Thyra-style numbering (default for most applications) or
203 // Xpetra style numbering
204 bool bRangeUseThyraStyleNumbering = false;
205 size_t numAllElements = 0;
206 for(size_t v = 0; v < subRangeMaps.size(); ++v) {
207 numAllElements += subRangeMaps[v]->getGlobalNumElements();
208 }
209 if ( fullRangeMap->getGlobalNumElements() != numAllElements) bRangeUseThyraStyleNumbering = true;
210 rangemaps_ = Teuchos::rcp(new Xpetra::MapExtractor<Scalar,LocalOrdinal,GlobalOrdinal,Node>(fullRangeMap, subRangeMaps, bRangeUseThyraStyleNumbering));
211
212 // build domain map extractor from Thyra::BlockedLinearOpBase object
213 std::vector<Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > > subDomainMaps(numDomainBlocks);
214 for (size_t c=0; c<Teuchos::as<size_t>(numDomainBlocks); ++c) {
215 for (size_t r=0; r<Teuchos::as<size_t>(numRangeBlocks); ++r) {
216 if (thyraOp->blockExists(r,c)) {
217 // we only need at least one block in each block row to extract the range map
218 Teuchos::RCP<const Thyra::LinearOpBase<Scalar> > const_op = thyraOp->getBlock(r,c); // nonConst access is not allowed.
221 subDomainMaps[c] = xop->getDomainMap();
222 break;
223 }
224 }
225 }
226 Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > fullDomainMap = mergeMaps(subDomainMaps);
227 // plausibility check for numbering style (Xpetra or Thyra)
228 bool bDomainUseThyraStyleNumbering = false;
229 numAllElements = 0;
230 for(size_t v = 0; v < subDomainMaps.size(); ++v) {
231 numAllElements += subDomainMaps[v]->getGlobalNumElements();
232 }
233 if (fullDomainMap->getGlobalNumElements() != numAllElements) bDomainUseThyraStyleNumbering = true;
234 domainmaps_ = Teuchos::rcp(new Xpetra::MapExtractor<Scalar,LocalOrdinal,GlobalOrdinal,Node>(fullDomainMap, subDomainMaps, bDomainUseThyraStyleNumbering));
235
236 // store numbering mode
237 bRangeThyraMode_ = bRangeUseThyraStyleNumbering;
238 bDomainThyraMode_ = bDomainUseThyraStyleNumbering;
239
240 // add CrsMatrix objects in row,column order
241 blocks_.reserve(Rows()*Cols());
242 for (size_t r = 0; r < Rows(); ++r) {
243 for (size_t c = 0; c < Cols(); ++c) {
244 if(thyraOp->blockExists(r,c)) {
245 // TODO we do not support nested Thyra operators here!
246 Teuchos::RCP<const Thyra::LinearOpBase<Scalar> > const_op = thyraOp->getBlock(r,c); // nonConst access is not allowed.
247 Teuchos::RCP<Thyra::LinearOpBase<Scalar> > op = Teuchos::rcp_const_cast<Thyra::LinearOpBase<Scalar> >(const_op); // cast away const
252 blocks_.push_back(xwrap);
253 } else {
254 // add empty block
256 }
257 }
258 }
259 // Default view
261 }
262
263 private:
265
273 // TODO merging for Thyra mode is missing (similar to what we do in constructor of MapExtractor
274
275 // merge submaps to global map
276 std::vector<GlobalOrdinal> gids;
277 for(size_t tt = 0; tt<subMaps.size(); ++tt) {
279#if 1
280 Teuchos::ArrayView< const GlobalOrdinal > subMapGids = subMap->getNodeElementList();
281 gids.insert(gids.end(), subMapGids.begin(), subMapGids.end());
282#else
283 size_t myNumElements = subMap->getNodeNumElements();
284 for(LocalOrdinal l = 0; l < Teuchos::as<LocalOrdinal>(myNumElements); ++l) {
285 GlobalOrdinal gid = subMap->getGlobalElement(l);
286 gids.push_back(gid);
287 }
288#endif
289 }
290
291 // we have to sort the matrix entries and get rid of the double entries
292 // since we use this to detect Thyra-style numbering or Xpetra-style
293 // numbering. In Thyra-style numbering mode, the Xpetra::MapExtractor builds
294 // the correct row maps.
296 std::sort(gids.begin(), gids.end());
297 gids.erase(std::unique(gids.begin(), gids.end()), gids.end());
298 Teuchos::ArrayView<GO> gidsView(&gids[0], gids.size());
299 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());
300 return fullMap;
301 }
302
303 public:
304#endif
305
307 virtual ~BlockedCrsMatrix() {}
308
310
311
313
314
316
341 void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal>& cols, const ArrayView<const Scalar>& vals) {
342 XPETRA_MONITOR("XpetraBlockedCrsMatrix::insertGlobalValues");
343 if (Rows() == 1 && Cols () == 1) {
344 getMatrix(0,0)->insertGlobalValues(globalRow, cols, vals);
345 return;
346 }
347 throw Xpetra::Exceptions::RuntimeError("insertGlobalValues not supported by BlockedCrsMatrix");
348 }
349
351
361 void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal>& cols, const ArrayView<const Scalar>& vals) {
362 XPETRA_MONITOR("XpetraBlockedCrsMatrix::insertLocalValues");
363 if (Rows() == 1 && Cols () == 1) {
364 getMatrix(0,0)->insertLocalValues(localRow, cols, vals);
365 return;
366 }
367 throw Xpetra::Exceptions::RuntimeError("insertLocalValues not supported by BlockedCrsMatrix");
368 }
369
371 XPETRA_MONITOR("XpetraBlockedCrsMatrix::removeEmptyProcessesInPlace");
372 if (Rows() == 1 && Cols () == 1) {
373 getMatrix(0,0)->removeEmptyProcessesInPlace(newMap);
374 return;
375 }
376 throw Xpetra::Exceptions::RuntimeError("removeEmptyProcesses not supported by BlockedCrsMatrix");
377 }
378
380
388 void replaceGlobalValues(GlobalOrdinal globalRow,
390 const ArrayView<const Scalar> &vals) {
391 XPETRA_MONITOR("XpetraBlockedCrsMatrix::replaceGlobalValues");
392 if (Rows() == 1 && Cols () == 1) {
393 getMatrix(0,0)->replaceGlobalValues(globalRow,cols,vals);
394 return;
395 }
396 throw Xpetra::Exceptions::RuntimeError("replaceGlobalValues not supported by BlockedCrsMatrix");
397 }
398
400
404 void replaceLocalValues(LocalOrdinal localRow,
406 const ArrayView<const Scalar> &vals) {
407 XPETRA_MONITOR("XpetraBlockedCrsMatrix::replaceLocalValues");
408 if (Rows() == 1 && Cols () == 1) {
409 getMatrix(0,0)->replaceLocalValues(localRow,cols,vals);
410 return;
411 }
412 throw Xpetra::Exceptions::RuntimeError("replaceLocalValues not supported by BlockedCrsMatrix");
413 }
414
416 virtual void setAllToScalar(const Scalar& alpha) {
417 XPETRA_MONITOR("XpetraBlockedCrsMatrix::setAllToScalar");
418 for (size_t row = 0; row < Rows(); row++) {
419 for (size_t col = 0; col < Cols(); col++) {
420 if (!getMatrix(row,col).is_null()) {
421 getMatrix(row,col)->setAllToScalar(alpha);
422 }
423 }
424 }
425 }
426
428 void scale(const Scalar& alpha) {
429 XPETRA_MONITOR("XpetraBlockedCrsMatrix::scale");
430 for (size_t row = 0; row < Rows(); row++) {
431 for (size_t col = 0; col < Cols(); col++) {
432 if (!getMatrix(row,col).is_null()) {
433 getMatrix(row,col)->scale(alpha);
434 }
435 }
436 }
437 }
438
440
442
443
451 void resumeFill(const RCP< ParameterList >& params = null) {
452 XPETRA_MONITOR("XpetraBlockedCrsMatrix::resumeFill");
453 for (size_t row = 0; row < Rows(); row++) {
454 for (size_t col = 0; col < Cols(); col++) {
455 if (!getMatrix(row,col).is_null()) {
456 getMatrix(row,col)->resumeFill(params);
457 }
458 }
459 }
460 }
461
467 void fillComplete(const RCP<const Map>& domainMap, const RCP<const Map>& rangeMap, const RCP<ParameterList>& params = null) {
468 XPETRA_MONITOR("XpetraBlockedCrsMatrix::fillComplete");
469 if (Rows() == 1 && Cols () == 1) {
470 getMatrix(0,0)->fillComplete(domainMap, rangeMap, params);
471 return;
472 }
473 fillComplete(params);
474 }
475
490 void fillComplete(const RCP<ParameterList>& params = null) {
491 XPETRA_MONITOR("XpetraBlockedCrsMatrix::fillComplete");
492 TEUCHOS_TEST_FOR_EXCEPTION(rangemaps_==Teuchos::null, Xpetra::Exceptions::RuntimeError,"BlockedCrsMatrix::fillComplete: rangemaps_ is not set. Error.");
493
494 for (size_t r = 0; r < Rows(); ++r)
495 for (size_t c = 0; c < Cols(); ++c) {
496 if(getMatrix(r,c) != Teuchos::null) {
497 Teuchos::RCP<Matrix> Ablock = getMatrix(r,c);
498 if(r!=c) is_diagonal_ = false;
499 if (Ablock != Teuchos::null && !Ablock->isFillComplete())
500 Ablock->fillComplete(getDomainMap(c, bDomainThyraMode_), getRangeMap(r, bRangeThyraMode_), params);
501 }
502 }
503
504#if 0
505 // get full row map
506 RCP<const Map> rangeMap = rangemaps_->getFullMap();
507 fullrowmap_ = MapFactory::Build(rangeMap()->lib(), rangeMap()->getGlobalNumElements(), rangeMap()->getNodeElementList(), rangeMap()->getIndexBase(), rangeMap()->getComm());
508
509 // build full col map
510 fullcolmap_ = Teuchos::null; // delete old full column map
511
512 std::vector<GO> colmapentries;
513 for (size_t c = 0; c < Cols(); ++c) {
514 // copy all local column lids of all block rows to colset
515 std::set<GO> colset;
516 for (size_t r = 0; r < Rows(); ++r) {
518
519 if (Ablock != Teuchos::null) {
520 Teuchos::ArrayView<const GO> colElements = Ablock->getColMap()->getNodeElementList();
521 Teuchos::RCP<const Map> colmap = Ablock->getColMap();
522 copy(colElements.getRawPtr(), colElements.getRawPtr() + colElements.size(), inserter(colset, colset.begin()));
523 }
524 }
525
526 // remove duplicates (entries which are in column maps of more than one block row)
527 colmapentries.reserve(colmapentries.size() + colset.size());
528 copy(colset.begin(), colset.end(), back_inserter(colmapentries));
529 sort(colmapentries.begin(), colmapentries.end());
530 typename std::vector<GO>::iterator gendLocation;
531 gendLocation = std::unique(colmapentries.begin(), colmapentries.end());
532 colmapentries.erase(gendLocation,colmapentries.end());
533 }
534
535 // sum up number of local elements
536 size_t numGlobalElements = 0;
537 Teuchos::reduceAll(*(rangeMap->getComm()), Teuchos::REDUCE_SUM, colmapentries.size(), Teuchos::outArg(numGlobalElements));
538
539 // store global full column map
541 fullcolmap_ = MapFactory::Build(rangeMap->lib(), numGlobalElements, aView, 0, rangeMap->getComm());
542#endif
543 }
544
546
548
551 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalNumRows");
552 global_size_t globalNumRows = 0;
553
554 for (size_t row = 0; row < Rows(); row++)
555 for (size_t col = 0; col < Cols(); col++)
556 if (!getMatrix(row,col).is_null()) {
557 globalNumRows += getMatrix(row,col)->getGlobalNumRows();
558 break; // we need only one non-null matrix in a row
559 }
560
561 return globalNumRows;
562 }
563
565
568 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalNumCols");
569 global_size_t globalNumCols = 0;
570
571 for (size_t col = 0; col < Cols(); col++)
572 for (size_t row = 0; row < Rows(); row++)
573 if (!getMatrix(row,col).is_null()) {
574 globalNumCols += getMatrix(row,col)->getGlobalNumCols();
575 break; // we need only one non-null matrix in a col
576 }
577
578 return globalNumCols;
579 }
580
582 size_t getNodeNumRows() const {
583 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNodeNumRows");
584 global_size_t nodeNumRows = 0;
585
586 for (size_t row = 0; row < Rows(); ++row)
587 for (size_t col = 0; col < Cols(); col++)
588 if (!getMatrix(row,col).is_null()) {
589 nodeNumRows += getMatrix(row,col)->getNodeNumRows();
590 break; // we need only one non-null matrix in a row
591 }
592
593 return nodeNumRows;
594 }
595
598 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalNumEntries");
599 global_size_t globalNumEntries = 0;
600
601 for (size_t row = 0; row < Rows(); ++row)
602 for (size_t col = 0; col < Cols(); ++col)
603 if (!getMatrix(row,col).is_null())
604 globalNumEntries += getMatrix(row,col)->getGlobalNumEntries();
605
606 return globalNumEntries;
607 }
608
610 size_t getNodeNumEntries() const {
611 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNodeNumEntries");
612 global_size_t nodeNumEntries = 0;
613
614 for (size_t row = 0; row < Rows(); ++row)
615 for (size_t col = 0; col < Cols(); ++col)
616 if (!getMatrix(row,col).is_null())
617 nodeNumEntries += getMatrix(row,col)->getNodeNumEntries();
618
619 return nodeNumEntries;
620 }
621
623
624 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const {
625 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNumEntriesInLocalRow");
626 GlobalOrdinal gid = this->getRowMap()->getGlobalElement(localRow);
627 size_t row = getBlockedRangeMap()->getMapIndexForGID(gid);
628 LocalOrdinal lid = getBlockedRangeMap()->getMap(row)->getLocalElement(gid);
629 size_t numEntriesInLocalRow = 0;
630 for (size_t col = 0; col < Cols(); ++col)
631 if (!getMatrix(row,col).is_null())
632 numEntriesInLocalRow += getMatrix(row,col)->getNumEntriesInLocalRow(lid);
633 return numEntriesInLocalRow;
634 }
635
637
638 size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const {
639 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNumEntriesInGlobalRow");
640 size_t row = getBlockedRangeMap()->getMapIndexForGID(globalRow);
641 size_t numEntriesInGlobalRow = 0;
642 for (size_t col = 0; col < Cols(); ++col)
643 if (!getMatrix(row,col).is_null())
644 numEntriesInGlobalRow += getMatrix(row,col)->getNumEntriesInGlobalRow(globalRow);
645 return numEntriesInGlobalRow;
646 }
647
649
652 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalMaxNumRowEntries");
653 global_size_t globalMaxEntries = 0;
654
655 for (size_t row = 0; row < Rows(); row++) {
656 global_size_t globalMaxEntriesBlockRows = 0;
657 for (size_t col = 0; col < Cols(); col++) {
658 if (!getMatrix(row,col).is_null()) {
659 globalMaxEntriesBlockRows += getMatrix(row,col)->getGlobalMaxNumRowEntries();
660 }
661 }
662 if(globalMaxEntriesBlockRows > globalMaxEntries)
663 globalMaxEntries = globalMaxEntriesBlockRows;
664 }
665 return globalMaxEntries;
666 }
667
669
671 size_t getNodeMaxNumRowEntries() const {
672 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNodeMaxNumRowEntries");
673 size_t localMaxEntries = 0;
674
675 for (size_t row = 0; row < Rows(); row++) {
676 size_t localMaxEntriesBlockRows = 0;
677 for (size_t col = 0; col < Cols(); col++) {
678 if (!getMatrix(row,col).is_null()) {
679 localMaxEntriesBlockRows += getMatrix(row,col)->getNodeMaxNumRowEntries();
680 }
681 }
682 if(localMaxEntriesBlockRows > localMaxEntries)
683 localMaxEntries = localMaxEntriesBlockRows;
684 }
685 return localMaxEntries;
686 }
687
689
692 bool isLocallyIndexed() const {
693 XPETRA_MONITOR("XpetraBlockedCrsMatrix::isLocallyIndexed");
694 for (size_t i = 0; i < blocks_.size(); ++i)
695 if (blocks_[i] != Teuchos::null && !blocks_[i]->isLocallyIndexed())
696 return false;
697 return true;
698 }
699
701
704 bool isGloballyIndexed() const {
705 XPETRA_MONITOR("XpetraBlockedCrsMatrix::isGloballyIndexed");
706 for (size_t i = 0; i < blocks_.size(); i++)
707 if (blocks_[i] != Teuchos::null && !blocks_[i]->isGloballyIndexed())
708 return false;
709 return true;
710 }
711
713 bool isFillComplete() const {
714 XPETRA_MONITOR("XpetraBlockedCrsMatrix::isFillComplete");
715 for (size_t i = 0; i < blocks_.size(); i++)
716 if (blocks_[i] != Teuchos::null && !blocks_[i]->isFillComplete())
717 return false;
718 return true;
719 }
720
722
736 virtual void getLocalRowCopy(LocalOrdinal LocalRow,
737 const ArrayView<LocalOrdinal>& Indices,
738 const ArrayView<Scalar>& Values,
739 size_t &NumEntries) const {
740 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalRowCopy");
741 if (Rows() == 1 && Cols () == 1) {
742 getMatrix(0,0)->getLocalRowCopy(LocalRow, Indices, Values, NumEntries);
743 return;
744 }
745 throw Xpetra::Exceptions::RuntimeError("getLocalRowCopy not supported by BlockedCrsMatrix" );
746 }
747
749
758 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal>& indices, ArrayView<const Scalar>& values) const {
759 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalRowView");
760 if (Rows() == 1 && Cols () == 1) {
761 getMatrix(0,0)->getGlobalRowView(GlobalRow, indices, values);
762 return;
763 }
764 throw Xpetra::Exceptions::RuntimeError("getGlobalRowView not supported by BlockedCrsMatrix");
765 }
766
768
777 void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal>& indices, ArrayView<const Scalar>& values) const {
778 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalRowView");
779 if (Rows() == 1 && Cols () == 1) {
780 getMatrix(0,0)->getLocalRowView(LocalRow, indices, values);
781 return;
782 }
783 else if(is_diagonal_) {
784 GlobalOrdinal gid = this->getRowMap()->getGlobalElement(LocalRow);
785 size_t row = getBlockedRangeMap()->getMapIndexForGID(gid);
786 getMatrix(row,row)->getLocalRowView(getMatrix(row,row)->getRowMap()->getLocalElement(gid),indices,values);
787 return;
788 }
789 throw Xpetra::Exceptions::RuntimeError("getLocalRowView not supported by BlockedCrsMatrix");
790 }
791
793
796 void getLocalDiagCopy(Vector& diag) const {
797 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalDiagCopy");
798
799 //RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
800
801 Teuchos::RCP<Vector> rcpdiag = Teuchos::rcpFromRef(diag);
802 Teuchos::RCP<BlockedVector> bdiag = Teuchos::rcp_dynamic_cast<BlockedVector>(rcpdiag);
803
804 // special treatment for 1x1 block matrices
805 // ReorderedBlockedCrsMatrix object encapsulate single blocks in ReorderedBlocks
806 // BlockedVectors have Vector objects as Leaf objects.
807 if(Rows() == 1 && Cols() == 1 && bdiag.is_null() == true) {
809 rm->getLocalDiagCopy(diag);
810 return;
811 }
812
813 TEUCHOS_TEST_FOR_EXCEPTION(bdiag.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::getLocalDiagCopy: diag must be a Blocked(Multi)Vector.");
814 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());
815 TEUCHOS_TEST_FOR_EXCEPTION(bdiag->getBlockedMap()->getNumMaps() != this->Rows(), Xpetra::Exceptions::RuntimeError,
816 "BlockedCrsMatrix::getLocalDiagCopy(): the number of blocks in diag differ from the number of blocks in this operator." );
817 //XPETRA_TEST_FOR_EXCEPTION(bdiag->getMap()->isSameAs(*(getMap())) == false, Xpetra::Exceptions::RuntimeError,
818 // "BlockedCrsMatrix::getLocalDiagCopy(): the map of the vector diag is not compatible with the map of the blocked operator." );
819
820 for (size_t row = 0; row < Rows(); row++) {
822 if (!rm.is_null()) {
823 Teuchos::RCP<Vector> rv = VectorFactory::Build(bdiag->getBlockedMap()->getMap(row,bdiag->getBlockedMap()->getThyraMode()));
824 rm->getLocalDiagCopy(*rv);
825 bdiag->setMultiVector(row,rv,bdiag->getBlockedMap()->getThyraMode());
826 }
827 }
828 }
829
831 void leftScale (const Vector& x) {
832 XPETRA_MONITOR("XpetraBlockedCrsMatrix::leftScale");
833
834 Teuchos::RCP<const Vector> rcpx = Teuchos::rcpFromRef(x);
835 Teuchos::RCP<const BlockedVector> bx = Teuchos::rcp_dynamic_cast<const BlockedVector>(rcpx);
836
837 //RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
838
839 // special treatment for 1xn block matrices
840 // ReorderedBlockedCrsMatrix object encapsulate single blocks in ReorderedBlocks
841 // BlockedVectors have Vector objects as Leaf objects.
842 if(Rows() == 1 && bx.is_null() == true) {
843 for (size_t col = 0; col < Cols(); ++col) {
845 rm->leftScale(x);
846 }
847 return;
848 }
849
850 TEUCHOS_TEST_FOR_EXCEPTION(bx.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::leftScale: x must be a Blocked(Multi)Vector.");
851 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());
852 TEUCHOS_TEST_FOR_EXCEPTION(bx->getBlockedMap()->getNumMaps() != this->Rows(), Xpetra::Exceptions::RuntimeError,
853 "BlockedCrsMatrix::leftScale(): the number of blocks in diag differ from the number of blocks in this operator." );
854
855 for (size_t row = 0; row < Rows(); row++) {
856 Teuchos::RCP<const MultiVector> rmv = bx->getMultiVector(row);
857 Teuchos::RCP<const Vector> rscale = rmv->getVector(0);
858 XPETRA_TEST_FOR_EXCEPTION(rscale.is_null()==true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::leftScale: x must be a Vector.");
859 for (size_t col = 0; col < Cols(); ++col) {
860 Teuchos::RCP<Matrix> rm = getMatrix(row,col);
861 if (!rm.is_null()) {
862 rm->leftScale(*rscale);
863 }
864 }
865 }
866 }
867
869 void rightScale (const Vector& x) {
870 XPETRA_MONITOR("XpetraBlockedCrsMatrix::rightScale");
871
872 Teuchos::RCP<const Vector> rcpx = Teuchos::rcpFromRef(x);
873 Teuchos::RCP<const BlockedVector> bx = Teuchos::rcp_dynamic_cast<const BlockedVector>(rcpx);
874
875 //RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
876
877 // special treatment for nx1 block matrices
878 // ReorderedBlockedCrsMatrix object encapsulate single blocks in ReorderedBlocks
879 // BlockedVectors have Vector objects as Leaf objects.
880 if(Cols() == 1 && bx.is_null() == true) {
881 for (size_t row = 0; row < Rows(); ++row) {
883 rm->rightScale(x);
884 }
885 return;
886 }
887
888 TEUCHOS_TEST_FOR_EXCEPTION(bx.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::rightScale: x must be a Blocked(Multi)Vector.");
889 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());
890 TEUCHOS_TEST_FOR_EXCEPTION(bx->getBlockedMap()->getNumMaps() != this->Cols(), Xpetra::Exceptions::RuntimeError,
891 "BlockedCrsMatrix::rightScale(): the number of blocks in diag differ from the number of blocks in this operator." );
892
893 for (size_t col = 0; col < Cols(); ++col) {
894 Teuchos::RCP<const MultiVector> rmv = bx->getMultiVector(col);
895 Teuchos::RCP<const Vector> rscale = rmv->getVector(0);
896 XPETRA_TEST_FOR_EXCEPTION(rscale.is_null()==true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::leftScale: x must be a Vector.");
897 for (size_t row = 0; row < Rows(); row++) {
898 Teuchos::RCP<Matrix> rm = getMatrix(row,col);
899 if (!rm.is_null()) {
900 rm->rightScale(*rscale);
901 }
902 }
903 }
904 }
905
906
909 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getFrobeniusNorm");
911 for (size_t col = 0; col < Cols(); ++col) {
912 for (size_t row = 0; row < Rows(); ++row) {
913 if(getMatrix(row,col)!=Teuchos::null) {
914 typename ScalarTraits<Scalar>::magnitudeType n = getMatrix(row,col)->getFrobeniusNorm();
915 ret += n * n;
916 }
917 }
918 }
920 }
921
922
924 virtual bool haveGlobalConstants() const {return true;}
925
926
928
930
931
933
953
955
956
958 virtual void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues,
959 const RCP<Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node> >& regionInterfaceImporter,
960 const Teuchos::ArrayRCP<LocalOrdinal>& regionInterfaceLIDs) const
961 { }
962
964
967 virtual void apply(const MultiVector& X, MultiVector& Y,
969 Scalar alpha = ScalarTraits<Scalar>::one(),
970 Scalar beta = ScalarTraits<Scalar>::zero()) const
971 {
972 XPETRA_MONITOR("XpetraBlockedCrsMatrix::apply");
973 //using Teuchos::RCP;
974
976 "apply() only supports the following modes: NO_TRANS and TRANS." );
977
978 // check whether input parameters are blocked or not
979 RCP<const MultiVector> refX = rcpFromRef(X);
980 RCP<const BlockedMultiVector> refbX = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(refX);
981 //RCP<MultiVector> tmpY = rcpFromRef(Y);
982 //RCP<BlockedMultiVector> tmpbY = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(tmpY);
983
984 // TODO get rid of me: adapt MapExtractor
985 bool bBlockedX = (refbX != Teuchos::null) ? true : false;
986
987 // create (temporary) vectors for output
988 // In the end we call Y.update(alpha, *tmpY, beta). Therefore we need a new vector storing the temporary results
990
991 //RCP<Teuchos::FancyOStream> out = rcp(new Teuchos::FancyOStream(rcp(&std::cout,false)));
992
994
995 if (mode == Teuchos::NO_TRANS) {
996
997 for (size_t row = 0; row < Rows(); row++) {
998 RCP<MultiVector> Yblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, true);
999 for (size_t col = 0; col < Cols(); col++) {
1000
1001 // extract matrix block
1002 RCP<Matrix> Ablock = getMatrix(row, col);
1003
1004 if (Ablock.is_null())
1005 continue;
1006
1007 // check whether Ablock is itself a blocked operator
1008 // If it is a blocked operator we have to provide Xpetra style GIDs, i.e. we have to transform GIDs
1009 bool bBlockedSubMatrix = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Ablock) == Teuchos::null ? false : true;
1010
1011 // input/output vectors for local block operation
1012 RCP<const MultiVector> Xblock = Teuchos::null; // subpart of X vector to be applied to subblock of A
1013
1014 // extract sub part of X using Xpetra or Thyra GIDs
1015 // if submatrix is again blocked, we extract it using Xpetra style gids. If it is a single
1016 // block matrix we use the Thyra or Xpetra style GIDs that are used to store the matrix
1017 if(bBlockedX) Xblock = domainmaps_->ExtractVector(refbX, col, bDomainThyraMode_);
1018 else Xblock = domainmaps_->ExtractVector(refX, col, bBlockedSubMatrix == true ? false : bDomainThyraMode_);
1019
1020 RCP<MultiVector> tmpYblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, false); // subpart of Y vector containing part of solution of Xblock applied to Ablock
1021 Ablock->apply(*Xblock, *tmpYblock);
1022
1023 Yblock->update(one, *tmpYblock, one);
1024 }
1025 rangemaps_->InsertVector(Yblock, row, tmpY, bRangeThyraMode_);
1026 }
1027
1028 } else if (mode == Teuchos::TRANS) {
1029 // TODO: test me!
1030 for (size_t col = 0; col < Cols(); col++) {
1031 RCP<MultiVector> Yblock = domainmaps_->getVector(col, Y.getNumVectors(), bDomainThyraMode_, true);
1032
1033 for (size_t row = 0; row<Rows(); row++) {
1034 RCP<Matrix> Ablock = getMatrix(row, col);
1035
1036 if (Ablock.is_null())
1037 continue;
1038
1039 // check whether Ablock is itself a blocked operator
1040 // If it is a blocked operator we have to provide Xpetra style GIDs, i.e. we have to transform GIDs
1041 bool bBlockedSubMatrix = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Ablock) == Teuchos::null ? false : true;
1042
1043 RCP<const MultiVector> Xblock = Teuchos::null;
1044
1045 // extract sub part of X using Xpetra or Thyra GIDs
1046 if(bBlockedX) Xblock = rangemaps_->ExtractVector(refbX, row, bRangeThyraMode_);
1047 else Xblock = rangemaps_->ExtractVector(refX, row, bBlockedSubMatrix == true ? false : bRangeThyraMode_);
1048 RCP<MultiVector> tmpYblock = domainmaps_->getVector(col, Y.getNumVectors(), bDomainThyraMode_, false);
1049 Ablock->apply(*Xblock, *tmpYblock, Teuchos::TRANS);
1050
1051 Yblock->update(one, *tmpYblock, one);
1052 }
1053 domainmaps_->InsertVector(Yblock, col, tmpY, bDomainThyraMode_);
1054 }
1055 }
1056 Y.update(alpha, *tmpY, beta);
1057 }
1058
1060 RCP<const Map > getFullDomainMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getFullDomainMap()"); return domainmaps_->getFullMap(); }
1061
1063 RCP<const BlockedMap > getBlockedDomainMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getBlockedDomainMap()"); return domainmaps_->getBlockedMap(); }
1064
1066 RCP<const Map > getDomainMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMap()"); return domainmaps_->getMap(); /*domainmaps_->getFullMap();*/ }
1067
1069 RCP<const Map > getDomainMap(size_t i) const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMap(size_t)"); return domainmaps_->getMap(i, bDomainThyraMode_); }
1070
1072 RCP<const Map > getDomainMap(size_t i, bool bThyraMode) const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMap(size_t,bool)"); return domainmaps_->getMap(i, bThyraMode); }
1073
1075 RCP<const Map > getFullRangeMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap()"); return rangemaps_->getFullMap(); }
1076
1078 RCP<const BlockedMap > getBlockedRangeMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getBlockedRangeMap()"); return rangemaps_->getBlockedMap(); }
1079
1081 RCP<const Map > getRangeMap() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap()"); return rangemaps_->getMap(); /*rangemaps_->getFullMap();*/ }
1082
1084 RCP<const Map > getRangeMap(size_t i) const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap(size_t)"); return rangemaps_->getMap(i, bRangeThyraMode_); }
1085
1087 RCP<const Map > getRangeMap(size_t i, bool bThyraMode) const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap(size_t,bool)"); return rangemaps_->getMap(i, bThyraMode); }
1088
1090 RCP<const MapExtractor> getRangeMapExtractor() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMapExtractor()"); return rangemaps_; }
1091
1093 RCP<const MapExtractor> getDomainMapExtractor() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMapExtractor()"); return domainmaps_; }
1094
1096
1098 //{@
1099
1108 virtual void bgs_apply(
1109 const MultiVector& X,
1110 MultiVector& Y,
1111 size_t row,
1113 Scalar alpha = ScalarTraits<Scalar>::one(),
1114 Scalar beta = ScalarTraits<Scalar>::zero()
1115 ) const
1116 {
1117 XPETRA_MONITOR("XpetraBlockedCrsMatrix::bgs_apply");
1118 //using Teuchos::RCP;
1119
1121 "apply() only supports the following modes: NO_TRANS and TRANS." );
1122
1123 // check whether input parameters are blocked or not
1124 RCP<const MultiVector> refX = rcpFromRef(X);
1125 RCP<const BlockedMultiVector> refbX = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(refX);
1126 //RCP<MultiVector> tmpY = rcpFromRef(Y);
1127 //RCP<BlockedMultiVector> tmpbY = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(tmpY);
1128
1129 bool bBlockedX = (refbX != Teuchos::null) ? true : false;
1130
1131 // create (temporary) vectors for output
1132 // In the end we call Y.update(alpha, *tmpY, beta). Therefore we need a new vector storing the temporary results
1134
1135 SC one = ScalarTraits<SC>::one();
1136
1137 if (mode == Teuchos::NO_TRANS) {
1138 RCP<MultiVector> Yblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, true);
1139 for (size_t col = 0; col < Cols(); col++) {
1140
1141 // extract matrix block
1142 RCP<Matrix> Ablock = getMatrix(row, col);
1143
1144 if (Ablock.is_null())
1145 continue;
1146
1147 // check whether Ablock is itself a blocked operator
1148 // If it is a blocked operator we have to provide Xpetra style GIDs, i.e. we have to transform GIDs
1149 bool bBlockedSubMatrix = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Ablock) == Teuchos::null ? false : true;
1150
1151 // input/output vectors for local block operation
1152 RCP<const MultiVector> Xblock = Teuchos::null; // subpart of X vector to be applied to subblock of A
1153
1154 // extract sub part of X using Xpetra or Thyra GIDs
1155 // if submatrix is again blocked, we extract it using Xpetra style gids. If it is a single
1156 // block matrix we use the Thyra or Xpetra style GIDs that are used to store the matrix
1157 if(bBlockedX) Xblock = domainmaps_->ExtractVector(refbX, col, bDomainThyraMode_);
1158 else Xblock = domainmaps_->ExtractVector(refX, col, bBlockedSubMatrix == true ? false : bDomainThyraMode_);
1159
1160 RCP<MultiVector> tmpYblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, false); // subpart of Y vector containing part of solution of Xblock applied to Ablock
1161 Ablock->apply(*Xblock, *tmpYblock);
1162
1163 Yblock->update(one, *tmpYblock, one);
1164 }
1165 rangemaps_->InsertVector(Yblock, row, tmpY, bRangeThyraMode_);
1166 } else {
1167 TEUCHOS_TEST_FOR_EXCEPTION(true,Xpetra::Exceptions::NotImplemented,"Xpetar::BlockedCrsMatrix::bgs_apply: not implemented for transpose case.");
1168 }
1169 Y.update(alpha, *tmpY, beta);
1170 }
1171
1172
1174
1175
1177 //{@
1178
1181 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getMap");
1182 if (Rows() == 1 && Cols () == 1) {
1183 return getMatrix(0,0)->getMap();
1184 }
1185 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::getMap(): operation not supported.");
1186 }
1187
1189 void doImport(const Matrix &source, const Import& importer, CombineMode CM) {
1190 XPETRA_MONITOR("XpetraBlockedCrsMatrix::doImport");
1191 if (Rows() == 1 && Cols () == 1) {
1192 getMatrix(0,0)->doImport(source, importer, CM);
1193 return;
1194 }
1195 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doImport(): operation not supported.");
1196 }
1197
1199 void doExport(const Matrix& dest, const Import& importer, CombineMode CM) {
1200 XPETRA_MONITOR("XpetraBlockedCrsMatrix::doExport");
1201 if (Rows() == 1 && Cols () == 1) {
1202 getMatrix(0,0)->doExport(dest, importer, CM);
1203 return;
1204 }
1205 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doExport(): operation not supported.");
1206 }
1207
1209 void doImport(const Matrix& source, const Export& exporter, CombineMode CM) {
1210 XPETRA_MONITOR("XpetraBlockedCrsMatrix::doImport");
1211 if (Rows() == 1 && Cols () == 1) {
1212 getMatrix(0,0)->doImport(source, exporter, CM);
1213 return;
1214 }
1215 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doImport(): operation not supported.");
1216 }
1217
1219 void doExport(const Matrix& dest, const Export& exporter, CombineMode CM) {
1220 XPETRA_MONITOR("XpetraBlockedCrsMatrix::doExport");
1221 if (Rows() == 1 && Cols () == 1) {
1222 getMatrix(0,0)->doExport(dest, exporter, CM);
1223 return;
1224 }
1225 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doExport(): operation not supported.");
1226 }
1227
1228 // @}
1229
1231
1232
1234 std::string description() const { return "Xpetra_BlockedCrsMatrix.description()"; }
1235
1238 out << "Xpetra::BlockedCrsMatrix: " << Rows() << " x " << Cols() << std::endl;
1239
1240 if (isFillComplete()) {
1241 out << "BlockMatrix is fillComplete" << std::endl;
1242
1243 /*if(fullrowmap_ != Teuchos::null) {
1244 out << "fullRowMap" << std::endl;
1245 fullrowmap_->describe(out,verbLevel);
1246 } else {
1247 out << "fullRowMap not set. Check whether block matrix is properly fillCompleted!" << std::endl;
1248 }*/
1249
1250 //out << "fullColMap" << std::endl;
1251 //fullcolmap_->describe(out,verbLevel);
1252
1253 } else {
1254 out << "BlockMatrix is NOT fillComplete" << std::endl;
1255 }
1256
1257 for (size_t r = 0; r < Rows(); ++r)
1258 for (size_t c = 0; c < Cols(); ++c) {
1259 if(getMatrix(r,c)!=Teuchos::null) {
1260 out << "Block(" << r << "," << c << ")" << std::endl;
1261 getMatrix(r,c)->describe(out,verbLevel);
1262 } else out << "Block(" << r << "," << c << ") = null" << std::endl;
1263 }
1264 }
1265
1267
1268 void setObjectLabel( const std::string &objectLabel ) {
1269 XPETRA_MONITOR("TpetraBlockedCrsMatrix::setObjectLabel");
1270 for (size_t r = 0; r < Rows(); ++r)
1271 for (size_t c = 0; c < Cols(); ++c) {
1272 if(getMatrix(r,c)!=Teuchos::null) {
1273 std::ostringstream oss; oss<< objectLabel << "(" << r << "," << c << ")";
1274 getMatrix(r,c)->setObjectLabel(oss.str());
1275 }
1276 }
1277 }
1279
1280
1282 bool hasCrsGraph() const {
1283 if (Rows() == 1 && Cols () == 1) return true;
1284 else return false;
1285 }
1286
1289 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getCrsGraph");
1290 if (Rows() == 1 && Cols () == 1) {
1291 return getMatrix(0,0)->getCrsGraph();
1292 }
1293 throw Xpetra::Exceptions::RuntimeError("getCrsGraph() not supported by BlockedCrsMatrix");
1294 }
1295
1297
1299
1300
1301 virtual bool isDiagonal() const {return is_diagonal_;}
1302
1304 virtual size_t Rows() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::Rows"); return rangemaps_->NumMaps(); }
1305
1307 virtual size_t Cols() const { XPETRA_MONITOR("XpetraBlockedCrsMatrix::Cols"); return domainmaps_->NumMaps(); }
1308
1311 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getCrsMatrix");
1312 TEUCHOS_TEST_FOR_EXCEPTION(Rows()!=1, std::out_of_range, "Can only unwrap a 1x1 blocked matrix. The matrix has " << Rows() << " block rows, though.");
1313 TEUCHOS_TEST_FOR_EXCEPTION(Cols()!=1, std::out_of_range, "Can only unwrap a 1x1 blocked matrix. The matrix has " << Cols() << " block columns, though.");
1314
1315 RCP<Matrix> mat = getMatrix(0,0);
1316 RCP<BlockedCrsMatrix> bmat = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(mat);
1317 if (bmat == Teuchos::null) return mat;
1318 return bmat->getCrsMatrix();
1319 }
1320
1323 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getInnermostCrsMatrix");
1324 size_t row = Rows()+1, col = Cols()+1;
1325 for (size_t r = 0; r < Rows(); ++r)
1326 for(size_t c = 0; c < Cols(); ++c)
1327 if (getMatrix(r,c) != Teuchos::null) {
1328 row = r;
1329 col = c;
1330 break;
1331 }
1332 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.")
1333 RCP<Matrix> mm = getMatrix(row,col);
1334 RCP<BlockedCrsMatrix> bmat = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(mm);
1335 if (bmat == Teuchos::null) return mm;
1336 return bmat->getInnermostCrsMatrix();
1337 }
1338
1340 Teuchos::RCP<Matrix> getMatrix(size_t r, size_t c) const {
1341 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getMatrix");
1342 TEUCHOS_TEST_FOR_EXCEPTION(r > Rows(), std::out_of_range, "Error, r = " << Rows() << " is too big");
1343 TEUCHOS_TEST_FOR_EXCEPTION(c > Cols(), std::out_of_range, "Error, c = " << Cols() << " is too big");
1344
1345 // transfer strided/blocked map information
1346 /* if (blocks_[r*Cols()+c] != Teuchos::null &&
1347 blocks_[r*Cols()+c]->IsView("stridedMaps") == false)
1348 blocks_[r*Cols()+c]->CreateView("stridedMaps", getRangeMap(r,bRangeThyraMode_), getDomainMap(c,bDomainThyraMode_));*/
1349 return blocks_[r*Cols()+c];
1350 }
1351
1353 //void setMatrix(size_t r, size_t c, Teuchos::RCP<CrsMatrix> mat) {
1354 void setMatrix(size_t r, size_t c, Teuchos::RCP<Matrix> mat) {
1355 XPETRA_MONITOR("XpetraBlockedCrsMatrix::setMatrix");
1356 // TODO: if filled -> return error
1357
1358 TEUCHOS_TEST_FOR_EXCEPTION(r > Rows(), std::out_of_range, "Error, r = " << Rows() << " is too big");
1359 TEUCHOS_TEST_FOR_EXCEPTION(c > Cols(), std::out_of_range, "Error, c = " << Cols() << " is too big");
1360 if(!mat.is_null() && r!=c) is_diagonal_ = false;
1361 // set matrix
1362 blocks_[r*Cols() + c] = mat;
1363 }
1364
1366 // NOTE: This is a rather expensive operation, since all blocks are copied
1367 // into a new big CrsMatrix
1369 XPETRA_MONITOR("XpetraBlockedCrsMatrix::Merge");
1370 using Teuchos::RCP;
1371 using Teuchos::rcp_dynamic_cast;
1372 Scalar one = ScalarTraits<SC>::one();
1373
1375 "BlockedCrsMatrix::Merge: only implemented for Xpetra-style or Thyra-style numbering. No mixup allowed!" );
1376
1378 "BlockedCrsMatrix::Merge: BlockMatrix must be fill-completed." );
1379
1380 LocalOrdinal lclNumRows = getFullRangeMap()->getNodeNumElements();
1381 Teuchos::ArrayRCP<size_t> numEntPerRow (lclNumRows);
1382 for (LocalOrdinal lclRow = 0; lclRow < lclNumRows; ++lclRow)
1383 numEntPerRow[lclRow] = getNumEntriesInLocalRow(lclRow);
1384
1385 RCP<Matrix> sparse = MatrixFactory::Build(getFullRangeMap(), numEntPerRow);
1386
1387 if(bRangeThyraMode_ == false) {
1388 // Xpetra mode
1389 for (size_t i = 0; i < Rows(); i++) {
1390 for (size_t j = 0; j < Cols(); j++) {
1391 if (getMatrix(i,j) != Teuchos::null) {
1392 RCP<const Matrix> mat = getMatrix(i,j);
1393
1394 // recursively call Merge routine
1395 RCP<const BlockedCrsMatrix> bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1396 if (bMat != Teuchos::null) mat = bMat->Merge();
1397
1398 bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1400 "BlockedCrsMatrix::Merge: Merging of blocked sub-operators failed?!" );
1401
1402 // jump over empty blocks
1403 if(mat->getNodeNumEntries() == 0) continue;
1404
1405 this->Add(*mat, one, *sparse, one);
1406 }
1407 }
1408 }
1409 } else {
1410 // Thyra mode
1411 for (size_t i = 0; i < Rows(); i++) {
1412 for (size_t j = 0; j < Cols(); j++) {
1413 if (getMatrix(i,j) != Teuchos::null) {
1414 RCP<const Matrix> mat = getMatrix(i,j);
1415 // recursively call Merge routine
1416 RCP<const BlockedCrsMatrix> bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1417 if (bMat != Teuchos::null) mat = bMat->Merge();
1418
1419 bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1421 "BlockedCrsMatrix::Merge: Merging of blocked sub-operators failed?!" );
1422
1423 // check whether we have a CrsMatrix block (no blocked operator)
1424 RCP<const CrsMatrixWrap> crsMat = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(mat);
1425 TEUCHOS_ASSERT(crsMat != Teuchos::null);
1426
1427 // these are the thyra style maps of the matrix
1428 RCP<const Map> trowMap = mat->getRowMap();
1429 RCP<const Map> tcolMap = mat->getColMap();
1430 RCP<const Map> tdomMap = mat->getDomainMap();
1431
1432 // get Xpetra maps
1433 RCP<const Map> xrowMap = getRangeMapExtractor()->getMap(i,false);
1434 RCP<const Map> xdomMap = getDomainMapExtractor()->getMap(j,false);
1435
1436 // generate column map with Xpetra GIDs
1437 // We have to do this separately for each block since the column
1438 // map of each block might be different in the same block column
1440 *tcolMap,
1441 *tdomMap,
1442 *xdomMap);
1443
1444 // jump over empty blocks
1445 if(mat->getNodeNumEntries() == 0) continue;
1446
1447 size_t maxNumEntries = mat->getNodeMaxNumRowEntries();
1448
1449 size_t numEntries;
1450 Array<GO> inds (maxNumEntries);
1451 Array<GO> inds2(maxNumEntries);
1452 Array<SC> vals (maxNumEntries);
1453
1454 // loop over all rows and add entries
1455 for (size_t k = 0; k < mat->getNodeNumRows(); k++) {
1456 GlobalOrdinal rowTGID = trowMap->getGlobalElement(k);
1457 crsMat->getCrsMatrix()->getGlobalRowCopy(rowTGID, inds(), vals(), numEntries);
1458
1459 // create new indices array
1460 for (size_t l = 0; l < numEntries; ++l) {
1461 LocalOrdinal lid = tcolMap->getLocalElement(inds[l]);
1462 inds2[l] = xcolMap->getGlobalElement(lid);
1463 }
1464
1465 GlobalOrdinal rowXGID = xrowMap->getGlobalElement(k);
1466 sparse->insertGlobalValues(
1467 rowXGID, inds2(0, numEntries),
1468 vals(0, numEntries));
1469 }
1470 }
1471 }
1472 }
1473 }
1474
1475 sparse->fillComplete(getFullDomainMap(), getFullRangeMap());
1476
1478 "BlockedCrsMatrix::Merge: Local number of entries of merged matrix does not coincide with local number of entries of blocked operator." );
1479
1481 "BlockedCrsMatrix::Merge: Global number of entries of merged matrix does not coincide with global number of entries of blocked operator." );
1482
1483 return sparse;
1484 }
1486
1487#ifdef HAVE_XPETRA_KOKKOS_REFACTOR
1488 typedef typename CrsMatrix::local_matrix_type local_matrix_type;
1489#ifdef TPETRA_ENABLE_DEPRECATED_CODE
1490 local_matrix_type getLocalMatrix () const {
1491 return getLocalMatrixDevice();
1492 }
1493#endif
1495 local_matrix_type getLocalMatrixDevice () const {
1496 if (Rows() == 1 && Cols () == 1) {
1497 return getMatrix(0,0)->getLocalMatrixDevice();
1498 }
1499 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::getLocalMatrix(): operation not supported.");
1500 }
1502 typename local_matrix_type::HostMirror getLocalMatrixHost () const {
1503 if (Rows() == 1 && Cols () == 1) {
1504 return getMatrix(0,0)->getLocalMatrixHost();
1505 }
1506 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::getLocalMatrix(): operation not supported.");
1507 }
1508
1509
1510#endif
1511
1512#ifdef HAVE_XPETRA_THYRA
1515 Xpetra::ThyraUtils<Scalar,LO,GO,Node>::toThyra(Teuchos::rcpFromRef(*this));
1517
1519 Teuchos::rcp_dynamic_cast<Thyra::BlockedLinearOpBase<Scalar> >(thOp);
1521 return thbOp;
1522 }
1523#endif
1524
1526 void residual(const MultiVector & X,
1527 const MultiVector & B,
1528 MultiVector & R) const {
1530 R.update(STS::one(),B,STS::zero());
1531 this->apply (X, R, Teuchos::NO_TRANS, -STS::one(), STS::one());
1532 }
1533
1534 private:
1535
1538
1540
1550 void Add(const Matrix& A, const Scalar scalarA, Matrix& B, const Scalar scalarB) const {
1551 XPETRA_MONITOR("XpetraBlockedCrsMatrix::Add");
1553 "Matrix A is not completed");
1554 using Teuchos::Array;
1555 using Teuchos::ArrayView;
1556
1557 B.scale(scalarB);
1558
1559 Scalar one = ScalarTraits<SC>::one();
1560 Scalar zero = ScalarTraits<SC>::zero();
1561
1562 if (scalarA == zero)
1563 return;
1564
1565 Teuchos::RCP<const Matrix> rcpA = Teuchos::rcpFromRef(A);
1566 Teuchos::RCP<const CrsMatrixWrap> rcpAwrap = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(rcpA);
1568 "BlockedCrsMatrix::Add: matrix A must be of type CrsMatrixWrap.");
1569 Teuchos::RCP<const CrsMatrix> crsA = rcpAwrap->getCrsMatrix();
1570
1571 size_t maxNumEntries = crsA->getNodeMaxNumRowEntries();
1572
1573 size_t numEntries;
1574 Array<GO> inds(maxNumEntries);
1575 Array<SC> vals(maxNumEntries);
1576
1577 RCP<const Map> rowMap = crsA->getRowMap();
1578 RCP<const Map> colMap = crsA->getColMap();
1579
1580 ArrayView<const GO> rowGIDs = crsA->getRowMap()->getNodeElementList();
1581 for (size_t i = 0; i < crsA->getNodeNumRows(); i++) {
1582 GO row = rowGIDs[i];
1583 crsA->getGlobalRowCopy(row, inds(), vals(), numEntries);
1584
1585 if (scalarA != one)
1586 for (size_t j = 0; j < numEntries; ++j)
1587 vals[j] *= scalarA;
1588
1589 B.insertGlobalValues(row, inds(0, numEntries), vals(0, numEntries)); // insert should be ok, since blocks in BlockedCrsOpeartor do not overlap!
1590 }
1591 }
1592
1594
1595 // Default view is created after fillComplete()
1596 // Because ColMap might not be available before fillComplete().
1598 XPETRA_MONITOR("XpetraBlockedCrsMatrix::CreateDefaultView");
1599
1600 // Create default view
1601 this->defaultViewLabel_ = "point";
1603
1604 // Set current view
1605 this->currentViewLabel_ = this->GetDefaultViewLabel();
1606 }
1607
1608 private:
1612
1613 std::vector<Teuchos::RCP<Matrix> > blocks_;
1614#ifdef HAVE_XPETRA_THYRA
1616#endif
1619
1620};
1621
1622} //namespace Xpetra
1623
1624#define XPETRA_BLOCKEDCRSMATRIX_SHORT
1625#endif /* XPETRA_BLOCKEDCRSMATRIX_HPP */
#define XPETRA_MONITOR(funcName)
#define XPETRA_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
GlobalOrdinal GO
iterator end() const
iterator begin() const
size_type size() const
T * getRawPtr() const
static const EVerbosityLevel verbLevel_default
bool is_null() const
Teuchos::RCP< const MapExtractor > domainmaps_
full domain map together with all partial domain maps
void doImport(const Matrix &source, const Export &exporter, CombineMode CM)
Import (using an Exporter).
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.
virtual ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
Teuchos::RCP< Matrix > getMatrix(size_t r, size_t c) const
return block (r,c)
Teuchos::RCP< Matrix > getInnermostCrsMatrix()
helper routine recursively returns the first inner-most non-null matrix block from a (nested) blocked...
bool isLocallyIndexed() const
If matrix indices of all matrix blocks are in the local range, this function returns true....
Teuchos::RCP< const MapExtractor > rangemaps_
full range map together with all partial domain maps
RCP< const BlockedMap > getBlockedRangeMap() const
Returns the BlockedMap associated with the range of this operator.
void setMatrix(size_t r, size_t c, Teuchos::RCP< Matrix > mat)
set matrix block
RCP< const BlockedMap > getBlockedDomainMap() const
Returns the BlockedMap associated with the domain of this operator.
virtual size_t Rows() const
number of row blocks
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
Teuchos::RCP< Matrix > Merge() const
merge BlockedCrsMatrix blocks in a CrsMatrix
BlockedCrsMatrix(const Teuchos::RCP< const BlockedMap > &rangeMaps, const Teuchos::RCP< const BlockedMap > &domainMaps, size_t numEntriesPerRow)
Constructor.
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete.
bool bDomainThyraMode_
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub...
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
RCP< const Map > getRangeMap(size_t i) const
Returns the Map associated with the i'th block range of this operator.
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/Jacobi smoother)
void doExport(const Matrix &dest, const Export &exporter, CombineMode CM)
Export (using an Importer).
RCP< const Map > getDomainMap() const
Returns the Map associated with the domain of this operator.
bool hasCrsGraph() const
Supports the getCrsGraph() call.
std::vector< Teuchos::RCP< Matrix > > blocks_
row major matrix block storage
virtual ~BlockedCrsMatrix()
Destructor.
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...
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
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.
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode.
Teuchos::RCP< Matrix > getCrsMatrix() const
return unwrap 1x1 blocked operators
RCP< const MapExtractor > getRangeMapExtractor() const
Returns map extractor class for range map.
RCP< const Map > getRangeMap() const
Returns the Map associated with the range of this operator.
virtual void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP< Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > &regionInterfaceImporter, const Teuchos::ArrayRCP< LocalOrdinal > &regionInterfaceLIDs) const
sparse matrix-multivector multiplication for the region layout matrices (currently no blocked impleme...
RCP< const Map > getFullRangeMap() const
Returns the Map associated with the full range of this operator.
void leftScale(const Vector &x)
Left scale matrix using the given vector entries.
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
RCP< const Map > getRangeMap(size_t i, bool bThyraMode) const
Returns the Map associated with the i'th block range of this operator.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
void doImport(const Matrix &source, const Import &importer, CombineMode CM)
Import.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
bool is_diagonal_
If we're diagonal, a bunch of the extraction stuff should work.
const Teuchos::RCP< const Map > getMap() const
Implements DistObject interface.
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.
bool bRangeThyraMode_
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub...
BlockedCrsMatrix(Teuchos::RCP< const MapExtractor > &rangeMapExtractor, Teuchos::RCP< const MapExtractor > &domainMapExtractor, size_t numEntriesPerRow)
Constructor.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
virtual bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
RCP< const Map > getDomainMap(size_t i) const
Returns the Map associated with the i'th block domain of this operator.
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.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the specified (locally owned) global row.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
void rightScale(const Vector &x)
Right scale matrix using the given vector entries.
void setObjectLabel(const std::string &objectLabel)
void fillComplete(const RCP< ParameterList > &params=null)
Signal that data entry is complete.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
void resumeFill(const RCP< ParameterList > &params=null)
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
void getLocalDiagCopy(Vector &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
global_size_t getGlobalNumRows() const
Returns the number of global rows.
RCP< const MapExtractor > getDomainMapExtractor() const
Returns map extractor for domain map.
virtual size_t Cols() const
number of column blocks
RCP< const Map > getDomainMap(size_t i, bool bThyraMode) const
Returns the Map associated with the i'th block domain of this operator.
RCP< const Map > getFullDomainMap() const
Returns the Map associated with the full domain of this operator.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this 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.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
void doExport(const Matrix &dest, const Import &importer, CombineMode CM)
Export.
std::string description() const
Return a simple one-line description of this object.
Concrete implementation of Xpetra::Matrix.
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const =0
The Map describing the parallel distribution of this object.
Exception indicating invalid cast attempted.
Exception throws to report incompatible objects (like maps).
Exception throws when you call an unimplemented method of Xpetra.
Exception throws to report errors in the internal logical of the program.
static Teuchos::RCP< Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const Teuchos::RCP< const Map > &fullmap, const std::vector< Teuchos::RCP< const Map > > &maps, bool bThyraMode=false)
Build MapExtrasctor from a given set of partial maps.
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)
Map constructor with Xpetra-defined contiguous uniform distribution.
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
static RCP< Matrix > Build(const RCP< const Map > &rowMap)
Xpetra-specific matrix class.
void CreateView(viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap)
virtual bool isFillComplete() const =0
Returns true if fillComplete() has been called and the matrix is in compute mode.
virtual const RCP< const Map > & getRowMap() const
Returns the Map that describes the row distribution in this matrix.
virtual void scale(const Scalar &alpha)=0
Scale the current values of a matrix, this = alpha*this.
virtual void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Insert matrix entries, using global IDs.
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.
virtual size_t getNumVectors() const =0
Number of columns in the multivector.
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.
static Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &map, bool zeroOut=true)
Constructor specifying the number of non-zeros for all rows.
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
bool is_null(const std::shared_ptr< T > &p)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Xpetra namespace
std::string viewLabel_t
size_t global_size_t
Global size_t object.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
CombineMode
Xpetra::Combine Mode enumerable type.
static magnitudeType magnitude(T a)