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;
126 static RCP<Matrix>
Crs2Op(RCP<CrsMatrix> Op) {
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) {
144 A.getLocalRowView(i, cols, vals);
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();
168 RCP<const Matrix> rcpA = Teuchos::rcpFromRef(A);
170 RCP<const Map> rowMap = rcpA->getRowMap();
171 RCP<Vector> diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,
true);
173 rcpA->getLocalDiagCopy(*diag);
190 size_t numRows = A.getRowMap()->getNodeNumElements();
191 Teuchos::ArrayRCP<Scalar> diag(numRows);
192 Teuchos::ArrayView<const LocalOrdinal> cols;
193 Teuchos::ArrayView<const Scalar> vals;
194 for (
size_t i = 0; i < numRows; ++i) {
195 A.getLocalRowView(i, cols, vals);
196 diag[i] = Teuchos::ScalarTraits<Scalar>::zero();
198 diag[i] += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
212 RCP<Vector> diag = Teuchos::null;
214 RCP<const Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > bA =
215 Teuchos::rcp_dynamic_cast<const Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(rcpA);
216 if(bA == Teuchos::null) {
217 RCP<const Map> rowMap = rcpA->getRowMap();
218 diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,
true);
219 ArrayRCP<Scalar> diagVals = diag->getDataNonConst(0);
220 Teuchos::ArrayView<const LocalOrdinal> cols;
221 Teuchos::ArrayView<const Scalar> vals;
222 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
223 rcpA->getLocalRowView(i, cols, vals);
224 diagVals[i] = Teuchos::ScalarTraits<Scalar>::zero();
226 diagVals[i] += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
234 diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(bA->getRangeMapExtractor()->getFullMap(),
true);
236 for (
size_t row = 0; row < bA->Rows(); ++row) {
237 for (
size_t col = 0; col < bA->Cols(); ++col) {
238 if (!bA->getMatrix(row,col).is_null()) {
240 bool bThyraMode = bA->getRangeMapExtractor()->getThyraMode() && (Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(bA->getMatrix(row,col)) == Teuchos::null);
241 RCP<Vector> ddtemp = bA->getRangeMapExtractor()->ExtractVector(diag,row,bThyraMode);
243 ddtemp->update(Teuchos::as<Scalar>(1.0),*dd,Teuchos::as<Scalar>(1.0));
244 bA->getRangeMapExtractor()->InsertVector(ddtemp,row,diag,bThyraMode);
262 static Teuchos::RCP<Vector>
GetInverse(Teuchos::RCP<const Vector> v,
Magnitude tol = Teuchos::ScalarTraits<Scalar>::eps()*100,
Scalar tolReplacement = Teuchos::ScalarTraits<Scalar>::zero()) {
264 RCP<Vector> ret = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(v->getMap(),
true);
267 RCP<const BlockedVector> bv = Teuchos::rcp_dynamic_cast<const BlockedVector>(v);
268 if(bv.is_null() ==
false) {
269 RCP<BlockedVector> bret = Teuchos::rcp_dynamic_cast<BlockedVector>(ret);
270 TEUCHOS_TEST_FOR_EXCEPTION(bret.is_null() ==
true,
MueLu::Exceptions::RuntimeError,
"MueLu::UtilitiesBase::GetInverse: return vector should be of type BlockedVector");
271 RCP<const BlockedMap> bmap = bv->getBlockedMap();
272 for(
size_t r = 0; r < bmap->getNumMaps(); ++r) {
273 RCP<const MultiVector> submvec = bv->getMultiVector(r,bmap->getThyraMode());
274 RCP<const Vector> subvec = submvec->getVector(0);
276 bret->setMultiVector(r, subvecinf, bmap->getThyraMode());
282 ArrayRCP<Scalar> retVals = ret->getDataNonConst(0);
283 ArrayRCP<const Scalar> inputVals = v->getData(0);
284 for (
size_t i = 0; i < v->getMap()->getNodeNumElements(); ++i) {
285 if(Teuchos::ScalarTraits<Scalar>::magnitude(inputVals[i]) > tol)
286 retVals[i] = Teuchos::ScalarTraits<Scalar>::one() / inputVals[i];
288 retVals[i] = tolReplacement;
301 RCP<const Map> rowMap = A.getRowMap(), colMap = A.getColMap();
304 RCP<const BlockedMap> browMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(rowMap);
305 if(!browMap.is_null()) rowMap = browMap->getMap();
307 RCP<Vector> localDiag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap);
313 Teuchos::ArrayRCP<size_t> offsets;
314 crsOp->getLocalDiagOffsets(offsets);
315 crsOp->getLocalDiagCopy(*localDiag,offsets());
318 ArrayRCP<Scalar> localDiagVals = localDiag->getDataNonConst(0);
321 localDiagVals[i] = diagVals[i];
322 localDiagVals = diagVals =
null;
325 RCP<Vector> diagonal = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(colMap);
326 RCP< const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > importer;
327 importer = A.getCrsGraph()->getImporter();
328 if (importer == Teuchos::null) {
329 importer = Xpetra::ImportFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap, colMap);
331 diagonal->doImport(*localDiag, *(importer), Xpetra::INSERT);
341 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(),
Exceptions::RuntimeError,
"Number of solution vectors != number of right-hand sides")
342 const size_t numVecs = X.getNumVectors();
343 RCP<MultiVector> RES =
Residual(Op, X, RHS);
344 Teuchos::Array<Magnitude> norms(numVecs);
350 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(),
Exceptions::RuntimeError,
"Number of solution vectors != number of right-hand sides")
351 const size_t numVecs = X.getNumVectors();
353 Teuchos::Array<Magnitude> norms(numVecs);
359 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(),
Exceptions::RuntimeError,
"Number of solution vectors != number of right-hand sides")
360 const size_t numVecs = X.getNumVectors();
361 Scalar one = Teuchos::ScalarTraits<Scalar>::one(), negone = -one, zero = Teuchos::ScalarTraits<Scalar>::zero();
363 RCP<MultiVector> RES = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(RHS.getMap(), numVecs,
false);
364 Op.apply(X, *RES, Teuchos::NO_TRANS, one, zero);
365 RES->update(one, RHS, negone);
371 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(),
Exceptions::RuntimeError,
"Number of solution vectors != number of right-hand sides");
372 TEUCHOS_TEST_FOR_EXCEPTION(Resid.getNumVectors() != RHS.getNumVectors(),
Exceptions::RuntimeError,
"Number of residual vectors != number of right-hand sides");
373 Scalar one = Teuchos::ScalarTraits<Scalar>::one(), negone = -one, zero = Teuchos::ScalarTraits<Scalar>::zero();
375 Op.apply(X, Resid, Teuchos::NO_TRANS, one, zero);
376 Resid.update(one, RHS, negone);
382 RCP<const Teuchos::Comm<int> > comm = Teuchos::DefaultComm<int>::getComm();
383 int myPID = comm->getRank();
386 for (
int i = 0; i <comm->getSize(); i++) {
388 gethostname(hostname,
sizeof(hostname));
389 std::cout <<
"Host: " << hostname <<
"\tMPI rank: " << myPID <<
",\tPID: " << pid <<
"\n\tattach " << pid << std::endl;
394 std::cout <<
"** Enter a character to continue > " << std::endl;
396 int r = scanf(
"%c", &go);
424 LocalOrdinal niters = 10,
Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
426 "Utils::PowerMethod: operator must have domain and range maps that are equivalent.");
429 RCP<Vector> q = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getDomainMap());
430 RCP<Vector> r = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRangeMap());
431 RCP<Vector> z = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRangeMap());
436 Teuchos::Array<Magnitude> norms(1);
438 typedef Teuchos::ScalarTraits<Scalar> STS;
440 const Scalar zero = STS::zero(), one = STS::one();
443 Magnitude residual = STS::magnitude(zero);
446 RCP<Vector> diagInvVec;
448 RCP<Vector> diagVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
449 A.getLocalDiagCopy(*diagVec);
450 diagInvVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
451 diagInvVec->reciprocal(*diagVec);
454 for (
int iter = 0; iter < niters; ++iter) {
456 q->update(one/norms[0], *z, zero);
459 z->elementWiseMultiply(one, *diagInvVec, *z, zero);
462 if (iter % 100 == 0 || iter + 1 == niters) {
463 r->update(1.0, *z, -lambda, *q, zero);
465 residual = STS::magnitude(norms[0] / lambda);
467 std::cout <<
"Iter = " << iter
468 <<
" Lambda = " << lambda
469 <<
" Residual of A*q - lambda*q = " << residual
481 static RCP<Teuchos::FancyOStream>
MakeFancy(std::ostream& os) {
482 RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(os));
491 const size_t numVectors = v.size();
493 Scalar d = Teuchos::ScalarTraits<Scalar>::zero();
494 for (
size_t j = 0; j < numVectors; j++) {
495 d += (v[j][i0] - v[j][i1])*(v[j][i0] - v[j][i1]);
497 return Teuchos::ScalarTraits<Scalar>::magnitude(d);
512 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) {
514 typedef Teuchos::ScalarTraits<Scalar> STS;
515 ArrayRCP<bool> boundaryNodes(numRows,
true);
516 if (count_twos_as_dirichlet) {
518 ArrayView<const LocalOrdinal> indices;
519 ArrayView<const Scalar> vals;
520 A.getLocalRowView(row, indices, vals);
521 size_t nnz = A.getNumEntriesInLocalRow(row);
524 for (col = 0; col < nnz; col++)
525 if ( (indices[col] != row) && STS::magnitude(vals[col]) > tol) {
526 if (!boundaryNodes[row])
528 boundaryNodes[row] =
false;
531 boundaryNodes[row] =
true;
536 ArrayView<const LocalOrdinal> indices;
537 ArrayView<const Scalar> vals;
538 A.getLocalRowView(row, indices, vals);
539 size_t nnz = A.getNumEntriesInLocalRow(row);
541 for (
size_t col = 0; col < nnz; col++)
542 if ( (indices[col] != row) && STS::magnitude(vals[col]) > tol) {
543 boundaryNodes[row] =
false;
548 return boundaryNodes;
564 static Teuchos::ArrayRCP<const bool>
DetectDirichletRowsExt(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
bool & bHasZeroDiagonal,
const Magnitude& tol = Teuchos::ScalarTraits<Scalar>::zero()) {
567 bHasZeroDiagonal =
false;
569 Teuchos::RCP<Vector> diagVec = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A.getRowMap());
570 A.getLocalDiagCopy(*diagVec);
571 Teuchos::ArrayRCP< const Scalar > diagVecData = diagVec->getData(0);
574 typedef Teuchos::ScalarTraits<Scalar> STS;
575 ArrayRCP<bool> boundaryNodes(numRows,
false);
577 ArrayView<const LocalOrdinal> indices;
578 ArrayView<const Scalar> vals;
579 A.getLocalRowView(row, indices, vals);
581 bool bHasDiag =
false;
582 for (decltype(indices.size()) col = 0; col < indices.size(); col++) {
583 if ( indices[col] != row) {
584 if (STS::magnitude(vals[col] / STS::magnitude(sqrt(STS::magnitude(diagVecData[row]) * STS::magnitude(diagVecData[col]))) ) > tol) {
587 }
else bHasDiag =
true;
589 if (bHasDiag ==
false) bHasZeroDiagonal =
true;
590 else if(nnz == 0) boundaryNodes[row] =
true;
592 return boundaryNodes;
606 static Teuchos::ArrayRCP<const bool>
DetectDirichletCols(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
607 const Teuchos::ArrayRCP<const bool>& dirichletRows) {
608 Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
609 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
610 Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > domMap = A.getDomainMap();
611 Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > colMap = A.getColMap();
612 Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > myColsToZero = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(colMap,1);
613 myColsToZero->putScalar(zero);
615 for(
size_t i=0; i<(size_t) dirichletRows.size(); i++) {
616 if (dirichletRows[i]) {
617 Teuchos::ArrayView<const LocalOrdinal> indices;
618 Teuchos::ArrayView<const Scalar> values;
619 A.getLocalRowView(i,indices,values);
620 for(
size_t j=0; j<static_cast<size_t>(indices.size()); j++)
621 myColsToZero->replaceLocalValue(indices[j],0,one);
625 Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > globalColsToZero = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(domMap,1);
626 globalColsToZero->putScalar(zero);
627 Teuchos::RCP<Xpetra::Export<LocalOrdinal,GlobalOrdinal,Node> > exporter = Xpetra::ExportFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(colMap,domMap);
629 globalColsToZero->doExport(*myColsToZero,*exporter,Xpetra::ADD);
631 myColsToZero->doImport(*globalColsToZero,*exporter,Xpetra::INSERT);
632 Teuchos::ArrayRCP<const Scalar> myCols = myColsToZero->getData(0);
633 Teuchos::ArrayRCP<bool> dirichletCols(colMap->getNodeNumElements(),
true);
634 Magnitude eps = Teuchos::ScalarTraits<Magnitude>::eps();
635 for(
size_t i=0; i<colMap->getNodeNumElements(); i++) {
636 dirichletCols[i] = Teuchos::ScalarTraits<Scalar>::magnitude(myCols[i])>2.0*eps;
638 return dirichletCols;
646 static Scalar Frobenius(
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A,
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>& B) {
651 TEUCHOS_TEST_FOR_EXCEPTION(!A.getRowMap()->isSameAs(*B.getRowMap()),
Exceptions::Incompatible,
"MueLu::CGSolver::Frobenius: row maps are incompatible");
652 TEUCHOS_TEST_FOR_EXCEPTION(!A.isFillComplete() || !B.isFillComplete(),
Exceptions::RuntimeError,
"Matrices must be fill completed");
654 const Map& AColMap = *A.getColMap();
655 const Map& BColMap = *B.getColMap();
657 Teuchos::ArrayView<const LocalOrdinal> indA, indB;
658 Teuchos::ArrayView<const Scalar> valA, valB;
659 size_t nnzA = 0, nnzB = 0;
671 Teuchos::Array<Scalar> valBAll(BColMap.getNodeNumElements());
673 LocalOrdinal invalid = Teuchos::OrdinalTraits<LocalOrdinal>::invalid();
674 Scalar zero = Teuchos::ScalarTraits<Scalar> ::zero(), f = zero, gf;
675 size_t numRows = A.getNodeNumRows();
676 for (
size_t i = 0; i < numRows; i++) {
677 A.getLocalRowView(i, indA, valA);
678 B.getLocalRowView(i, indB, valB);
683 for (
size_t j = 0; j < nnzB; j++)
684 valBAll[indB[j]] = valB[j];
686 for (
size_t j = 0; j < nnzA; j++) {
689 LocalOrdinal ind = BColMap.getLocalElement(AColMap.getGlobalElement(indA[j]));
691 f += valBAll[ind] * valA[j];
695 for (
size_t j = 0; j < nnzB; j++)
696 valBAll[indB[j]] = zero;
719 int maxint = INT_MAX;
720 int mySeed = Teuchos::as<int>((maxint-1) * (one -(comm.getRank()+1)/(comm.getSize()+one)) );
721 if (mySeed < 1 || mySeed == maxint) {
722 std::ostringstream errStr;
723 errStr <<
"Error detected with random seed = " << mySeed <<
". It should be in the interval [1,2^31-2].";
728 Teuchos::ScalarTraits<Scalar>::seedrandom(mySeed);
740 static void FindDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &A,
741 std::vector<LocalOrdinal>& dirichletRows,
bool count_twos_as_dirichlet=
false) {
742 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
743 dirichletRows.resize(0);
744 for(
size_t i=0; i<A->getNodeNumRows(); i++) {
745 Teuchos::ArrayView<const LocalOrdinal> indices;
746 Teuchos::ArrayView<const Scalar> values;
747 A->getLocalRowView(i,indices,values);
749 for (
size_t j=0; j<(size_t)indices.size(); j++) {
750 if (Teuchos::ScalarTraits<Scalar>::magnitude(values[j]) > Teuchos::ScalarTraits<MT>::eps()) {
754 if (nnz == 1 || (count_twos_as_dirichlet && nnz == 2)) {
755 dirichletRows.push_back(i);
763 const std::vector<LocalOrdinal>& dirichletRows) {
764 RCP<const Map> Rmap = A->getColMap();
765 RCP<const Map> Cmap = A->getColMap();
766 Scalar one =Teuchos::ScalarTraits<Scalar>::one();
767 Scalar zero =Teuchos::ScalarTraits<Scalar>::zero();
769 for(
size_t i=0; i<dirichletRows.size(); i++) {
770 GlobalOrdinal row_gid = Rmap->getGlobalElement(dirichletRows[i]);
772 Teuchos::ArrayView<const LocalOrdinal> indices;
773 Teuchos::ArrayView<const Scalar> values;
774 A->getLocalRowView(dirichletRows[i],indices,values);
776 Scalar* valuesNC =
const_cast<Scalar*
>(values.getRawPtr());
777 for(
size_t j=0; j<(size_t)indices.size(); j++) {
778 if(Cmap->getGlobalElement(indices[j])==row_gid)
789 const Teuchos::ArrayRCP<const bool>& dirichletRows) {
790 RCP<const Map> Rmap = A->getColMap();
791 RCP<const Map> Cmap = A->getColMap();
792 Scalar one =Teuchos::ScalarTraits<Scalar>::one();
793 Scalar zero =Teuchos::ScalarTraits<Scalar>::zero();
795 for(
size_t i=0; i<(size_t) dirichletRows.size(); i++) {
796 if (dirichletRows[i]){
799 Teuchos::ArrayView<const LocalOrdinal> indices;
800 Teuchos::ArrayView<const Scalar> values;
801 A->getLocalRowView(i,indices,values);
803 Scalar* valuesNC =
const_cast<Scalar*
>(values.getRawPtr());
804 for(
size_t j=0; j<(size_t)indices.size(); j++) {
805 if(Cmap->getGlobalElement(indices[j])==row_gid)
816 static void ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
817 const std::vector<LocalOrdinal>& dirichletRows,
818 Scalar replaceWith=Teuchos::ScalarTraits<Scalar>::zero()) {
819 for(
size_t i=0; i<dirichletRows.size(); i++) {
820 Teuchos::ArrayView<const LocalOrdinal> indices;
821 Teuchos::ArrayView<const Scalar> values;
822 A->getLocalRowView(dirichletRows[i],indices,values);
824 Scalar* valuesNC =
const_cast<Scalar*
>(values.getRawPtr());
825 for(
size_t j=0; j<(size_t)indices.size(); j++)
826 valuesNC[j]=replaceWith;
832 static void ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& A,
833 const Teuchos::ArrayRCP<const bool>& dirichletRows,
834 Scalar replaceWith=Teuchos::ScalarTraits<Scalar>::zero()) {
835 for(
size_t i=0; i<(size_t) dirichletRows.size(); i++) {
836 if (dirichletRows[i]) {
837 Teuchos::ArrayView<const LocalOrdinal> indices;
838 Teuchos::ArrayView<const Scalar> values;
839 A->getLocalRowView(i,indices,values);
841 Scalar* valuesNC =
const_cast<Scalar*
>(values.getRawPtr());
842 for(
size_t j=0; j<(size_t)indices.size(); j++)
843 valuesNC[j]=replaceWith;
850 static void ZeroDirichletRows(Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& X,
851 const Teuchos::ArrayRCP<const bool>& dirichletRows,
852 Scalar replaceWith=Teuchos::ScalarTraits<Scalar>::zero()) {
853 for(
size_t i=0; i<(size_t) dirichletRows.size(); i++) {
854 if (dirichletRows[i]) {
855 for(
size_t j=0; j<X->getNumVectors(); j++)
856 X->replaceLocalValue(i,j,replaceWith);
864 const Teuchos::ArrayRCP<const bool>& dirichletCols,
865 Scalar replaceWith=Teuchos::ScalarTraits<Scalar>::zero()) {
866 for(
size_t i=0; i<A->getNodeNumRows(); i++) {
867 Teuchos::ArrayView<const LocalOrdinal> indices;
868 Teuchos::ArrayView<const Scalar> values;
869 A->getLocalRowView(i,indices,values);
871 Scalar* valuesNC =
const_cast<Scalar*
>(values.getRawPtr());
872 for(
size_t j=0; j<static_cast<size_t>(indices.size()); j++)
873 if (dirichletCols[indices[j]])
874 valuesNC[j] = replaceWith;
880 Teuchos::RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> >& isDirichletRow,
881 Teuchos::RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> >& isDirichletCol) {
884 if(!A->getRowMap()->isSameAs(*A->getDomainMap())) {
885 throw std::runtime_error(
"UtilitiesBase::FindDirichletRowsAndPropagateToCols row and domain maps must match.");
887 RCP<const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > importer = A->getCrsGraph()->getImporter();
888 bool has_import = !importer.is_null();
891 std::vector<LocalOrdinal> dirichletRows;
895 printf(
"[%d] DirichletRow Ids = ",A->getRowMap()->getComm()->getRank());
896 for(
size_t i=0; i<(size_t) dirichletRows.size(); i++)
897 printf(
"%d ",dirichletRows[i]);
902 isDirichletRow = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getRowMap(),
true);
903 isDirichletCol = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getColMap(),
true);
905 Teuchos::ArrayRCP<int> dr_rcp = isDirichletRow->getDataNonConst(0);
906 Teuchos::ArrayView<int> dr = dr_rcp();
907 Teuchos::ArrayRCP<int> dc_rcp = isDirichletCol->getDataNonConst(0);
908 Teuchos::ArrayView<int> dc = dc_rcp();
909 for(
size_t i=0; i<(size_t) dirichletRows.size(); i++) {
910 dr[dirichletRows[i]] = 1;
911 if(!has_import) dc[dirichletRows[i]] = 1;
915 isDirichletCol->doImport(*isDirichletRow,*importer,Xpetra::CombineMode::ADD);
922 static RCP<const Xpetra::BlockedMap<LocalOrdinal,GlobalOrdinal,Node> >
GeneratedBlockedTargetMap(
const Xpetra::BlockedMap<LocalOrdinal,GlobalOrdinal,Node> & sourceBlockedMap,
923 const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> & Importer) {
924 typedef Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> IntVector;
925 Xpetra::UnderlyingLib lib = sourceBlockedMap.lib();
928 RCP<const Map> fullMap = sourceBlockedMap.getMap();
929 RCP<const Map> stridedMap = Teuchos::rcp_dynamic_cast<const Xpetra::StridedMap<LocalOrdinal,GlobalOrdinal,Node> >(fullMap);
930 if(!stridedMap.is_null()) fullMap = stridedMap->getMap();
933 const size_t numSubMaps = sourceBlockedMap.getNumMaps();
934 if(!Importer.getSourceMap()->isCompatible(*fullMap))
935 throw std::runtime_error(
"GenerateBlockedTargetMap(): Map compatibility error");
938 RCP<IntVector> block_ids = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(fullMap);
940 for(
size_t i=0; i<numSubMaps; i++) {
941 RCP<const Map> map = sourceBlockedMap.getMap(i);
943 for(
size_t j=0; j<map->getNodeNumElements(); j++) {
944 LocalOrdinal jj = fullMap->getLocalElement(map->getGlobalElement(j));
945 block_ids->replaceLocalValue(jj,(
int)i);
950 RCP<const Map> targetMap = Importer.getTargetMap();
951 RCP<IntVector> new_block_ids = Xpetra::VectorFactory<int,LocalOrdinal,GlobalOrdinal,Node>::Build(targetMap);
952 new_block_ids->doImport(*block_ids,Importer,Xpetra::CombineMode::ADD);
953 Teuchos::ArrayRCP<const int> dataRCP = new_block_ids->getData(0);
954 Teuchos::ArrayView<const int> data = dataRCP();
958 Teuchos::Array<Teuchos::Array<GlobalOrdinal> > elementsInSubMap(numSubMaps);
959 for(
size_t i=0; i<targetMap->getNodeNumElements(); i++) {
960 elementsInSubMap[data[i]].push_back(targetMap->getGlobalElement(i));
964 std::vector<RCP<const Map> > subMaps(numSubMaps);
965 for(
size_t i=0; i<numSubMaps; i++) {
966 subMaps[i] = Xpetra::MapFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(lib,Teuchos::OrdinalTraits<GlobalOrdinal>::invalid(),elementsInSubMap[i](),targetMap->getIndexBase(),targetMap->getComm());
970 return rcp(
new BlockedMap(targetMap,subMaps));
981 #define MUELU_UTILITIESBASE_SHORT
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
#define MueLu_sumAll(rcpComm, in, out)
Exception throws to report incompatible objects (like maps).
Exception throws to report errors in the internal logical of the program.
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 Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS, MultiVector &Resid)
static void ZeroDirichletCols(Teuchos::RCP< Matrix > &A, const Teuchos::ArrayRCP< const bool > &dirichletCols, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
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 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< Matrix > Crs2Op(RCP< CrsMatrix > Op)
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 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::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > Vector
static void ApplyOAZToMatrixRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const std::vector< LocalOrdinal > &dirichletRows)
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixWrap
static Teuchos::RCP< Vector > GetLumpedMatrixDiagonal(Teuchos::RCP< const Matrix > rcpA)
Extract Matrix Diagonal of lumped matrix.
static void FindDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, std::vector< LocalOrdinal > &dirichletRows, bool count_twos_as_dirichlet=false)
Xpetra::BlockedVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > BlockedVector
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
static RCP< Vector > GetMatrixDiagonalInverse(const Matrix &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100)
Extract Matrix Diagonal.
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 Teuchos::ArrayRCP< Scalar > GetLumpedMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal of lumped matrix.
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.
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static void PauseForDebugger()
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::ArrayRCP< Scalar > GetMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal.
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MultiVector
static RCP< const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > > GeneratedBlockedTargetMap(const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > &sourceBlockedMap, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &Importer)
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > BlockedMultiVector
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 void ApplyOAZToMatrixRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Teuchos::ArrayRCP< const bool > &dirichletRows)
Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > BlockedMap
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
Extract Overlapped Matrix Diagonal.
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())
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.
Namespace for MueLu classes and methods.
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode