46#ifndef MUELU_UTILITIESBASE_DECL_HPP
47#define MUELU_UTILITIESBASE_DECL_HPP
57#include <Teuchos_DefaultComm.hpp>
58#include <Teuchos_ScalarTraits.hpp>
59#include <Teuchos_ParameterList.hpp>
61#include <Xpetra_BlockedCrsMatrix_fwd.hpp>
62#include <Xpetra_CrsMatrix_fwd.hpp>
63#include <Xpetra_CrsMatrixWrap_fwd.hpp>
64#include <Xpetra_Map_fwd.hpp>
65#include <Xpetra_BlockedMap_fwd.hpp>
66#include <Xpetra_MapFactory_fwd.hpp>
67#include <Xpetra_Matrix_fwd.hpp>
68#include <Xpetra_MatrixFactory_fwd.hpp>
69#include <Xpetra_MultiVector_fwd.hpp>
70#include <Xpetra_MultiVectorFactory_fwd.hpp>
71#include <Xpetra_Operator_fwd.hpp>
72#include <Xpetra_Vector_fwd.hpp>
73#include <Xpetra_BlockedMultiVector.hpp>
74#include <Xpetra_BlockedVector.hpp>
75#include <Xpetra_VectorFactory_fwd.hpp>
76#include <Xpetra_ExportFactory.hpp>
78#include <Xpetra_Import.hpp>
79#include <Xpetra_ImportFactory.hpp>
80#include <Xpetra_MatrixMatrix.hpp>
81#include <Xpetra_CrsMatrixWrap.hpp>
82#include <Xpetra_StridedMap.hpp>
90#define MueLu_sumAll(rcpComm, in, out) \
91 Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_SUM, in, Teuchos::outArg(out))
92#define MueLu_minAll(rcpComm, in, out) \
93 Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_MIN, in, Teuchos::outArg(out))
94#define MueLu_maxAll(rcpComm, in, out) \
95 Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_MAX, in, Teuchos::outArg(out))
110#undef MUELU_UTILITIESBASE_SHORT
113 typedef Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>
CrsMatrixWrap;
114 typedef Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
CrsMatrix;
115 typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
Matrix;
116 typedef Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
Vector;
117 typedef Xpetra::BlockedVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
BlockedVector;
118 typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
MultiVector;
120 typedef Xpetra::BlockedMap<LocalOrdinal,GlobalOrdinal,Node>
BlockedMap;
121 typedef Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node>
Map;
123 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType
Magnitude;
128 return Teuchos::null;
139 size_t numRows = A.getRowMap()->getNodeNumElements();
140 Teuchos::ArrayRCP<Scalar> diag(numRows);
141 Teuchos::ArrayView<const LocalOrdinal>
cols;
142 Teuchos::ArrayView<const Scalar>
vals;
143 for (
size_t i = 0; i < numRows; ++i) {
146 for (; j <
cols.size(); ++j) {
147 if (Teuchos::as<size_t>(
cols[j]) == i) {
152 if (j ==
cols.size()) {
154 diag[i] = Teuchos::ScalarTraits<Scalar>::zero();
167 Teuchos::TimeMonitor
MM = *Teuchos::TimeMonitor::getNewTimer(
"UtilitiesBase::GetMatrixDiagonalInverse");
170 RCP<Vector> diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,
true);
172 A.getLocalDiagCopy(*diag);
186 Magnitude tol = Teuchos::ScalarTraits<Scalar>::eps()*100,
190 typedef Teuchos::ScalarTraits<Scalar>
TST;
199 Teuchos::RCP<const Matrix>
rcpA = Teuchos::rcpFromRef(A);
202 Teuchos::rcp_dynamic_cast<const Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(
rcpA);
203 if(
bA == Teuchos::null) {
205 diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,
true);
207 Teuchos::Array<Scalar>
regSum(diag->getLocalLength());
208 Teuchos::ArrayView<const LocalOrdinal>
cols;
209 Teuchos::ArrayView<const Scalar>
vals;
211 std::vector<int>
nnzPerRow(rowMap->getNodeNumElements());
214 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
227 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
243 "UtilitiesBase::GetLumpedMatrixDiagonal(): extracting reciprocal of diagonal of a blocked matrix is not supported");
244 diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(
bA->getRangeMapExtractor()->getFullMap(),
true);
247 for (
size_t col = 0; col <
bA->Cols(); ++col) {
248 if (!
bA->getMatrix(
row,col).is_null()) {
250 bool bThyraMode =
bA->getRangeMapExtractor()->getThyraMode() && (Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(
bA->getMatrix(
row,col)) == Teuchos::null);
253 ddtemp->update(Teuchos::as<Scalar>(1.0),*
dd,Teuchos::as<Scalar>(1.0));
271 size_t numRows = A.getRowMap()->getNodeNumElements();
273 Teuchos::ArrayRCP<Magnitude>
maxvec(numRows);
274 Teuchos::ArrayView<const LocalOrdinal>
cols;
275 Teuchos::ArrayView<const Scalar>
vals;
276 for (
size_t i = 0; i < numRows; ++i) {
280 if (Teuchos::as<size_t>(
cols[j]) != i) {
281 mymax = std::max(
mymax,-Teuchos::ScalarTraits<Scalar>::real(
vals[j]));
289 static Teuchos::ArrayRCP<Magnitude>
GetMatrixMaxMinusOffDiagonal(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
const Xpetra::Vector<LocalOrdinal,LocalOrdinal,GlobalOrdinal,Node> &
BlockNumber) {
294 size_t numRows = A.getRowMap()->getNodeNumElements();
296 Teuchos::ArrayRCP<Magnitude>
maxvec(numRows);
297 Teuchos::ArrayView<const LocalOrdinal>
cols;
298 Teuchos::ArrayView<const Scalar>
vals;
299 for (
size_t i = 0; i < numRows; ++i) {
304 mymax = std::max(
mymax,-Teuchos::ScalarTraits<Scalar>::real(
vals[j]));
323 RCP<Vector> ret = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(v->getMap(),
true);
327 if(
bv.is_null() ==
false) {
331 for(
size_t r = 0;
r <
bmap->getNumMaps(); ++
r) {
344 if(Teuchos::ScalarTraits<Scalar>::magnitude(
inputVals[i]) >
tol)
366 RCP<Vector> localDiag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap);
372 Teuchos::ArrayRCP<size_t> offsets;
373 crsOp->getLocalDiagOffsets(offsets);
384 RCP<Vector> diagonal = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(colMap);
386 importer = A.getCrsGraph()->getImporter();
388 importer = Xpetra::ImportFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap, colMap);
406 using STS =
typename Teuchos::ScalarTraits<SC>;
418 size_t nnz = A.getNumEntriesInLocalRow(
row);
421 A.getLocalRowView(
row, indices,
vals);
434 importer = A.getCrsGraph()->getImporter();
436 importer = Xpetra::ImportFactory<LO,GO,Node>::Build(rowMap, colMap);
447 using STS =
typename Teuchos::ScalarTraits<Scalar>;
448 using MTS =
typename Teuchos::ScalarTraits<Magnitude>;
464 size_t nnz = A.getNumEntriesInLocalRow(
rowIdx);
480 importer = A.getCrsGraph()->getImporter();
482 importer = Xpetra::ImportFactory<LO,GO,Node>::Build(rowMap, colMap);
496 const size_t numVecs = X.getNumVectors();
505 const size_t numVecs = X.getNumVectors();
514 const size_t numVecs = X.getNumVectors();
516 RCP<MultiVector> RES = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(RHS.getMap(),
numVecs,
false);
517 Op.residual(X,RHS,*
RES);
525 Op.residual(X,RHS,
Resid);
543 "Utils::PowerMethod: operator must have domain and range maps that are equivalent.");
548 RCP<Vector> diagVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
550 diagInvVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
572 "Utils::PowerMethod: operator must have domain and range maps that are equivalent.");
575 RCP<Vector> q = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getDomainMap());
576 RCP<Vector> r = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRangeMap());
577 RCP<Vector> z = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRangeMap());
582 Teuchos::Array<Magnitude>
norms(1);
584 typedef Teuchos::ScalarTraits<Scalar> STS;
605 std::cout <<
"Iter = " <<
iter
607 <<
" Residual of A*q - lambda*q = " << residual
629 Scalar d = Teuchos::ScalarTraits<Scalar>::zero();
631 d += (v[j][
i0] - v[j][
i1])*(v[j][
i0] - v[j][
i1]);
633 return Teuchos::ScalarTraits<Scalar>::magnitude(
d);
643 Scalar d = Teuchos::ScalarTraits<Scalar>::zero();
647 return Teuchos::ScalarTraits<Scalar>::magnitude(
d);
665 typedef Teuchos::ScalarTraits<Scalar> STS;
671 A.getLocalRowView(
row, indices,
vals);
672 size_t nnz = A.getNumEntriesInLocalRow(
row);
675 for (col = 0; col <
nnz; col++)
676 if ( (indices[col] !=
row) && STS::magnitude(
vals[col]) >
tol) {
689 A.getLocalRowView(
row, indices,
vals);
690 size_t nnz = A.getNumEntriesInLocalRow(
row);
692 for (
size_t col = 0; col <
nnz; col++)
693 if ( (indices[col] !=
row) && STS::magnitude(
vals[col]) >
tol) {
720 Teuchos::RCP<Vector>
diagVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
725 typedef Teuchos::ScalarTraits<Scalar> STS;
730 A.getLocalRowView(
row, indices,
vals);
733 for (
decltype(indices.size()) col = 0; col < indices.size(); col++) {
734 if ( indices[col] !=
row) {
756 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitudeType;
757 const magnitudeType
eps = 2.0*Teuchos::ScalarTraits<magnitudeType>::eps();
759 nonzeros[i] = (Teuchos::ScalarTraits<Scalar>::magnitude(
vals[i]) >
eps);
775 const Scalar one = Teuchos::ScalarTraits<Scalar>::one();
788 A.getLocalRowView(i,indices,values);
797 globalColsToZero = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(
domMap, 1,
true);
824 typedef Teuchos::ScalarTraits<Scalar> STS;
825 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
826 typedef Teuchos::ScalarTraits<MT>
MTS;
829 size_t nnz = A.getNumEntriesInLocalRow(
row);
832 A.getLocalRowView(
row, indices,
vals);
852 typedef Teuchos::ScalarTraits<Scalar> STS;
853 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
854 typedef Teuchos::ScalarTraits<MT>
MTS;
861 size_t nnz = A.getNumEntriesInLocalRow(
row);
864 A.getLocalRowView(
row, indices,
vals);
896 static Teuchos::ArrayRCP<const bool>
DetectDirichletCols(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
898 Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
899 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
900 Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
domMap = A.getDomainMap();
901 Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > colMap = A.getColMap();
902 Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
myColsToZero = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(colMap,1);
907 Teuchos::ArrayView<const LocalOrdinal> indices;
908 Teuchos::ArrayView<const Scalar> values;
909 A.getLocalRowView(i,indices,values);
915 Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
globalColsToZero = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(
domMap,1);
917 Teuchos::RCP<Xpetra::Export<LocalOrdinal,GlobalOrdinal,Node> >
exporter = Xpetra::ExportFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(colMap,
domMap);
923 Teuchos::ArrayRCP<bool>
dirichletCols(colMap->getNodeNumElements(),
true);
924 Magnitude eps = Teuchos::ScalarTraits<Magnitude>::eps();
925 for(
size_t i=0; i<colMap->getNodeNumElements(); i++) {
936 static Scalar Frobenius(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& B) {
947 Teuchos::ArrayView<const LocalOrdinal>
indA,
indB;
948 Teuchos::ArrayView<const Scalar>
valA,
valB;
965 size_t numRows = A.getNodeNumRows();
966 for (
size_t i = 0; i < numRows; i++) {
973 for (
size_t j = 0; j <
nnzB; j++)
976 for (
size_t j = 0; j <
nnzA; j++) {
985 for (
size_t j = 0; j <
nnzB; j++)
1012 std::ostringstream
errStr;
1013 errStr <<
"Error detected with random seed = " <<
mySeed <<
". It should be in the interval [1,2^31-2].";
1018 Teuchos::ScalarTraits<Scalar>::seedrandom(
mySeed);
1030 static void FindDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &A,
1032 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
1034 for(
size_t i=0; i<A->getNodeNumRows(); i++) {
1035 Teuchos::ArrayView<const LocalOrdinal> indices;
1036 Teuchos::ArrayView<const Scalar> values;
1037 A->getLocalRowView(i,indices,values);
1039 for (
size_t j=0; j<(
size_t)indices.size(); j++) {
1040 if (Teuchos::ScalarTraits<Scalar>::magnitude(values[j]) > Teuchos::ScalarTraits<MT>::eps()) {
1056 Scalar one =Teuchos::ScalarTraits<Scalar>::one();
1057 Scalar zero =Teuchos::ScalarTraits<Scalar>::zero();
1062 Teuchos::ArrayView<const LocalOrdinal> indices;
1063 Teuchos::ArrayView<const Scalar> values;
1067 for(
size_t j=0; j<(
size_t)indices.size(); j++) {
1086 const Scalar one = Teuchos::ScalarTraits<Scalar>::one();
1087 const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
1093 Teuchos::ArrayView<const LocalOrdinal> indices;
1094 Teuchos::ArrayView<const Scalar> values;
1095 A->getLocalRowView(i,indices,values);
1097 Teuchos::ArrayRCP<Scalar>
valuesNC(values.size());
1098 for(
size_t j=0; j<(
size_t)indices.size(); j++) {
1104 A->replaceLocalValues(i,indices,
valuesNC());
1112 static void ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
1116 Teuchos::ArrayView<const LocalOrdinal> indices;
1117 Teuchos::ArrayView<const Scalar> values;
1121 for(
size_t j=0; j<(
size_t)indices.size(); j++)
1128 static void ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
1134 Teuchos::ArrayView<const LocalOrdinal> indices;
1135 Teuchos::ArrayView<const Scalar> values;
1136 A->getLocalRowView(i,indices,values);
1139 for(
size_t j=0; j<(
size_t)indices.size(); j++)
1147 static void ZeroDirichletRows(Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& X,
1153 for(
size_t j=0; j<X->getNumVectors(); j++)
1165 for(
size_t i=0; i<A->getNodeNumRows(); i++) {
1166 Teuchos::ArrayView<const LocalOrdinal> indices;
1167 Teuchos::ArrayView<const Scalar> values;
1168 A->getLocalRowView(i,indices,values);
1179 Teuchos::RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> >&
isDirichletRow,
1180 Teuchos::RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> >&
isDirichletCol) {
1183 if(!A->getRowMap()->isSameAs(*A->getDomainMap())) {
1184 throw std::runtime_error(
"UtilitiesBase::FindDirichletRowsAndPropagateToCols row and domain maps must match.");
1194 printf(
"[%d] DirichletRow Ids = ",A->getRowMap()->getComm()->getRank());
1201 isDirichletRow = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getRowMap(),
true);
1202 isDirichletCol = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getColMap(),
true);
1205 Teuchos::ArrayView<int>
dr =
dr_rcp();
1207 Teuchos::ArrayView<int>
dc =
dc_rcp();
1222 const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> & Importer) {
1223 typedef Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node>
IntVector;
1233 if(!Importer.getSourceMap()->isCompatible(*
fullMap))
1234 throw std::runtime_error(
"GenerateBlockedTargetMap(): Map compatibility error");
1242 for(
size_t j=0; j<
map->getNodeNumElements(); j++) {
1253 Teuchos::ArrayView<const int> data =
dataRCP();
1258 for(
size_t i=0; i<
targetMap->getNodeNumElements(); i++) {
1275 static bool MapsAreNested(
const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node>& rowMap,
const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node>& colMap) {
1304#define MUELU_UTILITIESBASE_SHORT
#define MueLu_sumAll(rcpComm, in, out)
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Exception throws to report incompatible objects (like maps).
Exception throws to report errors in the internal logical of the program.
static RCP< Vector > GetMatrixOverlappedDeletedRowsum(const Matrix &A)
Extract Overlapped Matrix Deleted Rowsum.
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &X, const Teuchos::ArrayRCP< const bool > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
static void ZeroDirichletCols(Teuchos::RCP< Matrix > &A, const Teuchos::ArrayRCP< const bool > &dirichletCols, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
static void DetectDirichletColsAndDomains(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayRCP< bool > &dirichletRows, Teuchos::ArrayRCP< bool > dirichletCols, Teuchos::ArrayRCP< bool > dirichletDomain)
Detects Dirichlet columns & domains from a list of Dirichlet rows.
static Teuchos::ScalarTraits< Scalar >::magnitudeType Distance2(const Teuchos::Array< Teuchos::ArrayRCP< const Scalar > > &v, LocalOrdinal i0, LocalOrdinal i1)
Squared distance between two rows in a multivector.
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
static Scalar PowerMethod(const Matrix &A, const RCP< Vector > &diagInvVec, LocalOrdinal niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
Power method.
static Teuchos::ArrayRCP< Scalar > GetMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal.
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
static void ApplyRowSumCriterion(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &BlockNumber, const Magnitude rowSumTol, Teuchos::ArrayRCP< bool > &dirichletRows)
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::zero(), bool count_twos_as_dirichlet=false)
Detect Dirichlet rows.
static Teuchos::ArrayRCP< Magnitude > GetMatrixMaxMinusOffDiagonal(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Return vector containing: max_{i\not=k}(-a_ik), for each for i in the matrix.
static RCP< Vector > GetMatrixDiagonalInverse(const Matrix &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero())
Extract Matrix Diagonal.
static Scalar Frobenius(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B)
Frobenius inner product of two matrices.
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
static Teuchos::RCP< Vector > GetLumpedMatrixDiagonal(Matrix const &A, const bool doReciprocal=false, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool replaceSingleEntryRowWithZero=false)
Extract Matrix Diagonal of lumped matrix.
static Scalar PowerMethod(const Matrix &A, bool scaleByDiag=true, LocalOrdinal niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
Power method.
Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > Vector
static void ApplyOAZToMatrixRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const std::vector< LocalOrdinal > &dirichletRows)
static RCP< Matrix > Crs2Op(RCP< CrsMatrix > Op)
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixWrap
static Teuchos::ArrayRCP< const bool > DetectDirichletRowsExt(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool &bHasZeroDiagonal, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::zero())
Detect Dirichlet rows (extended version)
static void FindDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, std::vector< LocalOrdinal > &dirichletRows, bool count_twos_as_dirichlet=false)
static void ApplyRowSumCriterion(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude rowSumTol, Teuchos::ArrayRCP< bool > &dirichletRows)
Apply Rowsum Criterion.
Xpetra::BlockedVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > BlockedVector
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
static RCP< const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > > GeneratedBlockedTargetMap(const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > &sourceBlockedMap, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &Importer)
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Teuchos::ArrayRCP< const bool > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
static void FindNonZeros(const Teuchos::ArrayRCP< const Scalar > vals, Teuchos::ArrayRCP< bool > nonzeros)
Find non-zero values in an ArrayRCP Compares the value to 2 * machine epsilon.
static bool MapsAreNested(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &rowMap, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &colMap)
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
static void Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS, MultiVector &Resid)
static Teuchos::ScalarTraits< Scalar >::magnitudeType Distance2(const Teuchos::ArrayView< double > &weight, const Teuchos::Array< Teuchos::ArrayRCP< const Scalar > > &v, LocalOrdinal i0, LocalOrdinal i1)
Weighted squared distance between two rows in a multivector.
static Teuchos::ArrayRCP< const bool > DetectDirichletCols(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayRCP< const bool > &dirichletRows)
Detect Dirichlet columns based on Dirichlet rows.
static RCP< Xpetra::Vector< Magnitude, LocalOrdinal, GlobalOrdinal, Node > > GetMatrixOverlappedAbsDeletedRowsum(const Matrix &A)
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS, MultiVector &Resid)
static Teuchos::ArrayRCP< Magnitude > GetMatrixMaxMinusOffDiagonal(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &BlockNumber)
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MultiVector
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > BlockedMultiVector
static void ApplyOAZToMatrixRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Teuchos::ArrayRCP< const bool > &dirichletRows)
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
Extract Overlapped Matrix Diagonal.
static Teuchos::RCP< Vector > GetInverse(Teuchos::RCP< const Vector > v, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero())
Return vector containing inverse of input vector.
Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > BlockedMap
static void FindDirichletRowsAndPropagateToCols(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, Teuchos::RCP< Xpetra::Vector< int, LocalOrdinal, GlobalOrdinal, Node > > &isDirichletRow, Teuchos::RCP< Xpetra::Vector< int, LocalOrdinal, GlobalOrdinal, Node > > &isDirichletCol)
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const std::vector< LocalOrdinal > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
Namespace for MueLu classes and methods.
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode