47#ifndef PACKAGES_XPETRA_SUP_MATRIX_UTILS_HPP_
48#define PACKAGES_XPETRA_SUP_MATRIX_UTILS_HPP_
53#include "Xpetra_Map.hpp"
55#include "Xpetra_StridedMap.hpp"
56#include "Xpetra_MapFactory.hpp"
57#include "Xpetra_MapExtractor.hpp"
66#ifdef HAVE_XPETRA_TPETRA
67#include "Xpetra_TpetraMultiVector.hpp"
68#include <Tpetra_RowMatrixTransposer.hpp>
69#include <Tpetra_Details_extractBlockDiagonal.hpp>
70#include <Tpetra_Details_scaleBlockDiagonal.hpp>
84template <
class Scalar,
89#undef XPETRA_MATRIXUTILS_SHORT
101 ret->replaceLocalValue(Teuchos::as<LocalOrdinal>(r), c, data[r]);
119 for(
size_t i = 0; i<input.
getColMap()->getNodeNumElements(); i++) {
120 GlobalOrdinal gcid = input.
getColMap()->getGlobalElement(i);
128 std::vector<int> myUnknownDofGIDs(comm->getSize(),0);
129 std::vector<int> numUnknownDofGIDs(comm->getSize(),0);
130 myUnknownDofGIDs[comm->getRank()] = ovlUnknownStatusGids.
size();
131 Teuchos::reduceAll(*comm,Teuchos::REDUCE_MAX,comm->getSize(),&myUnknownDofGIDs[0],&numUnknownDofGIDs[0]);
134 size_t cntUnknownDofGIDs = 0;
135 for(
int p = 0; p < comm->getSize(); p++) cntUnknownDofGIDs += numUnknownDofGIDs[p];
136 std::vector<GlobalOrdinal> lUnknownDofGIDs(cntUnknownDofGIDs,-1);
137 std::vector<GlobalOrdinal> gUnknownDofGIDs(cntUnknownDofGIDs,-1);
139 size_t cntUnknownOffset = 0;
140 for(
int p = 0; p < comm->getRank(); p++) cntUnknownOffset += numUnknownDofGIDs[p];
141 for(
size_t k=0; k < Teuchos::as<size_t>(ovlUnknownStatusGids.
size()); k++) {
142 lUnknownDofGIDs[k+cntUnknownOffset] = ovlUnknownStatusGids[k];
144 if(cntUnknownDofGIDs > 0)
145 Teuchos::reduceAll(*comm,Teuchos::REDUCE_MAX,Teuchos::as<int>(cntUnknownDofGIDs),&lUnknownDofGIDs[0],&gUnknownDofGIDs[0]);
146 std::sort(gUnknownDofGIDs.begin(), gUnknownDofGIDs.end());
147 gUnknownDofGIDs.erase(std::unique(gUnknownDofGIDs.begin(), gUnknownDofGIDs.end()), gUnknownDofGIDs.end());
150 for(
size_t k=0; k < gUnknownDofGIDs.size(); k++) {
151 GlobalOrdinal curgid = gUnknownDofGIDs[k];
157 std::vector<int> myFoundDofGIDs(comm->getSize(),0);
158 std::vector<int> numFoundDofGIDs(comm->getSize(),0);
159 myFoundDofGIDs[comm->getRank()] = ovlFoundStatusGids.
size();
160 Teuchos::reduceAll(*comm,Teuchos::REDUCE_MAX,comm->getSize(),&myFoundDofGIDs[0],&numFoundDofGIDs[0]);
163 size_t cntFoundDofGIDs = 0;
164 for(
int p = 0; p < comm->getSize(); p++) cntFoundDofGIDs += numFoundDofGIDs[p];
165 std::vector<GlobalOrdinal> lFoundDofGIDs(cntFoundDofGIDs,-1);
166 std::vector<GlobalOrdinal> gFoundDofGIDs(cntFoundDofGIDs,-1);
168 size_t cntFoundOffset = 0;
169 for(
int p = 0; p < comm->getRank(); p++) cntFoundOffset += numFoundDofGIDs[p];
170 for(
size_t k=0; k < Teuchos::as<size_t>(ovlFoundStatusGids.
size()); k++) {
171 lFoundDofGIDs[k+cntFoundOffset] = ovlFoundStatusGids[k];
173 if(cntFoundDofGIDs > 0)
174 Teuchos::reduceAll(*comm,Teuchos::REDUCE_MAX,Teuchos::as<int>(cntFoundDofGIDs),&lFoundDofGIDs[0],&gFoundDofGIDs[0]);
177 for(
size_t i = 0; i<input.
getColMap()->getNodeNumElements(); i++) {
178 GlobalOrdinal gcid = input.
getColMap()->getGlobalElement(i);
180 std::find(gFoundDofGIDs.begin(), gFoundDofGIDs.end(), gcid) != gFoundDofGIDs.end()) {
203 bool bThyraMode =
false) {
205 size_t numRows = rangeMapExtractor->NumMaps();
206 size_t numCols = domainMapExtractor->NumMaps();
208 TEUCHOS_TEST_FOR_EXCEPTION(rangeMapExtractor->getThyraMode() ==
true,
Xpetra::Exceptions::Incompatible,
"Xpetra::MatrixUtils::Split: RangeMapExtractor must not use Thyra style numbering of GIDs. The MapExtractor must contain all GIDs of the full range map in order to define a proper splitting.")
209 TEUCHOS_TEST_FOR_EXCEPTION(domainMapExtractor->getThyraMode() ==
true,
Xpetra::Exceptions::Incompatible,
"Xpetra::MatrixUtils::Split: DomainMapExtractor must not use Thyra style numbering of GIDs. The MapExtractor must contain all GIDs of the full domain map in order to define a proper splitting.")
221 TEUCHOS_TEST_FOR_EXCEPTION(fullDomainMap->getMaxAllGlobalIndex() != input.
getColMap()->getMaxAllGlobalIndex(),
Xpetra::Exceptions::Incompatible,
"Xpetra::MatrixUtils::Split: DomainMapExtractor incompatible to domain map of input matrix. fullDomainMap->getMaxAllGlobalIndex() = " << fullDomainMap->getMaxAllGlobalIndex() <<
" vs. input.getColMap()->getMaxAllGlobalIndex() = " << input.
getColMap()->getMaxAllGlobalIndex())
228 if(columnMapExtractor == Teuchos::null) {
231 std::vector<Teuchos::RCP<const Map> > ovlxmaps(numCols, Teuchos::null);
232 for(
size_t c = 0; c < numCols; c++) {
235 ovlxmaps[c] = colMap;
241 myColumnMapExtractor = columnMapExtractor;
248 if(bThyraMode ==
true) {
250 std::vector<Teuchos::RCP<const Map> > thyRgMapExtractorMaps(numRows, Teuchos::null);
251 for (
size_t r = 0; r < numRows; r++) {
255 if(strRangeMap != Teuchos::null) {
256 std::vector<size_t> strInfo = strRangeMap->getStridingData();
257 GlobalOrdinal offset = strRangeMap->getOffset();
258 LocalOrdinal stridedBlockId = strRangeMap->getStridedBlockId();
260 thyRgMapExtractorMaps[r] = strShrinkedMap;
262 thyRgMapExtractorMaps[r] = shrinkedMap;
268 std::vector<Teuchos::RCP<const Map> > thyDoMapExtractorMaps (numCols, Teuchos::null);
269 std::vector<Teuchos::RCP<const Map> > thyColMapExtractorMaps(numCols, Teuchos::null);
270 for (
size_t c = 0; c < numCols; c++) {
275 if(strDomainMap != Teuchos::null) {
276 std::vector<size_t> strInfo = strDomainMap->getStridingData();
277 GlobalOrdinal offset = strDomainMap->getOffset();
278 LocalOrdinal stridedBlockId = strDomainMap->getStridedBlockId();
280 thyDoMapExtractorMaps[c] = strShrinkedDomainMap;
282 thyDoMapExtractorMaps[c] = shrinkedDomainMap;
287 if(strColMap != Teuchos::null) {
288 std::vector<size_t> strInfo = strColMap->getStridingData();
289 GlobalOrdinal offset = strColMap->getOffset();
290 LocalOrdinal stridedBlockId = strColMap->getStridedBlockId();
292 thyColMapExtractorMaps[c] = strShrinkedColMap;
294 thyColMapExtractorMaps[c] = shrinkedColMap;
306 std::vector<Teuchos::RCP<Matrix> > subMatrices(numRows*numCols, Teuchos::null);
307 for (
size_t r = 0; r < numRows; r++) {
308 for (
size_t c = 0; c < numCols; c++) {
312 if(bThyraMode ==
true)
332 doCheck->putScalar(1.0);
338 doCheck->putScalar(-1.0);
339 coCheck->putScalar(-1.0);
342 for (
size_t rrr = 0; rrr < input.
getDomainMap()->getNodeNumElements(); rrr++) {
344 GlobalOrdinal
id = input.
getDomainMap()->getGlobalElement(rrr);
347 size_t BlockId = domainMapExtractor->getMapIndexForGID(
id);
349 doCheckData[rrr] = Teuchos::as<Scalar>(BlockId);
357 for (
size_t rr = 0; rr < input.
getRowMap()->getNodeNumElements(); rr++) {
360 GlobalOrdinal growid = input.
getRowMap()->getGlobalElement(rr);
369 size_t rowBlockId = rangeMapExtractor->getMapIndexForGID(growid);
372 GlobalOrdinal subblock_growid = growid;
373 if(bThyraMode ==
true) {
375 LocalOrdinal lrowid = rangeMapExtractor->getMap(rowBlockId)->getLocalElement(growid);
377 subblock_growid = thyRangeMapExtractor->getMap(rowBlockId,
true)->getGlobalElement(lrowid);
390 for(
size_t i=0; i<(size_t)indices.
size(); i++) {
392 GlobalOrdinal gcolid = input.
getColMap()->getGlobalElement(indices[i]);
394 size_t colBlockId = myColumnMapExtractor->getMapIndexForGID(gcolid);
398 GlobalOrdinal subblock_gcolid = gcolid;
399 if(bThyraMode ==
true) {
401 LocalOrdinal lcolid = myColumnMapExtractor->getMap(colBlockId)->getLocalElement(gcolid);
403 subblock_gcolid = thyColMapExtractor->getMap(colBlockId,
true)->getGlobalElement(lcolid);
405 blockColIdx [colBlockId].push_back(subblock_gcolid);
406 blockColVals[colBlockId].push_back(vals[i]);
409 for (
size_t c = 0; c < numCols; c++) {
410 subMatrices[rowBlockId*numCols+c]->insertGlobalValues(subblock_growid,blockColIdx[c].view(0,blockColIdx[c].size()),blockColVals[c].view(0,blockColVals[c].size()));
417 if(bThyraMode ==
true) {
418 for (
size_t r = 0; r < numRows; r++) {
419 for (
size_t c = 0; c < numCols; c++) {
420 subMatrices[r*numCols+c]->fillComplete(thyDomainMapExtractor->getMap(c,
true), thyRangeMapExtractor->getMap(r,
true));
425 for (
size_t r = 0; r < numRows; r++) {
426 for (
size_t c = 0; c < numCols; c++) {
427 subMatrices[r*numCols+c]->fillComplete(domainMapExtractor->getMap(c), rangeMapExtractor->getMap(r));
433 for (
size_t r = 0; r < numRows; r++) {
434 for (
size_t c = 0; c < numCols; c++) {
435 bA->setMatrix(r,c,subMatrices[r*numCols+c]);
449 using Teuchos::rcp_dynamic_cast;
451 GlobalOrdinal gZeroDiags;
452 bool usedEfficientPath =
false;
454#ifdef HAVE_MUELU_TPETRA
458 tpCrsAc = rcp_dynamic_cast<TpetraCrsMatrix>(crsWrapAc->getCrsMatrix());
461 auto tpCrsGraph = tpCrsAc->getTpetra_CrsMatrix()->getCrsGraph();
462 size_t numRows = Ac->getRowMap()->getNodeNumElements();
463 typedef typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::offset_device_view_type offset_type;
464 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
465 auto offsets = offset_type(Kokkos::ViewAllocateWithoutInitializing(
"offsets"), numRows);
466 tpCrsGraph->getLocalDiagOffsets(offsets);
469 Tpetra::Details::OrdinalTraits<typename offset_type::value_type>::invalid ();
471 if (repairZeroDiagonals) {
475 LO numMissingDiagonalEntries = 0;
477 Kokkos::parallel_reduce(
"countMissingDiagonalEntries",
478 range_type(0, numRows),
479 KOKKOS_LAMBDA (
const LO i,
LO& missing) {
480 missing += (offsets(i) == STINV);
481 }, numMissingDiagonalEntries);
483 GlobalOrdinal gNumMissingDiagonalEntries;
484 Teuchos::reduceAll(*(Ac->getRowMap()->getComm()), Teuchos::REDUCE_SUM, Teuchos::as<GlobalOrdinal>(numMissingDiagonalEntries),
485 Teuchos::outArg(gNumMissingDiagonalEntries));
487 if (gNumMissingDiagonalEntries == 0) {
490 auto lclA = tpCrsAc->getTpetra_CrsMatrix()->getLocalMatrixDevice();
492 using ATS = Kokkos::ArithTraits<Scalar>;
493 using impl_ATS = Kokkos::ArithTraits<typename ATS::val_type>;
497 Kokkos::parallel_reduce(
"fixSmallDiagonalEntries",
498 range_type(0, numRows),
499 KOKKOS_LAMBDA (
const LO i,
LO& fixed) {
500 const auto offset = offsets(i);
501 auto curRow = lclA.row (i);
502 if (impl_ATS::magnitude(curRow.value(offset)) <= threshold) {
503 curRow.value(offset) = replacementValue;
508 Teuchos::reduceAll(*(Ac->getRowMap()->getComm()), Teuchos::REDUCE_SUM, Teuchos::as<GlobalOrdinal>(lZeroDiags),
509 Teuchos::outArg(gZeroDiags));
511 usedEfficientPath =
true;
517 auto lclA = tpCrsAc->getTpetra_CrsMatrix()->getLocalMatrixDevice();
519 using ATS = Kokkos::ArithTraits<Scalar>;
520 using impl_ATS = Kokkos::ArithTraits<typename ATS::val_type>;
524 Kokkos::parallel_reduce(
"detectSmallDiagonalEntries",
525 range_type(0, numRows),
526 KOKKOS_LAMBDA (
const LO i,
LO& small) {
527 const auto offset = offsets(i);
531 auto curRow = lclA.row (i);
532 if (impl_ATS::magnitude(curRow.value(offset)) <= threshold) {
538 Teuchos::reduceAll(*(Ac->getRowMap()->getComm()), Teuchos::REDUCE_SUM, Teuchos::as<GlobalOrdinal>(lZeroDiags),
539 Teuchos::outArg(gZeroDiags));
541 usedEfficientPath =
true;
547 if (!usedEfficientPath) {
550 Ac->getLocalDiagCopy(*diagVec);
552 LocalOrdinal lZeroDiags = 0;
555 for (
size_t i = 0; i < rowMap->getNodeNumElements(); i++) {
556 if (TST::magnitude(diagVal[i]) <= threshold) {
561 Teuchos::reduceAll(*(rowMap->getComm()), Teuchos::REDUCE_SUM, Teuchos::as<GlobalOrdinal>(lZeroDiags),
562 Teuchos::outArg(gZeroDiags));
564 if (repairZeroDiagonals && gZeroDiags > 0) {
582 for (
size_t r = 0; r < rowMap->getNodeNumElements(); r++) {
583 if (TST::magnitude(diagVal[r]) <= threshold) {
584 GlobalOrdinal grid = rowMap->getGlobalElement(r);
586 valout[0] = replacementValue;
587 fixDiagMatrix->insertGlobalValues(grid,indout(), valout());
592 fixDiagMatrix->fillComplete(Ac->getDomainMap(),Ac->getRangeMap());
596 Xpetra::MatrixMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TwoMatrixAdd(*Ac,
false, 1.0, *fixDiagMatrix,
false, 1.0, newAc, fos);
597 if (Ac->IsView(
"stridedMaps"))
598 newAc->CreateView(
"stridedMaps", Ac);
601 fixDiagMatrix = Teuchos::null;
607 p->set(
"DoOptimizeStorage",
true);
616 fos <<
"CheckRepairMainDiagonal: " << (repairZeroDiagonals ?
"repaired " :
"found ")
617 << gZeroDiags <<
" too small entries (threshold = " << threshold <<
") on main diagonal of Ac." << std::endl;
619#ifdef HAVE_XPETRA_DEBUG
625 Ac->getLocalDiagCopy(*diagVec);
626 diagVal = diagVec->getData(0);
627 for (
size_t r = 0; r < Ac->getRowMap()->getNodeNumElements(); r++) {
628 if (TST::magnitude(diagVal[r]) <= threshold) {
629 fos <<
"Error: there are too small entries left on diagonal after repair..." << std::endl;
652 LocalOrdinal numPDEs = A->GetFixedBlockSize();
655 Scalar zero = TST::zero();
656 Scalar one = TST::one();
660 A->getLocalDiagCopy(*diag);
662 size_t N = A->getRowMap()->getNodeNumElements();
665 std::vector<MT> l_diagMax(numPDEs), g_diagMax(numPDEs);
666 for(
size_t i=0; i<N; i++) {
667 int pde = (int) (i % numPDEs);
669 l_diagMax[pde] = TST::magnitude(dataVal[i]);
671 l_diagMax[pde] = std::max(l_diagMax[pde],TST::magnitude(dataVal[i]));
673 Teuchos::reduceAll(*A->getRowMap()->getComm(), Teuchos::REDUCE_MAX, numPDEs, l_diagMax.data(), g_diagMax.data() );
679 for (
size_t i = 0; i<N; i++) {
680 GlobalOrdinal GRID = A->getRowMap()->getGlobalElement(i);
681 int pde = (int) (i % numPDEs);
683 if (TST::magnitude(dataVal[i]) < relativeThreshold[pde] * g_diagMax[pde])
684 value[0] = relativeThreshold[pde] * g_diagMax[pde] - TST::magnitude(dataVal[i]);
687 boostMatrix->insertGlobalValues(GRID,index(),value());
689 boostMatrix->fillComplete(A->getDomainMap(),A->getRangeMap());
693 Xpetra::MatrixMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TwoMatrixAdd(*A,
false,one, *boostMatrix,
false,one,newA,fos);
694 if (A->IsView(
"stridedMaps"))
695 newA->CreateView(
"stridedMaps", A);
707#if defined(HAVE_XPETRA_EPETRA)
711#ifdef HAVE_XPETRA_TPETRA
714 Tpetra::Details::extractBlockDiagonal(At,Dt);
727#if defined(HAVE_XPETRA_EPETRA)
731#ifdef HAVE_XPETRA_TPETRA
734 Tpetra::Details::inverseScaleBlockDiagonal(Dt,doTranspose,St);
743 LO numRows = Teuchos::as<LocalOrdinal>(rowMap->getNodeNumElements());
745 for (
LO rowLID = 0; rowLID < numRows; rowLID++) {
746 GO rowGID = rowMap->getGlobalElement(rowLID);
747 LO colLID = colMap->getLocalElement(rowGID);
748 if (rowLID != colLID) {
750 std::cerr <<
"On rank " << comm->getRank() <<
", GID " << rowGID <<
" is LID " << rowLID <<
"in the rowmap, but LID " << colLID <<
" in the column map.\n";
754 "Local parts of row and column map do not match!");
766 std::vector<size_t>& rangeStridingInfo, std::vector<size_t>& domainStridingInfo)
771 if (matrix->IsView(
"stridedMaps") ==
true) matrix->RemoveView(
"stridedMaps");
772 matrix->CreateView(
"stridedMaps", stridedRowMap, stridedColMap);
779#define XPETRA_MATRIXUTILS_SHORT
void push_back(const value_type &x)
static RCP< Time > getNewTimer(const std::string &name)
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const =0
The Map describing the parallel distribution of this object.
Exception throws to report incompatible objects (like maps).
Exception throws to report errors in the internal logical of the program.
static RCP< const Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2TpetraCrs(RCP< Matrix > Op)
static RCP< Import< LocalOrdinal, GlobalOrdinal, Node > > Build(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &source, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &target, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null)
Constructor specifying the number of non-zeros for all rows.
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< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > concatenateMaps(const std::vector< Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > > &subMaps)
Helper function to concatenate several maps.
static Teuchos::RCP< Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > shrinkMapGIDs(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &input, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &nonOvlInput)
Helper function to shrink the GIDs and generate a standard map whith GIDs starting at 0.
virtual UnderlyingLib lib() const =0
Get the library used by this object (Tpetra or Epetra?)
virtual bool isNodeGlobalElement(GlobalOrdinal globalIndex) const =0
Whether the given global index is valid for this Map on this process.
static RCP< Matrix > Build(const RCP< const Map > &rowMap)
static void TwoMatrixAdd(const Matrix &A, bool transposeA, SC alpha, Matrix &B, SC beta)
Helper function to calculate B = alpha*A + beta*B.
Xpetra utility class for common matrix-related routines.
static void checkLocalRowMapMatchesColMap(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
static void convertMatrixToStridedMaps(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > matrix, std::vector< size_t > &rangeStridingInfo, std::vector< size_t > &domainStridingInfo)
static Teuchos::RCP< Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > findColumnSubMap(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &input, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &domainMap)
static void CheckRepairMainDiagonal(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Ac, bool const &repairZeroDiagonals, Teuchos::FancyOStream &fos, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType threshold=Teuchos::ScalarTraits< typename Teuchos::ScalarTraits< Scalar >::magnitudeType >::zero(), const Scalar replacementValue=Teuchos::ScalarTraits< Scalar >::one())
static void RelativeDiagonalBoost(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Teuchos::ArrayView< const double > &relativeThreshold, Teuchos::FancyOStream &fos)
static void extractBlockDiagonal(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diagonal)
static Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > xpetraGidNumbering2ThyraGidNumbering(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &input)
static Teuchos::RCP< Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > SplitMatrix(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &input, Teuchos::RCP< const Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > rangeMapExtractor, Teuchos::RCP< const Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > domainMapExtractor, Teuchos::RCP< const Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > columnMapExtractor=Teuchos::null, bool bThyraMode=false)
static void inverseScaleBlockDiagonal(Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &blockDiagonal, bool doTranspose, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &toBeScaled)
Xpetra-specific matrix class.
virtual const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
virtual size_t getNodeMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on this node.
virtual void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const =0
Extract a const, non-persisting view of local indices in a specified row of the matrix.
virtual const RCP< const Map > & getRowMap() const
Returns the Map that describes the row distribution in this matrix.
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 Teuchos::ArrayRCP< const Scalar > getData(size_t j) const =0
Const view of the local values in a particular vector of this multivector.
virtual size_t getLocalLength() const =0
Local number of rows on the calling process.
virtual Teuchos::RCP< const Map > getRangeMap() const =0
The Map associated with the range of this operator, which must be compatible with Y....
virtual Teuchos::RCP< const Map > getDomainMap() const =0
The Map associated with the domain of this operator, which must be compatible with X....
static RCP< Xpetra::StridedMap< LocalOrdinal, GlobalOrdinal, Node > > Build(UnderlyingLib lib, global_size_t numGlobalElements, GlobalOrdinal indexBase, std::vector< size_t > &stridingInfo, const Teuchos::RCP< const Teuchos::Comm< int > > &comm, LocalOrdinal stridedBlockId=-1, GlobalOrdinal offset=0, LocalGlobal lg=Xpetra::GloballyDistributed)
Map constructor with Xpetra-defined contiguous uniform distribution.
Class that stores a strided map.
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_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)