Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_BlockedCrsMatrix_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef XPETRA_BLOCKEDCRSMATRIX_DEF_HPP
11#define XPETRA_BLOCKEDCRSMATRIX_DEF_HPP
12
13#include <Tpetra_KokkosCompat_DefaultNode.hpp>
14
16#include <Teuchos_Hashtable.hpp>
17
18#include "Xpetra_ConfigDefs.hpp"
19#include "Xpetra_Exceptions.hpp"
20
22#include "Xpetra_MapFactory.hpp"
23#include "Xpetra_MultiVector.hpp"
24#include "Xpetra_BlockedMultiVector.hpp"
25#include "Xpetra_MultiVectorFactory.hpp"
26#include "Xpetra_BlockedVector.hpp"
27#include "Xpetra_CrsGraph.hpp"
28#include "Xpetra_CrsMatrix.hpp"
30
31#include "Xpetra_MapExtractor.hpp"
33
34#include "Xpetra_Matrix.hpp"
35#include "Xpetra_MatrixFactory.hpp"
36#include "Xpetra_CrsMatrixWrap.hpp"
37
38#ifdef HAVE_XPETRA_THYRA
39#include <Thyra_ProductVectorSpaceBase.hpp>
40#include <Thyra_VectorSpaceBase.hpp>
41#include <Thyra_LinearOpBase.hpp>
42#include <Thyra_BlockedLinearOpBase.hpp>
43#include <Thyra_PhysicallyBlockedLinearOpBase.hpp>
44#include "Xpetra_ThyraUtils.hpp"
45#endif
46
47#include "Xpetra_VectorFactory.hpp"
48
53namespace Xpetra {
54
55template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
57 const Teuchos::RCP<const BlockedMap>& domainMaps,
58 size_t numEntriesPerRow)
59 : is_diagonal_(true) {
62 bRangeThyraMode_ = rangeMaps->getThyraMode();
63 bDomainThyraMode_ = domainMaps->getThyraMode();
64
65 blocks_.reserve(Rows() * Cols());
66
67 // add CrsMatrix objects in row,column order
68 for (size_t r = 0; r < Rows(); ++r)
69 for (size_t c = 0; c < Cols(); ++c)
70 blocks_.push_back(MatrixFactory::Build(getRangeMap(r, bRangeThyraMode_), numEntriesPerRow));
71
72 // Default view
74}
75
76template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
79 size_t numEntriesPerRow)
80 : is_diagonal_(true)
81 , domainmaps_(domainMapExtractor)
82 , rangemaps_(rangeMapExtractor) {
83 bRangeThyraMode_ = rangeMapExtractor->getThyraMode();
84 bDomainThyraMode_ = domainMapExtractor->getThyraMode();
85
86 blocks_.reserve(Rows() * Cols());
87
88 // add CrsMatrix objects in row,column order
89 for (size_t r = 0; r < Rows(); ++r)
90 for (size_t c = 0; c < Cols(); ++c)
91 blocks_.push_back(MatrixFactory::Build(getRangeMap(r, bRangeThyraMode_), numEntriesPerRow));
92
93 // Default view
95}
96
97#ifdef HAVE_XPETRA_THYRA
98template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
99BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BlockedCrsMatrix(const Teuchos::RCP<const Thyra::BlockedLinearOpBase<Scalar>>& thyraOp, const Teuchos::RCP<const Teuchos::Comm<int>>& /* comm */)
100 : is_diagonal_(true)
101 , thyraOp_(thyraOp) {
102 // extract information from Thyra blocked operator and rebuilt information
103 const Teuchos::RCP<const Thyra::ProductVectorSpaceBase<Scalar>> productRangeSpace = thyraOp->productRange();
104 const Teuchos::RCP<const Thyra::ProductVectorSpaceBase<Scalar>> productDomainSpace = thyraOp->productDomain();
105
106 int numRangeBlocks = productRangeSpace->numBlocks();
107 int numDomainBlocks = productDomainSpace->numBlocks();
108
109 // build range map extractor from Thyra::BlockedLinearOpBase object
110 std::vector<Teuchos::RCP<const Map>> subRangeMaps(numRangeBlocks);
111 for (size_t r = 0; r < Teuchos::as<size_t>(numRangeBlocks); ++r) {
112 for (size_t c = 0; c < Teuchos::as<size_t>(numDomainBlocks); ++c) {
113 if (thyraOp->blockExists(r, c)) {
114 // we only need at least one block in each block row to extract the range map
115 Teuchos::RCP<const Thyra::LinearOpBase<Scalar>> const_op = thyraOp->getBlock(r, c); // nonConst access is not allowed.
117 Xpetra::ThyraUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::toXpetra(const_op);
118 subRangeMaps[r] = xop->getRangeMap();
119 if (r != c) is_diagonal_ = false;
120 break;
121 }
122 }
123 }
124 Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> fullRangeMap = mergeMaps(subRangeMaps);
125
126 // check whether the underlying Thyra operator uses Thyra-style numbering (default for most applications) or
127 // Xpetra style numbering
128 bool bRangeUseThyraStyleNumbering = false;
129 size_t numAllElements = 0;
130 for (size_t v = 0; v < subRangeMaps.size(); ++v) {
131 numAllElements += subRangeMaps[v]->getGlobalNumElements();
132 }
133 if (fullRangeMap->getGlobalNumElements() != numAllElements) bRangeUseThyraStyleNumbering = true;
134 rangemaps_ = Teuchos::rcp(new Xpetra::MapExtractor<Scalar, LocalOrdinal, GlobalOrdinal, Node>(fullRangeMap, subRangeMaps, bRangeUseThyraStyleNumbering));
135
136 // build domain map extractor from Thyra::BlockedLinearOpBase object
137 std::vector<Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>> subDomainMaps(numDomainBlocks);
138 for (size_t c = 0; c < Teuchos::as<size_t>(numDomainBlocks); ++c) {
139 for (size_t r = 0; r < Teuchos::as<size_t>(numRangeBlocks); ++r) {
140 if (thyraOp->blockExists(r, c)) {
141 // we only need at least one block in each block row to extract the range map
142 Teuchos::RCP<const Thyra::LinearOpBase<Scalar>> const_op = thyraOp->getBlock(r, c); // nonConst access is not allowed.
144 Xpetra::ThyraUtils<Scalar, LO, GO, Node>::toXpetra(const_op);
145 subDomainMaps[c] = xop->getDomainMap();
146 break;
147 }
148 }
149 }
150 Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> fullDomainMap = mergeMaps(subDomainMaps);
151 // plausibility check for numbering style (Xpetra or Thyra)
152 bool bDomainUseThyraStyleNumbering = false;
153 numAllElements = 0;
154 for (size_t v = 0; v < subDomainMaps.size(); ++v) {
155 numAllElements += subDomainMaps[v]->getGlobalNumElements();
156 }
157 if (fullDomainMap->getGlobalNumElements() != numAllElements) bDomainUseThyraStyleNumbering = true;
158 domainmaps_ = Teuchos::rcp(new Xpetra::MapExtractor<Scalar, LocalOrdinal, GlobalOrdinal, Node>(fullDomainMap, subDomainMaps, bDomainUseThyraStyleNumbering));
160 // store numbering mode
161 bRangeThyraMode_ = bRangeUseThyraStyleNumbering;
162 bDomainThyraMode_ = bDomainUseThyraStyleNumbering;
163
164 // add CrsMatrix objects in row,column order
165 blocks_.reserve(Rows() * Cols());
166 for (size_t r = 0; r < Rows(); ++r) {
167 for (size_t c = 0; c < Cols(); ++c) {
168 if (thyraOp->blockExists(r, c)) {
169 // TODO we do not support nested Thyra operators here!
170 Teuchos::RCP<const Thyra::LinearOpBase<Scalar>> const_op = thyraOp->getBlock(r, c); // nonConst access is not allowed.
171 Teuchos::RCP<Thyra::LinearOpBase<Scalar>> op = Teuchos::rcp_const_cast<Thyra::LinearOpBase<Scalar>>(const_op); // cast away const
173 Xpetra::ThyraUtils<Scalar, LO, GO, Node>::toXpetra(op);
175 Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LO, GO, Node>>(xop, true);
176 blocks_.push_back(xwrap);
177 } else {
178 // add empty block
179 blocks_.push_back(MatrixFactory::Build(getRangeMap(r, bRangeThyraMode_), 0));
180 }
181 }
182 }
183 // Default view
184 CreateDefaultView();
186
187template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
189 // TODO merging for Thyra mode is missing (similar to what we do in constructor of MapExtractor
190
191 // merge submaps to global map
192 std::vector<GlobalOrdinal> gids;
193 for (size_t tt = 0; tt < subMaps.size(); ++tt) {
195#if 1
196 Teuchos::ArrayView<const GlobalOrdinal> subMapGids = subMap->getLocalElementList();
197 gids.insert(gids.end(), subMapGids.begin(), subMapGids.end());
198#else
199 size_t myNumElements = subMap->getLocalNumElements();
200 for (LocalOrdinal l = 0; l < Teuchos::as<LocalOrdinal>(myNumElements); ++l) {
201 GlobalOrdinal gid = subMap->getGlobalElement(l);
202 gids.push_back(gid);
203 }
204#endif
205 }
206
207 // we have to sort the matrix entries and get rid of the double entries
208 // since we use this to detect Thyra-style numbering or Xpetra-style
209 // numbering. In Thyra-style numbering mode, the Xpetra::MapExtractor builds
210 // the correct row maps.
212 std::sort(gids.begin(), gids.end());
213 gids.erase(std::unique(gids.begin(), gids.end()), gids.end());
214 Teuchos::ArrayView<GO> gidsView(&gids[0], gids.size());
215 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());
216 return fullMap;
217}
218
219#endif
220
221template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
224template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
226 XPETRA_MONITOR("XpetraBlockedCrsMatrix::insertGlobalValues");
227 if (Rows() == 1 && Cols() == 1) {
228 getMatrix(0, 0)->insertGlobalValues(globalRow, cols, vals);
229 return;
230 }
231 throw Xpetra::Exceptions::RuntimeError("insertGlobalValues not supported by BlockedCrsMatrix");
232}
233
234template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
236 XPETRA_MONITOR("XpetraBlockedCrsMatrix::insertLocalValues");
237 if (Rows() == 1 && Cols() == 1) {
238 getMatrix(0, 0)->insertLocalValues(localRow, cols, vals);
239 return;
240 }
241 throw Xpetra::Exceptions::RuntimeError("insertLocalValues not supported by BlockedCrsMatrix");
242}
243
244template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246 XPETRA_MONITOR("XpetraBlockedCrsMatrix::removeEmptyProcessesInPlace");
247 if (Rows() == 1 && Cols() == 1) {
248 getMatrix(0, 0)->removeEmptyProcessesInPlace(newMap);
249 return;
250 }
251 throw Xpetra::Exceptions::RuntimeError("removeEmptyProcesses not supported by BlockedCrsMatrix");
252}
253
254template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
258 XPETRA_MONITOR("XpetraBlockedCrsMatrix::replaceGlobalValues");
259 if (Rows() == 1 && Cols() == 1) {
260 getMatrix(0, 0)->replaceGlobalValues(globalRow, cols, vals);
261 return;
262 }
263 throw Xpetra::Exceptions::RuntimeError("replaceGlobalValues not supported by BlockedCrsMatrix");
265
266template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
269 const ArrayView<const Scalar>& vals) {
270 XPETRA_MONITOR("XpetraBlockedCrsMatrix::replaceLocalValues");
271 if (Rows() == 1 && Cols() == 1) {
272 getMatrix(0, 0)->replaceLocalValues(localRow, cols, vals);
273 return;
274 }
275 throw Xpetra::Exceptions::RuntimeError("replaceLocalValues not supported by BlockedCrsMatrix");
276}
277
279template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
281 XPETRA_MONITOR("XpetraBlockedCrsMatrix::setAllToScalar");
282 for (size_t row = 0; row < Rows(); row++) {
283 for (size_t col = 0; col < Cols(); col++) {
284 if (!getMatrix(row, col).is_null()) {
285 getMatrix(row, col)->setAllToScalar(alpha);
286 }
287 }
288 }
289}
292template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
294 XPETRA_MONITOR("XpetraBlockedCrsMatrix::scale");
295 for (size_t row = 0; row < Rows(); row++) {
296 for (size_t col = 0; col < Cols(); col++) {
297 if (!getMatrix(row, col).is_null()) {
298 getMatrix(row, col)->scale(alpha);
299 }
300 }
301 }
302}
303
304template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
306 XPETRA_MONITOR("XpetraBlockedCrsMatrix::resumeFill");
307 for (size_t row = 0; row < Rows(); row++) {
308 for (size_t col = 0; col < Cols(); col++) {
309 if (!getMatrix(row, col).is_null()) {
310 getMatrix(row, col)->resumeFill(params);
311 }
312 }
313 }
314}
315
316template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
318 XPETRA_MONITOR("XpetraBlockedCrsMatrix::fillComplete");
319 if (Rows() == 1 && Cols() == 1) {
320 getMatrix(0, 0)->fillComplete(domainMap, rangeMap, params);
321 return;
322 }
323 fillComplete(params);
324}
326template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
328 XPETRA_MONITOR("XpetraBlockedCrsMatrix::fillComplete");
329 TEUCHOS_TEST_FOR_EXCEPTION(rangemaps_ == Teuchos::null, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::fillComplete: rangemaps_ is not set. Error.");
330
331 for (size_t r = 0; r < Rows(); ++r)
332 for (size_t c = 0; c < Cols(); ++c) {
333 if (getMatrix(r, c) != Teuchos::null) {
334 Teuchos::RCP<Matrix> Ablock = getMatrix(r, c);
335 if (r != c) is_diagonal_ = false;
336 if (Ablock != Teuchos::null && !Ablock->isFillComplete())
337 Ablock->fillComplete(getDomainMap(c, bDomainThyraMode_), getRangeMap(r, bRangeThyraMode_), params);
338 }
339 }
340
341#if 0
342 // get full row map
343 RCP<const Map> rangeMap = rangemaps_->getFullMap();
344 fullrowmap_ = MapFactory::Build(rangeMap()->lib(), rangeMap()->getGlobalNumElements(), rangeMap()->getLocalElementList(), rangeMap()->getIndexBase(), rangeMap()->getComm());
345
346 // build full col map
347 fullcolmap_ = Teuchos::null; // delete old full column map
348
349 std::vector<GO> colmapentries;
350 for (size_t c = 0; c < Cols(); ++c) {
351 // copy all local column lids of all block rows to colset
352 std::set<GO> colset;
353 for (size_t r = 0; r < Rows(); ++r) {
354 Teuchos::RCP<CrsMatrix> Ablock = getMatrix(r,c);
356 if (Ablock != Teuchos::null) {
357 Teuchos::ArrayView<const GO> colElements = Ablock->getColMap()->getLocalElementList();
358 Teuchos::RCP<const Map> colmap = Ablock->getColMap();
359 copy(colElements.getRawPtr(), colElements.getRawPtr() + colElements.size(), inserter(colset, colset.begin()));
360 }
361 }
362
363 // remove duplicates (entries which are in column maps of more than one block row)
364 colmapentries.reserve(colmapentries.size() + colset.size());
365 copy(colset.begin(), colset.end(), back_inserter(colmapentries));
366 sort(colmapentries.begin(), colmapentries.end());
367 typename std::vector<GO>::iterator gendLocation;
368 gendLocation = std::unique(colmapentries.begin(), colmapentries.end());
369 colmapentries.erase(gendLocation,colmapentries.end());
370 }
371
372 // sum up number of local elements
373 size_t numGlobalElements = 0;
374 Teuchos::reduceAll(*(rangeMap->getComm()), Teuchos::REDUCE_SUM, colmapentries.size(), Teuchos::outArg(numGlobalElements));
375
376 // store global full column map
378 fullcolmap_ = MapFactory::Build(rangeMap->lib(), numGlobalElements, aView, 0, rangeMap->getComm());
379#endif
380}
381
382template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
384 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalNumRows");
385 global_size_t globalNumRows = 0;
386
387 for (size_t row = 0; row < Rows(); row++)
388 for (size_t col = 0; col < Cols(); col++)
389 if (!getMatrix(row, col).is_null()) {
390 globalNumRows += getMatrix(row, col)->getGlobalNumRows();
391 break; // we need only one non-null matrix in a row
392 }
393
394 return globalNumRows;
395}
397template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
399 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalNumCols");
400 global_size_t globalNumCols = 0;
401
402 for (size_t col = 0; col < Cols(); col++)
403 for (size_t row = 0; row < Rows(); row++)
404 if (!getMatrix(row, col).is_null()) {
405 globalNumCols += getMatrix(row, col)->getGlobalNumCols();
406 break; // we need only one non-null matrix in a col
407 }
409 return globalNumCols;
410}
412template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
414 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalNumRows");
415 global_size_t nodeNumRows = 0;
416
417 for (size_t row = 0; row < Rows(); ++row)
418 for (size_t col = 0; col < Cols(); col++)
419 if (!getMatrix(row, col).is_null()) {
420 nodeNumRows += getMatrix(row, col)->getLocalNumRows();
421 break; // we need only one non-null matrix in a row
422 }
424 return nodeNumRows;
425}
427template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
429 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalNumEntries");
430 global_size_t globalNumEntries = 0;
431
432 for (size_t row = 0; row < Rows(); ++row)
433 for (size_t col = 0; col < Cols(); ++col)
434 if (!getMatrix(row, col).is_null())
435 globalNumEntries += getMatrix(row, col)->getGlobalNumEntries();
436
437 return globalNumEntries;
438}
439
440template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
442 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalNumEntries");
443 global_size_t nodeNumEntries = 0;
444
445 for (size_t row = 0; row < Rows(); ++row)
446 for (size_t col = 0; col < Cols(); ++col)
447 if (!getMatrix(row, col).is_null())
448 nodeNumEntries += getMatrix(row, col)->getLocalNumEntries();
449
450 return nodeNumEntries;
451}
452
453template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
455 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNumEntriesInLocalRow");
456 GlobalOrdinal gid = this->getRowMap()->getGlobalElement(localRow);
457 size_t row = getBlockedRangeMap()->getMapIndexForGID(gid);
458 LocalOrdinal lid = getBlockedRangeMap()->getMap(row)->getLocalElement(gid);
459 size_t numEntriesInLocalRow = 0;
460 for (size_t col = 0; col < Cols(); ++col)
461 if (!getMatrix(row, col).is_null())
462 numEntriesInLocalRow += getMatrix(row, col)->getNumEntriesInLocalRow(lid);
463 return numEntriesInLocalRow;
464}
466template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
468 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getNumEntriesInGlobalRow");
469 size_t row = getBlockedRangeMap()->getMapIndexForGID(globalRow);
470 size_t numEntriesInGlobalRow = 0;
471 for (size_t col = 0; col < Cols(); ++col)
472 if (!getMatrix(row, col).is_null())
473 numEntriesInGlobalRow += getMatrix(row, col)->getNumEntriesInGlobalRow(globalRow);
474 return numEntriesInGlobalRow;
475}
476
477template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
479 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalMaxNumRowEntries");
480 global_size_t globalMaxEntries = 0;
481
482 for (size_t row = 0; row < Rows(); row++) {
483 global_size_t globalMaxEntriesBlockRows = 0;
484 for (size_t col = 0; col < Cols(); col++) {
485 if (!getMatrix(row, col).is_null()) {
486 globalMaxEntriesBlockRows += getMatrix(row, col)->getGlobalMaxNumRowEntries();
487 }
489 if (globalMaxEntriesBlockRows > globalMaxEntries)
490 globalMaxEntries = globalMaxEntriesBlockRows;
491 }
492 return globalMaxEntries;
493}
494
495template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
497 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalMaxNumRowEntries");
498 size_t localMaxEntries = 0;
500 for (size_t row = 0; row < Rows(); row++) {
501 size_t localMaxEntriesBlockRows = 0;
502 for (size_t col = 0; col < Cols(); col++) {
503 if (!getMatrix(row, col).is_null()) {
504 localMaxEntriesBlockRows += getMatrix(row, col)->getLocalMaxNumRowEntries();
505 }
507 if (localMaxEntriesBlockRows > localMaxEntries)
508 localMaxEntries = localMaxEntriesBlockRows;
510 return localMaxEntries;
511}
513template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
515 XPETRA_MONITOR("XpetraBlockedCrsMatrix::isLocallyIndexed");
516 for (size_t i = 0; i < blocks_.size(); ++i)
517 if (blocks_[i] != Teuchos::null && !blocks_[i]->isLocallyIndexed())
518 return false;
519 return true;
520}
522template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
524 XPETRA_MONITOR("XpetraBlockedCrsMatrix::isGloballyIndexed");
525 for (size_t i = 0; i < blocks_.size(); i++)
526 if (blocks_[i] != Teuchos::null && !blocks_[i]->isGloballyIndexed())
527 return false;
528 return true;
529}
531template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
533 XPETRA_MONITOR("XpetraBlockedCrsMatrix::isFillComplete");
534 for (size_t i = 0; i < blocks_.size(); i++)
535 if (blocks_[i] != Teuchos::null && !blocks_[i]->isFillComplete())
536 return false;
537 return true;
538}
539
540template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
542 const ArrayView<LocalOrdinal>& Indices,
543 const ArrayView<Scalar>& Values,
544 size_t& NumEntries) const {
545 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalRowCopy");
546 if (Rows() == 1 && Cols() == 1) {
547 getMatrix(0, 0)->getLocalRowCopy(LocalRow, Indices, Values, NumEntries);
548 return;
549 }
550 throw Xpetra::Exceptions::RuntimeError("getLocalRowCopy not supported by BlockedCrsMatrix");
551}
552
553template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
555 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getGlobalRowView");
556 if (Rows() == 1 && Cols() == 1) {
557 getMatrix(0, 0)->getGlobalRowView(GlobalRow, indices, values);
558 return;
559 }
560 throw Xpetra::Exceptions::RuntimeError("getGlobalRowView not supported by BlockedCrsMatrix");
561}
562
563template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
565 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalRowView");
566 if (Rows() == 1 && Cols() == 1) {
567 getMatrix(0, 0)->getLocalRowView(LocalRow, indices, values);
568 return;
569 } else if (is_diagonal_) {
570 GlobalOrdinal gid = this->getRowMap()->getGlobalElement(LocalRow);
571 size_t row = getBlockedRangeMap()->getMapIndexForGID(gid);
572 getMatrix(row, row)->getLocalRowView(getMatrix(row, row)->getRowMap()->getLocalElement(gid), indices, values);
573 return;
574 }
575 throw Xpetra::Exceptions::RuntimeError("getLocalRowView not supported by BlockedCrsMatrix");
576}
577
578template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
580 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getLocalDiagCopy");
581
582 // RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
583
584 Teuchos::RCP<Vector> rcpdiag = Teuchos::rcpFromRef(diag);
585 Teuchos::RCP<BlockedVector> bdiag = Teuchos::rcp_dynamic_cast<BlockedVector>(rcpdiag);
586
587 // special treatment for 1x1 block matrices
588 // ReorderedBlockedCrsMatrix object encapsulate single blocks in ReorderedBlocks
589 // BlockedVectors have Vector objects as Leaf objects.
590 if (Rows() == 1 && Cols() == 1 && bdiag.is_null() == true) {
591 Teuchos::RCP<const Matrix> rm = getMatrix(0, 0);
592 rm->getLocalDiagCopy(diag);
593 return;
594 }
595
596 TEUCHOS_TEST_FOR_EXCEPTION(bdiag.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::getLocalDiagCopy: diag must be a Blocked(Multi)Vector.");
597 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());
598 TEUCHOS_TEST_FOR_EXCEPTION(bdiag->getBlockedMap()->getNumMaps() != this->Rows(), Xpetra::Exceptions::RuntimeError,
599 "BlockedCrsMatrix::getLocalDiagCopy(): the number of blocks in diag differ from the number of blocks in this operator.");
600 // XPETRA_TEST_FOR_EXCEPTION(bdiag->getMap()->isSameAs(*(getMap())) == false, Xpetra::Exceptions::RuntimeError,
601 // "BlockedCrsMatrix::getLocalDiagCopy(): the map of the vector diag is not compatible with the map of the blocked operator." );
602
603 for (size_t row = 0; row < Rows(); row++) {
604 Teuchos::RCP<const Matrix> rm = getMatrix(row, row);
605 if (!rm.is_null()) {
606 Teuchos::RCP<Vector> rv = VectorFactory::Build(bdiag->getBlockedMap()->getMap(row, bdiag->getBlockedMap()->getThyraMode()));
607 rm->getLocalDiagCopy(*rv);
608 bdiag->setMultiVector(row, rv, bdiag->getBlockedMap()->getThyraMode());
609 }
610 }
611}
612
613template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
615 XPETRA_MONITOR("XpetraBlockedCrsMatrix::leftScale");
616
617 Teuchos::RCP<const Vector> rcpx = Teuchos::rcpFromRef(x);
618 Teuchos::RCP<const BlockedVector> bx = Teuchos::rcp_dynamic_cast<const BlockedVector>(rcpx);
619
620 // RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
621
622 // special treatment for 1xn block matrices
623 // ReorderedBlockedCrsMatrix object encapsulate single blocks in ReorderedBlocks
624 // BlockedVectors have Vector objects as Leaf objects.
625 if (Rows() == 1 && bx.is_null() == true) {
626 for (size_t col = 0; col < Cols(); ++col) {
627 Teuchos::RCP<Matrix> rm = getMatrix(0, col);
628 rm->leftScale(x);
629 }
630 return;
631 }
632
633 TEUCHOS_TEST_FOR_EXCEPTION(bx.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::leftScale: x must be a Blocked(Multi)Vector.");
634 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());
635 TEUCHOS_TEST_FOR_EXCEPTION(bx->getBlockedMap()->getNumMaps() != this->Rows(), Xpetra::Exceptions::RuntimeError,
636 "BlockedCrsMatrix::leftScale(): the number of blocks in diag differ from the number of blocks in this operator.");
637
638 for (size_t row = 0; row < Rows(); row++) {
639 Teuchos::RCP<const MultiVector> rmv = bx->getMultiVector(row);
640 Teuchos::RCP<const Vector> rscale = rmv->getVector(0);
641 XPETRA_TEST_FOR_EXCEPTION(rscale.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::leftScale: x must be a Vector.");
642 for (size_t col = 0; col < Cols(); ++col) {
643 Teuchos::RCP<Matrix> rm = getMatrix(row, col);
644 if (!rm.is_null()) {
645 rm->leftScale(*rscale);
646 }
647 }
648 }
649}
650
651template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
653 XPETRA_MONITOR("XpetraBlockedCrsMatrix::rightScale");
654
655 Teuchos::RCP<const Vector> rcpx = Teuchos::rcpFromRef(x);
656 Teuchos::RCP<const BlockedVector> bx = Teuchos::rcp_dynamic_cast<const BlockedVector>(rcpx);
657
658 // RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
659
660 // special treatment for nx1 block matrices
661 // ReorderedBlockedCrsMatrix object encapsulate single blocks in ReorderedBlocks
662 // BlockedVectors have Vector objects as Leaf objects.
663 if (Cols() == 1 && bx.is_null() == true) {
664 for (size_t row = 0; row < Rows(); ++row) {
665 Teuchos::RCP<Matrix> rm = getMatrix(row, 0);
666 rm->rightScale(x);
667 }
668 return;
669 }
670
671 TEUCHOS_TEST_FOR_EXCEPTION(bx.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::rightScale: x must be a Blocked(Multi)Vector.");
672 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());
673 TEUCHOS_TEST_FOR_EXCEPTION(bx->getBlockedMap()->getNumMaps() != this->Cols(), Xpetra::Exceptions::RuntimeError,
674 "BlockedCrsMatrix::rightScale(): the number of blocks in diag differ from the number of blocks in this operator.");
675
676 for (size_t col = 0; col < Cols(); ++col) {
677 Teuchos::RCP<const MultiVector> rmv = bx->getMultiVector(col);
678 Teuchos::RCP<const Vector> rscale = rmv->getVector(0);
679 XPETRA_TEST_FOR_EXCEPTION(rscale.is_null() == true, Xpetra::Exceptions::RuntimeError, "BlockedCrsMatrix::leftScale: x must be a Vector.");
680 for (size_t row = 0; row < Rows(); row++) {
681 Teuchos::RCP<Matrix> rm = getMatrix(row, col);
682 if (!rm.is_null()) {
683 rm->rightScale(*rscale);
684 }
685 }
686 }
687}
688
689template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
691 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getFrobeniusNorm");
693 for (size_t col = 0; col < Cols(); ++col) {
694 for (size_t row = 0; row < Rows(); ++row) {
695 if (getMatrix(row, col) != Teuchos::null) {
696 typename ScalarTraits<Scalar>::magnitudeType n = getMatrix(row, col)->getFrobeniusNorm();
697 ret += n * n;
698 }
699 }
700 }
702}
703
704template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
706
707template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
708void BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::apply(const MultiVector& X, MultiVector& Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues,
709 const RCP<Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node>>& regionInterfaceImporter,
710 const Teuchos::ArrayRCP<LocalOrdinal>& regionInterfaceLIDs) const {}
711
712template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
714 Teuchos::ETransp mode,
715 Scalar alpha,
716 Scalar beta) const {
717 XPETRA_MONITOR("XpetraBlockedCrsMatrix::apply");
718 // using Teuchos::RCP;
719
721 "apply() only supports the following modes: NO_TRANS and TRANS.");
722
723 // check whether input parameters are blocked or not
724 RCP<const MultiVector> refX = rcpFromRef(X);
725 RCP<const BlockedMultiVector> refbX = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(refX);
726 // RCP<MultiVector> tmpY = rcpFromRef(Y);
727 // RCP<BlockedMultiVector> tmpbY = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(tmpY);
728
729 // TODO get rid of me: adapt MapExtractor
730 bool bBlockedX = (refbX != Teuchos::null) ? true : false;
731
732 // create (temporary) vectors for output
733 // In the end we call Y.update(alpha, *tmpY, beta). Therefore we need a new vector storing the temporary results
734 RCP<MultiVector> tmpY = MultiVectorFactory::Build(Y.getMap(), Y.getNumVectors(), true);
735
736 // RCP<Teuchos::FancyOStream> out = rcp(new Teuchos::FancyOStream(rcp(&std::cout,false)));
737
739
740 if (mode == Teuchos::NO_TRANS) {
741 for (size_t row = 0; row < Rows(); row++) {
742 RCP<MultiVector> Yblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, true);
743 for (size_t col = 0; col < Cols(); col++) {
744 // extract matrix block
745 RCP<Matrix> Ablock = getMatrix(row, col);
746
747 if (Ablock.is_null())
748 continue;
749
750 // check whether Ablock is itself a blocked operator
751 // If it is a blocked operator we have to provide Xpetra style GIDs, i.e. we have to transform GIDs
752 bool bBlockedSubMatrix = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Ablock) == Teuchos::null ? false : true;
753
754 // input/output vectors for local block operation
755 RCP<const MultiVector> Xblock = Teuchos::null; // subpart of X vector to be applied to subblock of A
756
757 // extract sub part of X using Xpetra or Thyra GIDs
758 // if submatrix is again blocked, we extract it using Xpetra style gids. If it is a single
759 // block matrix we use the Thyra or Xpetra style GIDs that are used to store the matrix
760 if (bBlockedX)
761 Xblock = domainmaps_->ExtractVector(refbX, col, bDomainThyraMode_);
762 else
763 Xblock = domainmaps_->ExtractVector(refX, col, bBlockedSubMatrix == true ? false : bDomainThyraMode_);
764
765 RCP<MultiVector> tmpYblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, false); // subpart of Y vector containing part of solution of Xblock applied to Ablock
766 Ablock->apply(*Xblock, *tmpYblock);
767
768 Yblock->update(one, *tmpYblock, one);
769 }
770 rangemaps_->InsertVector(Yblock, row, tmpY, bRangeThyraMode_);
771 }
772
773 } else if (mode == Teuchos::TRANS) {
774 // TODO: test me!
775 for (size_t col = 0; col < Cols(); col++) {
776 RCP<MultiVector> Yblock = domainmaps_->getVector(col, Y.getNumVectors(), bDomainThyraMode_, true);
777
778 for (size_t row = 0; row < Rows(); row++) {
779 RCP<Matrix> Ablock = getMatrix(row, col);
780
781 if (Ablock.is_null())
782 continue;
783
784 // check whether Ablock is itself a blocked operator
785 // If it is a blocked operator we have to provide Xpetra style GIDs, i.e. we have to transform GIDs
786 bool bBlockedSubMatrix = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Ablock) == Teuchos::null ? false : true;
787
788 RCP<const MultiVector> Xblock = Teuchos::null;
789
790 // extract sub part of X using Xpetra or Thyra GIDs
791 if (bBlockedX)
792 Xblock = rangemaps_->ExtractVector(refbX, row, bRangeThyraMode_);
793 else
794 Xblock = rangemaps_->ExtractVector(refX, row, bBlockedSubMatrix == true ? false : bRangeThyraMode_);
795 RCP<MultiVector> tmpYblock = domainmaps_->getVector(col, Y.getNumVectors(), bDomainThyraMode_, false);
796 Ablock->apply(*Xblock, *tmpYblock, Teuchos::TRANS);
797
798 Yblock->update(one, *tmpYblock, one);
799 }
800 domainmaps_->InsertVector(Yblock, col, tmpY, bDomainThyraMode_);
801 }
802 }
803 Y.update(alpha, *tmpY, beta);
804}
805
806template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
808 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getFullDomainMap()");
809 return domainmaps_->getFullMap();
810}
811
812template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
814 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getBlockedDomainMap()");
815 return domainmaps_->getBlockedMap();
816}
817
818template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
820 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMap()");
821 return domainmaps_->getMap(); /*domainmaps_->getFullMap();*/
822}
823
824template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
826 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMap(size_t)");
827 return domainmaps_->getMap(i, bDomainThyraMode_);
828}
829
830template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
832 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getDomainMap(size_t,bool)");
833 return domainmaps_->getMap(i, bThyraMode);
834}
835
836template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
838 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap()");
839 return rangemaps_->getFullMap();
840}
841
842template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
844 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getBlockedRangeMap()");
845 return rangemaps_->getBlockedMap();
846}
847
848template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
850 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap()");
851 return rangemaps_->getMap(); /*rangemaps_->getFullMap();*/
852}
853
854template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
856 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap(size_t)");
857 return rangemaps_->getMap(i, bRangeThyraMode_);
858}
859
860template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
862 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getRangeMap(size_t,bool)");
863 return rangemaps_->getMap(i, bThyraMode);
864}
865
866template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
871
872template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
877
878template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
880 const MultiVector& X,
881 MultiVector& Y,
882 size_t row,
883 Teuchos::ETransp mode,
884 Scalar alpha,
885 Scalar beta
886) const {
887 XPETRA_MONITOR("XpetraBlockedCrsMatrix::bgs_apply");
888 // using Teuchos::RCP;
889
891 "apply() only supports the following modes: NO_TRANS and TRANS.");
892
893 // check whether input parameters are blocked or not
894 RCP<const MultiVector> refX = rcpFromRef(X);
895 RCP<const BlockedMultiVector> refbX = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(refX);
896 // RCP<MultiVector> tmpY = rcpFromRef(Y);
897 // RCP<BlockedMultiVector> tmpbY = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(tmpY);
898
899 bool bBlockedX = (refbX != Teuchos::null) ? true : false;
900
901 // create (temporary) vectors for output
902 // In the end we call Y.update(alpha, *tmpY, beta). Therefore we need a new vector storing the temporary results
903 RCP<MultiVector> tmpY = MultiVectorFactory::Build(Y.getMap(), Y.getNumVectors(), true);
904
906
907 if (mode == Teuchos::NO_TRANS) {
908 RCP<MultiVector> Yblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, true);
909 for (size_t col = 0; col < Cols(); col++) {
910 // extract matrix block
911 RCP<Matrix> Ablock = getMatrix(row, col);
912
913 if (Ablock.is_null())
914 continue;
915
916 // check whether Ablock is itself a blocked operator
917 // If it is a blocked operator we have to provide Xpetra style GIDs, i.e. we have to transform GIDs
918 bool bBlockedSubMatrix = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Ablock) == Teuchos::null ? false : true;
919
920 // input/output vectors for local block operation
921 RCP<const MultiVector> Xblock = Teuchos::null; // subpart of X vector to be applied to subblock of A
922
923 // extract sub part of X using Xpetra or Thyra GIDs
924 // if submatrix is again blocked, we extract it using Xpetra style gids. If it is a single
925 // block matrix we use the Thyra or Xpetra style GIDs that are used to store the matrix
926 if (bBlockedX)
927 Xblock = domainmaps_->ExtractVector(refbX, col, bDomainThyraMode_);
928 else
929 Xblock = domainmaps_->ExtractVector(refX, col, bBlockedSubMatrix == true ? false : bDomainThyraMode_);
930
931 RCP<MultiVector> tmpYblock = rangemaps_->getVector(row, Y.getNumVectors(), bRangeThyraMode_, false); // subpart of Y vector containing part of solution of Xblock applied to Ablock
932 Ablock->apply(*Xblock, *tmpYblock);
933
934 Yblock->update(one, *tmpYblock, one);
935 }
936 rangemaps_->InsertVector(Yblock, row, tmpY, bRangeThyraMode_);
937 } else {
938 TEUCHOS_TEST_FOR_EXCEPTION(true, Xpetra::Exceptions::NotImplemented, "Xpetar::BlockedCrsMatrix::bgs_apply: not implemented for transpose case.");
939 }
940 Y.update(alpha, *tmpY, beta);
941}
942
943template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
945 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getMap");
946 if (Rows() == 1 && Cols() == 1) {
947 return getMatrix(0, 0)->getMap();
948 }
949 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::getMap(): operation not supported.");
950}
951
952template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
954 XPETRA_MONITOR("XpetraBlockedCrsMatrix::doImport");
955 if (Rows() == 1 && Cols() == 1) {
956 getMatrix(0, 0)->doImport(source, importer, CM);
957 return;
958 }
959 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doImport(): operation not supported.");
960}
961
962template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
964 XPETRA_MONITOR("XpetraBlockedCrsMatrix::doExport");
965 if (Rows() == 1 && Cols() == 1) {
966 getMatrix(0, 0)->doExport(dest, importer, CM);
967 return;
968 }
969 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doExport(): operation not supported.");
970}
971
972template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
974 XPETRA_MONITOR("XpetraBlockedCrsMatrix::doImport");
975 if (Rows() == 1 && Cols() == 1) {
976 getMatrix(0, 0)->doImport(source, exporter, CM);
977 return;
978 }
979 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doImport(): operation not supported.");
980}
981
983template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
985 XPETRA_MONITOR("XpetraBlockedCrsMatrix::doExport");
986 if (Rows() == 1 && Cols() == 1) {
987 getMatrix(0, 0)->doExport(dest, exporter, CM);
988 return;
989 }
990 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::doExport(): operation not supported.");
991}
992
993template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
994std::string BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::description() const { return "Xpetra_BlockedCrsMatrix.description()"; }
995
996template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
998 out << "Xpetra::BlockedCrsMatrix: " << Rows() << " x " << Cols() << std::endl;
999
1000 if (isFillComplete()) {
1001 out << "BlockMatrix is fillComplete" << std::endl;
1002
1003 /*if(fullrowmap_ != Teuchos::null) {
1004 out << "fullRowMap" << std::endl;
1005 fullrowmap_->describe(out,verbLevel);
1006 } else {
1007 out << "fullRowMap not set. Check whether block matrix is properly fillCompleted!" << std::endl;
1008 }*/
1009
1010 // out << "fullColMap" << std::endl;
1011 // fullcolmap_->describe(out,verbLevel);
1012
1013 } else {
1014 out << "BlockMatrix is NOT fillComplete" << std::endl;
1015 }
1016
1017 for (size_t r = 0; r < Rows(); ++r)
1018 for (size_t c = 0; c < Cols(); ++c) {
1019 if (getMatrix(r, c) != Teuchos::null) {
1020 out << "Block(" << r << "," << c << ")" << std::endl;
1021 getMatrix(r, c)->describe(out, verbLevel);
1022 } else
1023 out << "Block(" << r << "," << c << ") = null" << std::endl;
1024 }
1025}
1026
1027template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1029 XPETRA_MONITOR("TpetraBlockedCrsMatrix::setObjectLabel");
1030 for (size_t r = 0; r < Rows(); ++r)
1031 for (size_t c = 0; c < Cols(); ++c) {
1032 if (getMatrix(r, c) != Teuchos::null) {
1033 std::ostringstream oss;
1034 oss << objectLabel << "(" << r << "," << c << ")";
1035 getMatrix(r, c)->setObjectLabel(oss.str());
1036 }
1037 }
1038}
1039
1040template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1042 if (Rows() == 1 && Cols() == 1)
1043 return true;
1044 else
1045 return false;
1046}
1047
1048template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1050 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getCrsGraph");
1051 if (Rows() == 1 && Cols() == 1) {
1052 return getMatrix(0, 0)->getCrsGraph();
1053 }
1054 throw Xpetra::Exceptions::RuntimeError("getCrsGraph() not supported by BlockedCrsMatrix");
1055}
1056
1057template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1059
1060template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1062 XPETRA_MONITOR("XpetraBlockedCrsMatrix::Rows");
1063 return rangemaps_->NumMaps();
1064}
1065
1066template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1068 XPETRA_MONITOR("XpetraBlockedCrsMatrix::Cols");
1069 return domainmaps_->NumMaps();
1070}
1071
1072template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1074 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getCrsMatrix");
1075 TEUCHOS_TEST_FOR_EXCEPTION(Rows() != 1, std::out_of_range, "Can only unwrap a 1x1 blocked matrix. The matrix has " << Rows() << " block rows, though.");
1076 TEUCHOS_TEST_FOR_EXCEPTION(Cols() != 1, std::out_of_range, "Can only unwrap a 1x1 blocked matrix. The matrix has " << Cols() << " block columns, though.");
1077
1078 RCP<Matrix> mat = getMatrix(0, 0);
1079 RCP<BlockedCrsMatrix> bmat = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(mat);
1080 if (bmat == Teuchos::null) return mat;
1081 return bmat->getCrsMatrix();
1082}
1083
1084template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1086 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getInnermostCrsMatrix");
1087 size_t row = Rows() + 1, col = Cols() + 1;
1088 for (size_t r = 0; r < Rows(); ++r)
1089 for (size_t c = 0; c < Cols(); ++c)
1090 if (getMatrix(r, c) != Teuchos::null) {
1091 row = r;
1092 col = c;
1093 break;
1094 }
1095 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.")
1096 RCP<Matrix> mm = getMatrix(row, col);
1097 RCP<BlockedCrsMatrix> bmat = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(mm);
1098 if (bmat == Teuchos::null) return mm;
1099 return bmat->getInnermostCrsMatrix();
1100}
1101
1102template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1104 XPETRA_MONITOR("XpetraBlockedCrsMatrix::getMatrix");
1105 TEUCHOS_TEST_FOR_EXCEPTION(r > Rows(), std::out_of_range, "Error, r = " << Rows() << " is too big");
1106 TEUCHOS_TEST_FOR_EXCEPTION(c > Cols(), std::out_of_range, "Error, c = " << Cols() << " is too big");
1107
1108 // transfer strided/blocked map information
1109 /* if (blocks_[r*Cols()+c] != Teuchos::null &&
1110 blocks_[r*Cols()+c]->IsView("stridedMaps") == false)
1111 blocks_[r*Cols()+c]->CreateView("stridedMaps", getRangeMap(r,bRangeThyraMode_), getDomainMap(c,bDomainThyraMode_));*/
1112 return blocks_[r * Cols() + c];
1113}
1114
1115template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1117 XPETRA_MONITOR("XpetraBlockedCrsMatrix::setMatrix");
1118 // TODO: if filled -> return error
1119
1120 TEUCHOS_TEST_FOR_EXCEPTION(r > Rows(), std::out_of_range, "Error, r = " << Rows() << " is too big");
1121 TEUCHOS_TEST_FOR_EXCEPTION(c > Cols(), std::out_of_range, "Error, c = " << Cols() << " is too big");
1122 if (!mat.is_null() && r != c) is_diagonal_ = false;
1123 // set matrix
1124 blocks_[r * Cols() + c] = mat;
1125}
1126
1127template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1129 XPETRA_MONITOR("XpetraBlockedCrsMatrix::Merge");
1130 using Teuchos::RCP;
1131 using Teuchos::rcp_dynamic_cast;
1132 Scalar one = ScalarTraits<SC>::one();
1133
1134 TEUCHOS_TEST_FOR_EXCEPTION(bRangeThyraMode_ != bDomainThyraMode_, Xpetra::Exceptions::RuntimeError,
1135 "BlockedCrsMatrix::Merge: only implemented for Xpetra-style or Thyra-style numbering. No mixup allowed!");
1136
1138 "BlockedCrsMatrix::Merge: BlockMatrix must be fill-completed.");
1139
1140 LocalOrdinal lclNumRows = getFullRangeMap()->getLocalNumElements();
1141 Teuchos::ArrayRCP<size_t> numEntPerRow(lclNumRows);
1142 for (LocalOrdinal lclRow = 0; lclRow < lclNumRows; ++lclRow)
1143 numEntPerRow[lclRow] = getNumEntriesInLocalRow(lclRow);
1144
1145 RCP<Matrix> sparse = MatrixFactory::Build(getFullRangeMap(), numEntPerRow);
1146
1147 if (bRangeThyraMode_ == false) {
1148 // Xpetra mode
1149 for (size_t i = 0; i < Rows(); i++) {
1150 for (size_t j = 0; j < Cols(); j++) {
1151 if (getMatrix(i, j) != Teuchos::null) {
1152 RCP<const Matrix> mat = getMatrix(i, j);
1153
1154 // recursively call Merge routine
1155 RCP<const BlockedCrsMatrix> bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1156 if (bMat != Teuchos::null) mat = bMat->Merge();
1157
1158 bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1160 "BlockedCrsMatrix::Merge: Merging of blocked sub-operators failed?!");
1161
1162 // jump over empty blocks
1163 if (mat->getLocalNumEntries() == 0) continue;
1164
1165 this->Add(*mat, one, *sparse, one);
1166 }
1167 }
1168 }
1169 } else {
1170 // Thyra mode
1171 for (size_t i = 0; i < Rows(); i++) {
1172 for (size_t j = 0; j < Cols(); j++) {
1173 if (getMatrix(i, j) != Teuchos::null) {
1174 RCP<const Matrix> mat = getMatrix(i, j);
1175 // recursively call Merge routine
1176 RCP<const BlockedCrsMatrix> bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1177 if (bMat != Teuchos::null) mat = bMat->Merge();
1178
1179 bMat = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(mat);
1181 "BlockedCrsMatrix::Merge: Merging of blocked sub-operators failed?!");
1182
1183 // check whether we have a CrsMatrix block (no blocked operator)
1184 RCP<const CrsMatrixWrap> crsMat = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(mat);
1185 TEUCHOS_ASSERT(crsMat != Teuchos::null);
1186
1187 // these are the thyra style maps of the matrix
1188 RCP<const Map> trowMap = mat->getRowMap();
1189 RCP<const Map> tcolMap = mat->getColMap();
1190 RCP<const Map> tdomMap = mat->getDomainMap();
1191
1192 // get Xpetra maps
1193 RCP<const Map> xrowMap = getRangeMapExtractor()->getMap(i, false);
1194 RCP<const Map> xdomMap = getDomainMapExtractor()->getMap(j, false);
1195
1196 // generate column map with Xpetra GIDs
1197 // We have to do this separately for each block since the column
1198 // map of each block might be different in the same block column
1199 Teuchos::RCP<Map> xcolMap = MapUtils::transformThyra2XpetraGIDs(
1200 *tcolMap,
1201 *tdomMap,
1202 *xdomMap);
1203
1204 // jump over empty blocks
1205 if (mat->getLocalNumEntries() == 0) continue;
1206
1207 size_t maxNumEntries = mat->getLocalMaxNumRowEntries();
1208
1209 size_t numEntries;
1210 Array<GO> inds(maxNumEntries);
1211 Array<GO> inds2(maxNumEntries);
1212 Array<SC> vals(maxNumEntries);
1213
1214 // loop over all rows and add entries
1215 for (size_t k = 0; k < mat->getLocalNumRows(); k++) {
1216 GlobalOrdinal rowTGID = trowMap->getGlobalElement(k);
1217 crsMat->getCrsMatrix()->getGlobalRowCopy(rowTGID, inds(), vals(), numEntries);
1218
1219 // create new indices array
1220 for (size_t l = 0; l < numEntries; ++l) {
1221 LocalOrdinal lid = tcolMap->getLocalElement(inds[l]);
1222 inds2[l] = xcolMap->getGlobalElement(lid);
1223 }
1224
1225 GlobalOrdinal rowXGID = xrowMap->getGlobalElement(k);
1226 sparse->insertGlobalValues(
1227 rowXGID, inds2(0, numEntries),
1228 vals(0, numEntries));
1229 }
1230 }
1231 }
1232 }
1233 }
1234
1235 sparse->fillComplete(getFullDomainMap(), getFullRangeMap());
1236
1237 TEUCHOS_TEST_FOR_EXCEPTION(sparse->getLocalNumEntries() != getLocalNumEntries(), Xpetra::Exceptions::RuntimeError,
1238 "BlockedCrsMatrix::Merge: Local number of entries of merged matrix does not coincide with local number of entries of blocked operator.");
1239
1240 TEUCHOS_TEST_FOR_EXCEPTION(sparse->getGlobalNumEntries() != getGlobalNumEntries(), Xpetra::Exceptions::RuntimeError,
1241 "BlockedCrsMatrix::Merge: Global number of entries of merged matrix does not coincide with global number of entries of blocked operator.");
1242
1243 return sparse;
1244}
1245
1246template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1248 if (Rows() == 1 && Cols() == 1) {
1249 return getMatrix(0, 0)->getLocalMatrixDevice();
1250 }
1251 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::getLocalMatrix(): operation not supported.");
1252}
1253
1254template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1255#if KOKKOS_VERSION >= 40799
1257#else
1259#endif
1260 if (Rows() == 1 && Cols() == 1) {
1261 return getMatrix(0, 0)->getLocalMatrixHost();
1262 }
1263 throw Xpetra::Exceptions::RuntimeError("BlockedCrsMatrix::getLocalMatrix(): operation not supported.");
1264}
1265
1266#ifdef HAVE_XPETRA_THYRA
1267template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1270 Xpetra::ThyraUtils<Scalar, LO, GO, Node>::toThyra(Teuchos::rcpFromRef(*this));
1272
1274 Teuchos::rcp_dynamic_cast<Thyra::BlockedLinearOpBase<Scalar>>(thOp);
1276 return thbOp;
1277}
1278#endif
1279
1280template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1282
1283template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1285 const MultiVector& B,
1286 MultiVector& R) const {
1288 R.update(STS::one(), B, STS::zero());
1289 this->apply(X, R, Teuchos::NO_TRANS, -STS::one(), STS::one());
1290}
1291
1292template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1293void BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Add(const Matrix& A, const Scalar scalarA, Matrix& B, const Scalar scalarB) const {
1294 XPETRA_MONITOR("XpetraBlockedCrsMatrix::Add");
1296 "Matrix A is not completed");
1297 using Teuchos::Array;
1298 using Teuchos::ArrayView;
1299
1300 B.scale(scalarB);
1301
1302 Scalar one = ScalarTraits<SC>::one();
1303 Scalar zero = ScalarTraits<SC>::zero();
1304
1305 if (scalarA == zero)
1306 return;
1307
1308 Teuchos::RCP<const Matrix> rcpA = Teuchos::rcpFromRef(A);
1309 Teuchos::RCP<const CrsMatrixWrap> rcpAwrap = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(rcpA);
1311 "BlockedCrsMatrix::Add: matrix A must be of type CrsMatrixWrap.");
1312 Teuchos::RCP<const CrsMatrix> crsA = rcpAwrap->getCrsMatrix();
1313
1314 size_t maxNumEntries = crsA->getLocalMaxNumRowEntries();
1315
1316 size_t numEntries;
1317 Array<GO> inds(maxNumEntries);
1318 Array<SC> vals(maxNumEntries);
1319
1320 RCP<const Map> rowMap = crsA->getRowMap();
1321 RCP<const Map> colMap = crsA->getColMap();
1322
1323 ArrayView<const GO> rowGIDs = crsA->getRowMap()->getLocalElementList();
1324 for (size_t i = 0; i < crsA->getLocalNumRows(); i++) {
1325 GO row = rowGIDs[i];
1326 crsA->getGlobalRowCopy(row, inds(), vals(), numEntries);
1327
1328 if (scalarA != one)
1329 for (size_t j = 0; j < numEntries; ++j)
1330 vals[j] *= scalarA;
1331
1332 B.insertGlobalValues(row, inds(0, numEntries), vals(0, numEntries)); // insert should be ok, since blocks in BlockedCrsOpeartor do not overlap!
1333 }
1334}
1335
1336template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1338 XPETRA_MONITOR("XpetraBlockedCrsMatrix::CreateDefaultView");
1339
1340 // Create default view
1341 this->defaultViewLabel_ = "point";
1342 this->CreateView(this->GetDefaultViewLabel(), getRangeMap(), getDomainMap());
1343
1344 // Set current view
1345 this->currentViewLabel_ = this->GetDefaultViewLabel();
1346}
1347
1348} // namespace Xpetra
1349
1350#define XPETRA_BLOCKEDCRSMATRIX_SHORT
1351#endif /* XPETRA_BLOCKEDCRSMATRIX_DEF_HPP */
Add
#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
bool is_null() const
Teuchos::RCP< const MapExtractor > domainmaps_
full domain map together with all partial domain maps
virtual size_t Cols() const
number of column blocks
Teuchos::RCP< const MapExtractor > rangemaps_
full range map together with all partial domain maps
const RCP< const Map > getRangeMap() const
Returns the Map associated with the range of this operator.
BlockedCrsMatrix(const Teuchos::RCP< const BlockedMap > &rangeMaps, const Teuchos::RCP< const BlockedMap > &domainMaps, size_t numEntriesPerRow)
Constructor.
bool bDomainThyraMode_
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub...
std::vector< Teuchos::RCP< Matrix > > blocks_
row major matrix block storage
local_matrix_type::HostMirror getLocalMatrixHost() const
Access the underlying local Kokkos::CrsMatrix object.
bool bRangeThyraMode_
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub...
virtual size_t Rows() const
number of row blocks
local_matrix_type getLocalMatrixDevice() const
Access the underlying local Kokkos::CrsMatrix object.
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, execution_space, void, typename local_graph_type::size_type > local_matrix_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
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 RCP< Matrix > Build(const RCP< const Map > &rowMap)
Xpetra-specific matrix class.
virtual bool isFillComplete() const =0
Returns true if fillComplete() has been called and the matrix is in compute mode.
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.
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.
#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)
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.
static magnitudeType magnitude(T a)