46 #ifndef MUELU_UTILITIESBASE_DECL_HPP 47 #define MUELU_UTILITIESBASE_DECL_HPP 54 #include <Teuchos_DefaultComm.hpp> 55 #include <Teuchos_ScalarTraits.hpp> 56 #include <Teuchos_ParameterList.hpp> 58 #include <Xpetra_BlockedCrsMatrix_fwd.hpp> 59 #include <Xpetra_CrsMatrix_fwd.hpp> 60 #include <Xpetra_CrsMatrixWrap_fwd.hpp> 61 #include <Xpetra_Map_fwd.hpp> 62 #include <Xpetra_BlockedMap_fwd.hpp> 63 #include <Xpetra_MapFactory_fwd.hpp> 64 #include <Xpetra_Matrix_fwd.hpp> 65 #include <Xpetra_MatrixFactory_fwd.hpp> 66 #include <Xpetra_MultiVector_fwd.hpp> 67 #include <Xpetra_MultiVectorFactory_fwd.hpp> 68 #include <Xpetra_Operator_fwd.hpp> 69 #include <Xpetra_Vector_fwd.hpp> 70 #include <Xpetra_BlockedMultiVector.hpp> 71 #include <Xpetra_BlockedVector.hpp> 72 #include <Xpetra_VectorFactory_fwd.hpp> 73 #include <Xpetra_ExportFactory.hpp> 75 #include <Xpetra_Import.hpp> 76 #include <Xpetra_ImportFactory.hpp> 77 #include <Xpetra_MatrixMatrix.hpp> 78 #include <Xpetra_CrsMatrixWrap.hpp> 88 #define MueLu_sumAll(rcpComm, in, out) \ 89 Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_SUM, in, Teuchos::outArg(out)) 90 #define MueLu_minAll(rcpComm, in, out) \ 91 Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_MIN, in, Teuchos::outArg(out)) 92 #define MueLu_maxAll(rcpComm, in, out) \ 93 Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_MAX, in, Teuchos::outArg(out)) 102 template <
class Scalar,
103 class LocalOrdinal = int,
104 class GlobalOrdinal = LocalOrdinal,
105 class Node = KokkosClassic::DefaultNode::DefaultNodeType>
106 class UtilitiesBase {
108 #undef MUELU_UTILITIESBASE_SHORT 111 typedef Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>
CrsMatrixWrap;
112 typedef Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
CrsMatrix;
113 typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
Matrix;
114 typedef Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
Vector;
115 typedef Xpetra::BlockedVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
BlockedVector;
116 typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
MultiVector;
118 typedef Xpetra::BlockedMap<LocalOrdinal,GlobalOrdinal,Node>
BlockedMap;
119 typedef Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node>
Map;
121 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType
Magnitude;
124 static RCP<Matrix>
Crs2Op(RCP<CrsMatrix> Op) {
126 return Teuchos::null;
137 size_t numRows = A.getRowMap()->getNodeNumElements();
138 Teuchos::ArrayRCP<Scalar> diag(numRows);
139 Teuchos::ArrayView<const LocalOrdinal> cols;
140 Teuchos::ArrayView<const Scalar> vals;
141 for (
size_t i = 0; i < numRows; ++i) {
142 A.getLocalRowView(i, cols, vals);
144 for (; j < cols.size(); ++j) {
145 if (Teuchos::as<size_t>(cols[j]) == i) {
150 if (j == cols.size()) {
152 diag[i] = Teuchos::ScalarTraits<Scalar>::zero();
166 RCP<const Matrix> rcpA = Teuchos::rcpFromRef(A);
168 RCP<const Map> rowMap = rcpA->getRowMap();
169 RCP<Vector> diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,
true);
171 rcpA->getLocalDiagCopy(*diag);
188 size_t numRows = A.getRowMap()->getNodeNumElements();
189 Teuchos::ArrayRCP<Scalar> diag(numRows);
190 Teuchos::ArrayView<const LocalOrdinal> cols;
191 Teuchos::ArrayView<const Scalar> vals;
192 for (
size_t i = 0; i < numRows; ++i) {
193 A.getLocalRowView(i, cols, vals);
194 diag[i] = Teuchos::ScalarTraits<Scalar>::zero();
195 for (LocalOrdinal j = 0; j < cols.size(); ++j) {
196 diag[i] += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
210 RCP<Vector> diag = Teuchos::null;
212 RCP<const Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > bA =
213 Teuchos::rcp_dynamic_cast<
const Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(rcpA);
214 if(bA == Teuchos::null) {
215 RCP<const Map> rowMap = rcpA->getRowMap();
216 diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,
true);
217 ArrayRCP<Scalar> diagVals = diag->getDataNonConst(0);
218 Teuchos::ArrayView<const LocalOrdinal> cols;
219 Teuchos::ArrayView<const Scalar> vals;
220 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
221 rcpA->getLocalRowView(i, cols, vals);
222 diagVals[i] = Teuchos::ScalarTraits<Scalar>::zero();
223 for (LocalOrdinal j = 0; j < cols.size(); ++j) {
224 diagVals[i] += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
232 diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(bA->getRangeMapExtractor()->getFullMap(),
true);
234 for (
size_t row = 0; row < bA->Rows(); ++row) {
235 for (
size_t col = 0; col < bA->Cols(); ++col) {
236 if (!bA->getMatrix(row,col).is_null()) {
238 bool bThyraMode = bA->getRangeMapExtractor()->getThyraMode() && (Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(bA->getMatrix(row,col)) == Teuchos::null);
239 RCP<Vector> ddtemp = bA->getRangeMapExtractor()->ExtractVector(diag,row,bThyraMode);
241 ddtemp->update(Teuchos::as<Scalar>(1.0),*dd,Teuchos::as<Scalar>(1.0));
242 bA->getRangeMapExtractor()->InsertVector(ddtemp,row,diag,bThyraMode);
260 static Teuchos::RCP<Vector>
GetInverse(Teuchos::RCP<const Vector> v, Magnitude tol = Teuchos::ScalarTraits<Scalar>::eps()*100, Scalar tolReplacement = Teuchos::ScalarTraits<Scalar>::zero()) {
262 RCP<Vector> ret = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(v->getMap(),
true);
265 RCP<const BlockedVector> bv = Teuchos::rcp_dynamic_cast<
const BlockedVector>(v);
266 if(bv.is_null() ==
false) {
267 RCP<BlockedVector> bret = Teuchos::rcp_dynamic_cast<BlockedVector>(ret);
268 TEUCHOS_TEST_FOR_EXCEPTION(bret.is_null() ==
true,
MueLu::Exceptions::RuntimeError,
"MueLu::UtilitiesBase::GetInverse: return vector should be of type BlockedVector");
269 RCP<const BlockedMap> bmap = bv->getBlockedMap();
270 for(
size_t r = 0; r < bmap->getNumMaps(); ++r) {
271 RCP<const MultiVector> submvec = bv->getMultiVector(r,bmap->getThyraMode());
272 RCP<const Vector> subvec = submvec->getVector(0);
274 bret->setMultiVector(r, subvecinf, bmap->getThyraMode());
280 ArrayRCP<Scalar> retVals = ret->getDataNonConst(0);
281 ArrayRCP<const Scalar> inputVals = v->getData(0);
282 for (
size_t i = 0; i < v->getMap()->getNodeNumElements(); ++i) {
283 if(Teuchos::ScalarTraits<Scalar>::magnitude(inputVals[i]) > tol)
284 retVals[i] = Teuchos::ScalarTraits<Scalar>::one() / inputVals[i];
286 retVals[i] = tolReplacement;
299 RCP<const Map> rowMap = A.getRowMap(), colMap = A.getColMap();
300 RCP<Vector> localDiag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap);
303 const CrsMatrixWrap* crsOp =
dynamic_cast<const CrsMatrixWrap*
>(&A);
307 Teuchos::ArrayRCP<size_t> offsets;
308 crsOp->getLocalDiagOffsets(offsets);
309 crsOp->getLocalDiagCopy(*localDiag,offsets());
312 ArrayRCP<Scalar> localDiagVals = localDiag->getDataNonConst(0);
314 for (LocalOrdinal i = 0; i < localDiagVals.size(); i++)
315 localDiagVals[i] = diagVals[i];
316 localDiagVals = diagVals = null;
319 RCP<Vector> diagonal = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(colMap);
320 RCP< const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > importer;
321 importer = A.getCrsGraph()->getImporter();
322 if (importer == Teuchos::null) {
323 importer = Xpetra::ImportFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap, colMap);
325 diagonal->doImport(*localDiag, *(importer), Xpetra::INSERT);
334 static Teuchos::Array<Magnitude>
ResidualNorm(
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const MultiVector& X,
const MultiVector& RHS) {
335 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(),
Exceptions::RuntimeError,
"Number of solution vectors != number of right-hand sides")
336 const size_t numVecs = X.getNumVectors();
337 RCP<MultiVector> RES =
Residual(Op, X, RHS);
338 Teuchos::Array<Magnitude> norms(numVecs);
343 static RCP<MultiVector>
Residual(
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op,
const MultiVector& X,
const MultiVector& RHS) {
344 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(),
Exceptions::RuntimeError,
"Number of solution vectors != number of right-hand sides")
345 const size_t numVecs = X.getNumVectors();
346 Scalar one = Teuchos::ScalarTraits<Scalar>::one(), negone = -one, zero = Teuchos::ScalarTraits<Scalar>::zero();
348 RCP<MultiVector> RES = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(RHS.getMap(), numVecs,
false);
349 Op.apply(X, *RES, Teuchos::NO_TRANS, one, zero);
350 RES->update(one, RHS, negone);
357 RCP<const Teuchos::Comm<int> > comm = Teuchos::DefaultComm<int>::getComm();
358 int myPID = comm->getRank();
361 for (
int i = 0; i <comm->getSize(); i++) {
363 gethostname(hostname,
sizeof(hostname));
364 std::cout <<
"Host: " << hostname <<
"\tMPI rank: " << myPID <<
",\tPID: " << pid <<
"\n\tattach " << pid << std::endl;
369 std::cout <<
"** Enter a character to continue > " << std::endl;
371 int r = scanf(
"%c", &go);
398 static Scalar
PowerMethod(
const Matrix& A,
bool scaleByDiag =
true,
399 LocalOrdinal niters = 10, Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
401 "Utils::PowerMethod: operator must have domain and range maps that are equivalent.");
404 RCP<Vector> q = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getDomainMap());
405 RCP<Vector> r = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRangeMap());
406 RCP<Vector> z = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRangeMap());
411 Teuchos::Array<Magnitude> norms(1);
413 typedef Teuchos::ScalarTraits<Scalar> STS;
415 const Scalar zero = STS::zero(), one = STS::one();
417 Scalar lambda = zero;
418 Magnitude residual = STS::magnitude(zero);
421 RCP<Vector> diagInvVec;
423 RCP<Vector> diagVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
424 A.getLocalDiagCopy(*diagVec);
425 diagInvVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
426 diagInvVec->reciprocal(*diagVec);
429 for (
int iter = 0; iter < niters; ++iter) {
431 q->update(one/norms[0], *z, zero);
434 z->elementWiseMultiply(one, *diagInvVec, *z, zero);
437 if (iter % 100 == 0 || iter + 1 == niters) {
438 r->update(1.0, *z, -lambda, *q, zero);
440 residual = STS::magnitude(norms[0] / lambda);
442 std::cout <<
"Iter = " << iter
443 <<
" Lambda = " << lambda
444 <<
" Residual of A*q - lambda*q = " << residual
448 if (residual < tolerance)
456 static RCP<Teuchos::FancyOStream>
MakeFancy(std::ostream& os) {
457 RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(os));
465 static typename Teuchos::ScalarTraits<Scalar>::magnitudeType
Distance2(
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& v, LocalOrdinal i0, LocalOrdinal i1) {
466 size_t numVectors = v.getNumVectors();
468 Scalar d = Teuchos::ScalarTraits<Scalar>::zero();
469 for (
size_t j = 0; j < numVectors; j++) {
470 Teuchos::ArrayRCP<const Scalar> vv = v.getData(j);
471 d += (vv[i0] - vv[i1])*(vv[i0] - vv[i1]);
473 return Teuchos::ScalarTraits<Scalar>::magnitude(d);
488 static Teuchos::ArrayRCP<const bool>
DetectDirichletRows(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
const Magnitude& tol = Teuchos::ScalarTraits<Scalar>::zero()) {
489 LocalOrdinal numRows = A.getNodeNumRows();
490 typedef Teuchos::ScalarTraits<Scalar> STS;
491 ArrayRCP<bool> boundaryNodes(numRows,
true);
492 for (LocalOrdinal row = 0; row < numRows; row++) {
493 ArrayView<const LocalOrdinal> indices;
494 ArrayView<const Scalar> vals;
495 A.getLocalRowView(row, indices, vals);
496 size_t nnz = A.getNumEntriesInLocalRow(row);
498 for (
size_t col = 0; col < nnz; col++)
499 if ( (indices[col] != row) && STS::magnitude(vals[col]) > tol) {
500 boundaryNodes[row] =
false;
504 return boundaryNodes;
520 static Teuchos::ArrayRCP<const bool>
DetectDirichletRowsExt(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
bool & bHasZeroDiagonal,
const Magnitude& tol = Teuchos::ScalarTraits<Scalar>::zero()) {
523 bHasZeroDiagonal =
false;
525 Teuchos::RCP<Vector> diagVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
526 A.getLocalDiagCopy(*diagVec);
527 Teuchos::ArrayRCP< const Scalar > diagVecData = diagVec->getData(0);
529 LocalOrdinal numRows = A.getNodeNumRows();
530 typedef Teuchos::ScalarTraits<Scalar> STS;
531 ArrayRCP<bool> boundaryNodes(numRows,
false);
532 for (LocalOrdinal row = 0; row < numRows; row++) {
533 ArrayView<const LocalOrdinal> indices;
534 ArrayView<const Scalar> vals;
535 A.getLocalRowView(row, indices, vals);
537 bool bHasDiag =
false;
538 for (decltype(indices.size()) col = 0; col < indices.size(); col++) {
539 if ( indices[col] != row) {
540 if (STS::magnitude(vals[col] / sqrt(STS::magnitude(diagVecData[row]) * STS::magnitude(diagVecData[col])) ) > tol) {
543 }
else bHasDiag =
true;
545 if (bHasDiag ==
false) bHasZeroDiagonal =
true;
546 else if(nnz == 0) boundaryNodes[row] =
true;
548 return boundaryNodes;
555 static Scalar
Frobenius(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& B) {
560 TEUCHOS_TEST_FOR_EXCEPTION(!A.getRowMap()->isSameAs(*B.getRowMap()),
Exceptions::Incompatible,
"MueLu::CGSolver::Frobenius: row maps are incompatible");
561 TEUCHOS_TEST_FOR_EXCEPTION(!A.isFillComplete() || !B.isFillComplete(),
Exceptions::RuntimeError,
"Matrices must be fill completed");
563 const Map& AColMap = *A.getColMap();
564 const Map& BColMap = *B.getColMap();
566 Teuchos::ArrayView<const LocalOrdinal> indA, indB;
567 Teuchos::ArrayView<const Scalar> valA, valB;
568 size_t nnzA = 0, nnzB = 0;
580 Teuchos::Array<Scalar> valBAll(BColMap.getNodeNumElements());
582 LocalOrdinal invalid = Teuchos::OrdinalTraits<LocalOrdinal>::invalid();
583 Scalar zero = Teuchos::ScalarTraits<Scalar> ::zero(), f = zero, gf;
584 size_t numRows = A.getNodeNumRows();
585 for (
size_t i = 0; i < numRows; i++) {
586 A.getLocalRowView(i, indA, valA);
587 B.getLocalRowView(i, indB, valB);
592 for (
size_t j = 0; j < nnzB; j++)
593 valBAll[indB[j]] = valB[j];
595 for (
size_t j = 0; j < nnzA; j++) {
598 LocalOrdinal ind = BColMap.getLocalElement(AColMap.getGlobalElement(indA[j]));
600 f += valBAll[ind] * valA[j];
604 for (
size_t j = 0; j < nnzB; j++)
605 valBAll[indB[j]] = zero;
628 int maxint = INT_MAX;
629 int mySeed = Teuchos::as<int>((maxint-1) * (one -(comm.getRank()+1)/(comm.getSize()+one)) );
630 if (mySeed < 1 || mySeed == maxint) {
631 std::ostringstream errStr;
632 errStr <<
"Error detected with random seed = " << mySeed <<
". It should be in the interval [1,2^31-2].";
637 Teuchos::ScalarTraits<Scalar>::seedrandom(mySeed);
649 static void FindDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &A,
650 std::vector<LocalOrdinal>& dirichletRows,
bool count_twos_as_dirichlet=
false) {
651 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
652 dirichletRows.resize(0);
653 for(
size_t i=0; i<A->getNodeNumRows(); i++) {
654 Teuchos::ArrayView<const LocalOrdinal> indices;
655 Teuchos::ArrayView<const Scalar> values;
656 A->getLocalRowView(i,indices,values);
658 for (
size_t j=0; j<(size_t)indices.size(); j++) {
659 if (Teuchos::ScalarTraits<Scalar>::magnitude(values[j]) > Teuchos::ScalarTraits<MT>::eps()) {
663 if (nnz == 1 || (count_twos_as_dirichlet && nnz == 2)) {
664 dirichletRows.push_back(i);
671 const std::vector<LocalOrdinal>& dirichletRows) {
672 RCP<const Map> Rmap = A->getColMap();
673 RCP<const Map> Cmap = A->getColMap();
674 Scalar one =Teuchos::ScalarTraits<Scalar>::one();
675 Scalar zero =Teuchos::ScalarTraits<Scalar>::zero();
677 for(
size_t i=0; i<dirichletRows.size(); i++) {
678 GlobalOrdinal row_gid = Rmap->getGlobalElement(dirichletRows[i]);
680 Teuchos::ArrayView<const LocalOrdinal> indices;
681 Teuchos::ArrayView<const Scalar> values;
682 A->getLocalRowView(dirichletRows[i],indices,values);
684 Scalar* valuesNC =
const_cast<Scalar*
>(values.getRawPtr());
685 for(
size_t j=0; j<(size_t)indices.size(); j++) {
686 if(Cmap->getGlobalElement(indices[j])==row_gid)
695 static void ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
696 const std::vector<LocalOrdinal>& dirichletRows) {
697 Scalar zero =Teuchos::ScalarTraits<Scalar>::zero();
699 for(
size_t i=0; i<dirichletRows.size(); i++) {
700 Teuchos::ArrayView<const LocalOrdinal> indices;
701 Teuchos::ArrayView<const Scalar> values;
702 A->getLocalRowView(dirichletRows[i],indices,values);
704 Scalar* valuesNC =
const_cast<Scalar*
>(values.getRawPtr());
705 for(
size_t j=0; j<(size_t)indices.size(); j++)
712 Teuchos::RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> >& isDirichletRow,
713 Teuchos::RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> >& isDirichletCol) {
716 if(!A->getRowMap()->isSameAs(*A->getDomainMap())) {
717 throw std::runtime_error(
"UtilitiesBase::FindDirichletRowsAndPropagateToCols row and domain maps must match.");
719 RCP<const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > importer = A->getCrsGraph()->getImporter();
720 bool has_import = !importer.is_null();
723 std::vector<LocalOrdinal> dirichletRows;
727 printf(
"[%d] DirichletRow Ids = ",A->getRowMap()->getComm()->getRank());
728 for(
size_t i=0; i<(size_t) dirichletRows.size(); i++)
729 printf(
"%d ",dirichletRows[i]);
734 isDirichletRow = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getRowMap(),
true);
735 isDirichletCol = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getColMap(),
true);
737 Teuchos::ArrayRCP<int> dr_rcp = isDirichletRow->getDataNonConst(0);
738 Teuchos::ArrayView<int> dr = dr_rcp();
739 Teuchos::ArrayRCP<int> dc_rcp = isDirichletCol->getDataNonConst(0);
740 Teuchos::ArrayView<int> dc = dc_rcp();
741 for(
size_t i=0; i<(size_t) dirichletRows.size(); i++) {
742 dr[dirichletRows[i]] = 1;
743 if(!has_import) dc[dirichletRows[i]] = 1;
747 isDirichletCol->doImport(*isDirichletRow,*importer,Xpetra::CombineMode::ADD);
759 #define MUELU_UTILITIESBASE_SHORT 760 #endif // MUELU_UTILITIESBASE_DECL_HPP Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > BlockedMultiVector
static Scalar PowerMethod(const Matrix &A, bool scaleByDiag=true, LocalOrdinal niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
Simple transpose for Tpetra::CrsMatrix types.
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MultiVector
static void ApplyOAZToMatrixRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const std::vector< LocalOrdinal > &dirichletRows)
static Teuchos::RCP< Vector > GetInverse(Teuchos::RCP< const Vector > v, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar tolReplacement=Teuchos::ScalarTraits< Scalar >::zero())
Return vector containing inverse of input vector.
static void FindDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, std::vector< LocalOrdinal > &dirichletRows, bool count_twos_as_dirichlet=false)
static Teuchos::ScalarTraits< Scalar >::magnitudeType Distance2(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &v, LocalOrdinal i0, LocalOrdinal i1)
Squared distance between two rows in a multivector.
Namespace for MueLu classes and methods.
Exception throws to report incompatible objects (like maps).
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 PauseForDebugger()
#define MueLu_sumAll(rcpComm, in, out)
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)
Xpetra::BlockedVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > BlockedVector
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
static Teuchos::ArrayRCP< Scalar > GetLumpedMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal of lumped matrix.
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.
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > Vector
Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > BlockedMap
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< Matrix > Crs2Op(RCP< CrsMatrix > Op)
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixWrap
static RCP< Vector > GetMatrixDiagonalInverse(const Matrix &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100)
Extract Matrix Diagonal.
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::zero())
Detect Dirichlet rows.
Exception throws to report errors in the internal logical of the program.
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
Extract Overlapped Matrix Diagonal.
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const std::vector< LocalOrdinal > &dirichletRows)
static Teuchos::RCP< Vector > GetLumpedMatrixDiagonal(Teuchos::RCP< const Matrix > rcpA)
Extract Matrix Diagonal of lumped matrix.
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
static Teuchos::ArrayRCP< Scalar > GetMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal.