46 #ifndef MUELU_UTILITIES_KOKKOS_DECL_HPP 47 #define MUELU_UTILITIES_KOKKOS_DECL_HPP 50 #if defined(HAVE_MUELU_KOKKOS_REFACTOR) 54 #include <Teuchos_ScalarTraits.hpp> 55 #include <Teuchos_ParameterList.hpp> 57 #include <Xpetra_BlockedCrsMatrix_fwd.hpp> 58 #include <Xpetra_CrsMatrix_fwd.hpp> 59 #include <Xpetra_CrsMatrixWrap_fwd.hpp> 60 #include <Xpetra_ExportFactory.hpp> 61 #include <Xpetra_ImportFactory_fwd.hpp> 62 #include <Xpetra_MapFactory_fwd.hpp> 63 #include <Xpetra_Map_fwd.hpp> 64 #include <Xpetra_MatrixFactory_fwd.hpp> 65 #include <Xpetra_Matrix_fwd.hpp> 66 #include <Xpetra_MultiVectorFactory_fwd.hpp> 67 #include <Xpetra_MultiVector_fwd.hpp> 68 #include <Xpetra_Operator_fwd.hpp> 69 #include <Xpetra_VectorFactory_fwd.hpp> 70 #include <Xpetra_Vector_fwd.hpp> 72 #include <Xpetra_IO.hpp> 74 #ifdef HAVE_MUELU_EPETRA 75 #include <Epetra_MultiVector.h> 76 #include <Epetra_CrsMatrix.h> 77 #include <Xpetra_EpetraCrsMatrix_fwd.hpp> 78 #include <Xpetra_EpetraMultiVector_fwd.hpp> 82 #include "MueLu_Utilities.hpp" 83 #include "MueLu_UtilitiesBase.hpp" 85 #ifdef HAVE_MUELU_TPETRA 86 #include <Tpetra_CrsMatrix.hpp> 87 #include <Tpetra_Map.hpp> 88 #include <Tpetra_MultiVector.hpp> 89 #include <Xpetra_TpetraCrsMatrix_fwd.hpp> 90 #include <Xpetra_TpetraMultiVector_fwd.hpp> 103 template <
class Scalar,
104 class LocalOrdinal = int,
105 class GlobalOrdinal = LocalOrdinal,
106 class Node = KokkosClassic::DefaultNode::DefaultNodeType>
108 #undef MUELU_UTILITIES_KOKKOS_SHORT 112 typedef typename Teuchos::ScalarTraits<SC>::magnitudeType Magnitude;
113 typedef Xpetra::MultiVector<Magnitude,LO,GO,NO> RealValuedMultiVector;
115 #ifdef HAVE_MUELU_EPETRA 118 static RCP<const Epetra_MultiVector> MV2EpetraMV(RCP<MultiVector>
const vec) {
return Utilities::MV2EpetraMV(vec); }
121 static const Epetra_MultiVector& MV2EpetraMV(
const MultiVector& vec) {
return Utilities::MV2EpetraMV(vec); }
134 #ifdef HAVE_MUELU_TPETRA 137 static RCP<const Tpetra::MultiVector<SC,LO,GO,NO> > MV2TpetraMV(RCP<MultiVector>
const vec) {
return Utilities::MV2TpetraMV(vec); }
141 static const Tpetra::MultiVector<SC,LO,GO,NO>& MV2TpetraMV(
const MultiVector& vec) {
return Utilities::MV2TpetraMV(vec); }
144 static RCP<const Tpetra::CrsMatrix<SC,LO,GO,NO> > Op2TpetraCrs(RCP<const Matrix> Op) {
return Utilities::Op2TpetraCrs(Op); }
147 static const Tpetra::CrsMatrix<SC,LO,GO,NO>& Op2TpetraCrs(
const Matrix& Op) {
return Utilities::Op2TpetraCrs(Op); }
150 static RCP<const Tpetra::RowMatrix<SC,LO,GO,NO> > Op2TpetraRow(RCP<const Matrix> Op) {
return Utilities::Op2TpetraRow(Op); }
153 static const RCP<const Tpetra::Map<LO, GO, NO> > Map2TpetraMap(
const Map& map) {
return Utilities::Map2TpetraMap(map); }
156 static RCP<Xpetra::Matrix<SC,LO,GO,NO> > Crs2Op(RCP<CrsMatrix> Op) {
return Utilities::Crs2Op(Op); }
164 static Teuchos::ArrayRCP<SC> GetMatrixDiagonal(
const Matrix& A);
172 static RCP<Vector> GetMatrixDiagonalInverse(
const Matrix& A, Magnitude tol = Teuchos::ScalarTraits<SC>::eps()*100);
182 static Teuchos::ArrayRCP<SC> GetLumpedMatrixDiagonal(
const Matrix& A);
193 static RCP<Vector> GetMatrixOverlappedDiagonal(
const Matrix& A);
208 static Teuchos::Array<Magnitude> ResidualNorm(
const Operator& Op,
const MultiVector& X,
const MultiVector& RHS) {
212 static RCP<MultiVector> Residual(
const Operator& Op,
const MultiVector& X,
const MultiVector& RHS) {
216 static void PauseForDebugger();
233 static SC PowerMethod(
const Matrix& A,
bool scaleByDiag =
true,
234 LO niters = 10, Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
238 static void MyOldScaleMatrix(Matrix& Op,
const Teuchos::ArrayRCP<const SC>& scalingVector,
bool doInverse =
true,
239 bool doFillComplete =
true,
bool doOptimizeStorage =
true);
241 static void MyOldScaleMatrix_Tpetra(Matrix& Op,
const Teuchos::ArrayRCP<SC>& scalingVector,
242 bool doFillComplete,
bool doOptimizeStorage);
244 static void MyOldScaleMatrix_Epetra(Matrix& Op,
const Teuchos::ArrayRCP<SC>& scalingVector,
245 bool doFillComplete,
bool doOptimizeStorage);
247 static RCP<Teuchos::FancyOStream> MakeFancy(std::ostream& os) {
return Utilities::MakeFancy(os); }
256 static Kokkos::View<const bool*, typename NO::device_type>
DetectDirichletRows(
const Matrix& A,
const Magnitude& tol = Teuchos::ScalarTraits<SC>::zero(),
const bool count_twos_as_dirichlet=
false);
268 static Kokkos::View<const bool*, typename NO::device_type>
DetectDirichletCols(
const Matrix& A,
const Kokkos::View<const bool*, typename NO::device_type>& dirichletRows);
271 static void ZeroDirichletRows(RCP<Matrix>& A,
const Kokkos::View<const bool*, typename NO::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
273 static void ZeroDirichletRows(RCP<MultiVector>& X,
const Kokkos::View<const bool*, typename NO::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
275 static void ZeroDirichletCols(RCP<Matrix>& A,
const Kokkos::View<const bool*, typename NO::device_type>& dirichletCols, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
277 static RCP<MultiVector> RealValuedToScalarMultiVector(RCP<RealValuedMultiVector> X);
295 static RCP<Matrix> Transpose(Matrix& Op,
bool optimizeTranspose =
false,
const std::string & label = std::string()) {
299 static RCP<Xpetra::MultiVector<double,LocalOrdinal,GlobalOrdinal,Node> > ExtractCoordinatesFromParameterList(ParameterList& paramList) {
315 template <
class Node>
316 class Utilities_kokkos<double,int,int,Node> :
public UtilitiesBase<double,int,int,Node> {
318 typedef double Scalar;
319 typedef int LocalOrdinal;
320 typedef int GlobalOrdinal;
321 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
322 typedef Xpetra::MultiVector<Magnitude,LocalOrdinal,GlobalOrdinal,Node> RealValuedMultiVector;
325 #undef MUELU_UTILITIES_KOKKOS_SHORT 330 #ifdef HAVE_MUELU_EPETRA 333 static RCP<const Epetra_MultiVector> MV2EpetraMV(RCP<MultiVector>
const vec) {
return Utilities::MV2EpetraMV(vec); }
336 static const Epetra_MultiVector& MV2EpetraMV(
const MultiVector& vec) {
return Utilities::MV2EpetraMV(vec); }
349 #ifdef HAVE_MUELU_TPETRA 352 static RCP<const Tpetra::MultiVector<SC,LO,GO,NO> > MV2TpetraMV(RCP<MultiVector>
const vec) {
return Utilities::MV2TpetraMV(vec); }
356 static const Tpetra::MultiVector<SC,LO,GO,NO>& MV2TpetraMV(
const MultiVector& vec) {
return Utilities::MV2TpetraMV(vec); }
359 static RCP<const Tpetra::CrsMatrix<SC,LO,GO,NO> > Op2TpetraCrs(RCP<const Matrix> Op) {
return Utilities::Op2TpetraCrs(Op); }
362 static const Tpetra::CrsMatrix<SC,LO,GO,NO>& Op2TpetraCrs(
const Matrix& Op) {
return Utilities::Op2TpetraCrs(Op); }
365 static RCP<const Tpetra::RowMatrix<SC,LO,GO,NO> > Op2TpetraRow(RCP<const Matrix> Op) {
return Utilities::Op2TpetraRow(Op); }
368 static const RCP<const Tpetra::Map<LO, GO, NO> > Map2TpetraMap(
const Map& map) {
return Utilities::Map2TpetraMap(map); }
370 static RCP<Xpetra::Matrix<SC,LO,GO,NO> > Crs2Op(RCP<CrsMatrix> Op) {
return Utilities::Crs2Op(Op); }
372 static ArrayRCP<SC> GetMatrixDiagonal(
const Matrix& A) {
375 static RCP<Vector> GetMatrixDiagonalInverse(
const Matrix& A, Magnitude tol = Teuchos::ScalarTraits<SC>::eps()*100) {
378 static ArrayRCP<SC> GetLumpedMatrixDiagonal(
const Matrix& A) {
381 static RCP<Vector> GetLumpedMatrixDiagonal(RCP<const Matrix > A) {
384 static RCP<Vector> GetMatrixOverlappedDiagonal(
const Matrix& A) {
387 static RCP<Vector> GetInverse(RCP<const Vector> v, Magnitude tol = Teuchos::ScalarTraits<SC>::eps()*100, SC tolReplacement = Teuchos::ScalarTraits<SC>::zero()) {
390 static Teuchos::Array<Magnitude> ResidualNorm(
const Operator& Op,
const MultiVector& X,
const MultiVector& RHS) {
393 static RCP<MultiVector> Residual(
const Operator& Op,
const MultiVector& X,
const MultiVector& RHS) {
396 static void PauseForDebugger() {
399 static RCP<Teuchos::FancyOStream> MakeFancy(std::ostream& os) {
407 static Kokkos::View<const bool*, typename Node::device_type>
DetectDirichletRows(
const Matrix& A,
const Magnitude& tol = Teuchos::ScalarTraits<SC>::zero(),
const bool count_twos_as_dirichlet=
false);
409 static Kokkos::View<const bool*, typename Node::device_type>
DetectDirichletCols(
const Matrix& A,
const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows);
411 static void ZeroDirichletRows(RCP<Matrix>& A,
const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
413 static void ZeroDirichletRows(RCP<MultiVector>& X,
const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
415 static void ZeroDirichletCols(RCP<Matrix>& A,
const Kokkos::View<const bool*, typename Node::device_type>& dirichletCols, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
417 static RCP<MultiVector> RealValuedToScalarMultiVector(RCP<RealValuedMultiVector> X);
419 static Scalar PowerMethod(
const Matrix& A,
bool scaleByDiag =
true, LO niters = 10, Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
423 static void MyOldScaleMatrix(Matrix& Op,
const Teuchos::ArrayRCP<const SC>& scalingVector,
bool doInverse =
true,
424 bool doFillComplete =
true,
bool doOptimizeStorage =
true) {
425 SC one = Teuchos::ScalarTraits<SC>::one();
426 Teuchos::ArrayRCP<SC> sv(scalingVector.size());
428 for (
int i = 0; i < scalingVector.size(); ++i)
429 sv[i] = one / scalingVector[i];
431 for (
int i = 0; i < scalingVector.size(); ++i)
432 sv[i] = scalingVector[i];
435 switch (Op.getRowMap()->lib()) {
436 case Xpetra::UseTpetra:
437 MyOldScaleMatrix_Tpetra(Op, sv, doFillComplete, doOptimizeStorage);
440 case Xpetra::UseEpetra:
441 MyOldScaleMatrix_Epetra(Op, sv, doFillComplete, doOptimizeStorage);
445 throw Exceptions::RuntimeError(
"Only Epetra and Tpetra matrices can be scaled.");
446 #ifndef __NVCC__ //prevent nvcc warning 453 static void MyOldScaleMatrix_Tpetra(Matrix& Op,
const Teuchos::ArrayRCP<Scalar>& scalingVector,
454 bool doFillComplete,
bool doOptimizeStorage) {
455 #ifdef HAVE_MUELU_TPETRA 456 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT 458 Tpetra::CrsMatrix<SC,LO,GO,NO>& tpOp = Op2NonConstTpetraCrs(Op);
460 const RCP<const Tpetra::Map<LO,GO,NO> > rowMap = tpOp.getRowMap();
461 const RCP<const Tpetra::Map<LO,GO,NO> > domainMap = tpOp.getDomainMap();
462 const RCP<const Tpetra::Map<LO,GO,NO> > rangeMap = tpOp.getRangeMap();
464 size_t maxRowSize = tpOp.getNodeMaxNumRowEntries();
465 if (maxRowSize == Teuchos::as<size_t>(-1))
468 std::vector<SC> scaledVals(maxRowSize);
469 if (tpOp.isFillComplete())
472 if (Op.isLocallyIndexed() ==
true) {
473 Teuchos::ArrayView<const LO> cols;
474 Teuchos::ArrayView<const SC> vals;
476 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
477 tpOp.getLocalRowView(i, cols, vals);
478 size_t nnz = tpOp.getNumEntriesInLocalRow(i);
479 if (nnz > maxRowSize) {
481 scaledVals.resize(maxRowSize);
483 for (
size_t j = 0; j < nnz; ++j)
484 scaledVals[j] = vals[j]*scalingVector[i];
487 Teuchos::ArrayView<const SC> valview(&scaledVals[0], nnz);
488 tpOp.replaceLocalValues(i, cols, valview);
493 Teuchos::ArrayView<const GO> cols;
494 Teuchos::ArrayView<const SC> vals;
496 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
497 GO gid = rowMap->getGlobalElement(i);
498 tpOp.getGlobalRowView(gid, cols, vals);
499 size_t nnz = tpOp.getNumEntriesInGlobalRow(gid);
500 if (nnz > maxRowSize) {
502 scaledVals.resize(maxRowSize);
505 for (
size_t j = 0; j < nnz; ++j)
506 scaledVals[j] = vals[j]*scalingVector[i];
509 Teuchos::ArrayView<const SC> valview(&scaledVals[0], nnz);
510 tpOp.replaceGlobalValues(gid, cols, valview);
515 if (doFillComplete) {
516 if (domainMap == Teuchos::null || rangeMap == Teuchos::null)
517 throw Exceptions::RuntimeError(
"In Utilities::Scaling: cannot fillComplete because the domain and/or range map hasn't been defined");
519 RCP<Teuchos::ParameterList> params = rcp(
new Teuchos::ParameterList());
520 params->set(
"Optimize Storage", doOptimizeStorage);
521 params->set(
"No Nonlocal Changes",
true);
522 Op.fillComplete(Op.getDomainMap(), Op.getRangeMap(), params);
525 throw Exceptions::RuntimeError(
"Only Tpetra::CrsMatrix types can be scaled (Err.1)");
528 throw Exceptions::RuntimeError(
"Matrix scaling is not possible because Tpetra has not been compiled with support for LO=GO=int.");
531 throw Exceptions::RuntimeError(
"Matrix scaling is not possible because Tpetra has not been enabled.");
535 static void MyOldScaleMatrix_Epetra(Matrix& Op,
const Teuchos::ArrayRCP<Scalar>& scalingVector,
bool doFillComplete,
bool doOptimizeStorage) {
536 #ifdef HAVE_MUELU_EPETRA 539 const Epetra_CrsMatrix& epOp = Op2NonConstEpetraCrs(Op);
541 Epetra_Map
const &rowMap = epOp.RowMap();
546 for (
int i = 0; i < rowMap.NumMyElements(); ++i) {
547 epOp.ExtractMyRowView(i, nnz, vals, cols);
548 for (
int j = 0; j < nnz; ++j)
549 vals[j] *= scalingVector[i];
553 throw Exceptions::RuntimeError(
"Only Epetra_CrsMatrix types can be scaled");
556 throw Exceptions::RuntimeError(
"Matrix scaling is not possible because Epetra has not been enabled.");
557 #endif // HAVE_MUELU_EPETRA 565 static RCP<Matrix> Transpose(Matrix& Op,
bool optimizeTranspose =
false,
const std::string & label = std::string(),
const Teuchos::RCP<Teuchos::ParameterList> ¶ms=Teuchos::null) {
566 switch (Op.getRowMap()->lib()) {
567 case Xpetra::UseTpetra:
569 #ifdef HAVE_MUELU_TPETRA 570 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT 575 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A;
576 Tpetra::RowMatrixTransposer<SC,LO,GO,NO> transposer(rcpFromRef(tpetraOp),label);
577 A = transposer.createTranspose(params);
578 RCP<Xpetra::TpetraCrsMatrix<SC,LO,GO,NO> > AA = rcp(
new Xpetra::TpetraCrsMatrix<SC,LO,GO,NO>(A));
579 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
580 RCP<CrsMatrixWrap> AAAA = rcp(
new CrsMatrixWrap(AAA));
584 catch (std::exception& e) {
585 std::cout <<
"threw exception '" << e.what() <<
"'" << std::endl;
586 throw Exceptions::RuntimeError(
"Utilities::Transpose failed, perhaps because matrix is not a Crs matrix");
589 throw Exceptions::RuntimeError(
"Utilities::Transpose: Tpetra is not compiled with LO=GO=int. Add TPETRA_INST_INT_INT:BOOL=ON to your configuration!");
592 throw Exceptions::RuntimeError(
"Utilities::Transpose: Tpetra is not compiled!");
594 #ifndef __NVCC__ //prevent nvcc warning 598 case Xpetra::UseEpetra:
600 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT) 601 Teuchos::TimeMonitor tm(*Teuchos::TimeMonitor::getNewTimer(
"ZZ Entire Transpose"));
604 EpetraExt::RowMatrix_Transpose transposer;
605 Epetra_CrsMatrix * A =
dynamic_cast<Epetra_CrsMatrix*
>(&transposer(epetraOp));
606 transposer.ReleaseTranspose();
608 RCP<Epetra_CrsMatrix> rcpA(A);
609 RCP<Xpetra::EpetraCrsMatrixT<GO,NO> > AA = rcp(
new Xpetra::EpetraCrsMatrixT<GO,NO> (rcpA));
610 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
611 RCP<CrsMatrixWrap> AAAA = rcp(
new CrsMatrixWrap(AAA));
612 AAAA->fillComplete(Op.getRangeMap(), Op.getDomainMap());
616 throw Exceptions::RuntimeError(
"Epetra (Err. 2)");
618 #ifndef __NVCC__ //prevent nvcc warning 623 throw Exceptions::RuntimeError(
"Only Epetra and Tpetra matrices can be transposed.");
624 #ifndef __NVCC__ //prevent nvcc warning 629 #ifndef __NVCC__ //prevent nvcc warning 630 return Teuchos::null;
636 static RCP<Xpetra::MultiVector<double,LO,GO,NO> > ExtractCoordinatesFromParameterList(ParameterList& paramList) {
637 RCP<Xpetra::MultiVector<double,LO,GO,NO> > coordinates = Teuchos::null;
640 if(paramList.isParameter (
"Coordinates") ==
false)
643 #if defined(HAVE_MUELU_TPETRA) 644 #if ( defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_INT)) || \ 645 (!defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT)) 650 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT) 651 typedef Tpetra::MultiVector<float, LO,GO,NO> tfMV;
652 RCP<tfMV> floatCoords = Teuchos::null;
658 typedef Tpetra::MultiVector<SC,LO,GO,NO> tdMV;
659 RCP<tdMV> doubleCoords = Teuchos::null;
660 if (paramList.isType<RCP<tdMV> >(
"Coordinates")) {
662 doubleCoords = paramList.get<RCP<tdMV> >(
"Coordinates");
663 paramList.remove(
"Coordinates");
665 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT) 666 else if (paramList.isType<RCP<tfMV> >(
"Coordinates")) {
668 floatCoords = paramList.get<RCP<tfMV> >(
"Coordinates");
669 paramList.remove(
"Coordinates");
670 doubleCoords = rcp(
new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
671 deep_copy(*doubleCoords, *floatCoords);
675 if(doubleCoords != Teuchos::null) {
676 coordinates = Teuchos::rcp(
new Xpetra::TpetraMultiVector<SC,LO,GO,NO>(doubleCoords));
677 TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
679 #endif // Tpetra instantiated on GO=int and EpetraNode 680 #endif // endif HAVE_TPETRA 682 #if defined(HAVE_MUELU_EPETRA) 683 RCP<Epetra_MultiVector> doubleEpCoords;
684 if (paramList.isType<RCP<Epetra_MultiVector> >(
"Coordinates")) {
685 doubleEpCoords = paramList.get<RCP<Epetra_MultiVector> >(
"Coordinates");
686 paramList.remove(
"Coordinates");
687 RCP<Xpetra::EpetraMultiVectorT<GO,NO> > epCoordinates = Teuchos::rcp(
new Xpetra::EpetraMultiVectorT<GO,NO>(doubleEpCoords));
688 coordinates = rcp_dynamic_cast<Xpetra::MultiVector<double,LO,GO,NO> >(epCoordinates);
689 TEUCHOS_TEST_FOR_EXCEPT(doubleEpCoords->NumVectors() != Teuchos::as<int>(coordinates->getNumVectors()));
692 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(coordinates));
700 template <
class View,
unsigned AppendValue >
705 template <
class MT,
unsigned T >
706 struct CombineMemoryTraits {
710 template <
unsigned U,
unsigned T>
711 struct CombineMemoryTraits<Kokkos::MemoryTraits<U>, T> {
712 typedef Kokkos::MemoryTraits<U|T> type;
715 template <
class DataType,
unsigned T,
class... Pack >
716 struct AppendTrait< Kokkos::View< DataType, Pack... >, T> {
717 typedef Kokkos::View< DataType, Pack... > view_type;
718 using type = Kokkos::View< DataType, typename view_type::array_layout, typename view_type::device_type, typename CombineMemoryTraits<typename view_type::memory_traits,T>::type >;
723 #define MUELU_UTILITIES_KOKKOS_SHORT 725 #endif // #if defined(HAVE_MUELU_KOKKOS_REFACTOR) 727 #endif // MUELU_UTILITIES_KOKKOS_DECL_HPP static void SetRandomSeed(const Teuchos::Comm< int > &comm)
static const RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Map2TpetraMap(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map)
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2TpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
static RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static const Epetra_Map & Map2EpetraMap(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map)
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
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.
static RCP< const Epetra_MultiVector > MV2EpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Epetra objects from an Xpetra object.
void ZeroDirichletCols(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Kokkos::View< const bool *, typename Node::device_type > &dirichletCols, Scalar replaceWith)
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.
void ZeroDirichletRows(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Kokkos::View< const bool *, typename Node::device_type > &dirichletRows, Scalar replaceWith)
Namespace for MueLu classes and methods.
static RCP< Xpetra::MultiVector< double, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList ¶mList)
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
static void PauseForDebugger()
static RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS)
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
static RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV2(Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
static RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraRow(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static Teuchos::ArrayRCP< Scalar > GetLumpedMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal of lumped matrix.
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Crs2Op(RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS)
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)
Kokkos::View< const bool *, typename NO::device_type > DetectDirichletCols(const Xpetra::Matrix< SC, LO, GO, NO > &A, const Kokkos::View< const bool *, typename NO::device_type > &dirichletRows)
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
Extract Overlapped Matrix Diagonal.
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
static Scalar PowerMethod(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool scaleByDiag=true, LocalOrdinal niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
Kokkos::View< const bool *, typename NO::device_type > DetectDirichletRows(const Xpetra::Matrix< SC, LO, GO, NO > &A, const typename Teuchos::ScalarTraits< SC >::magnitudeType &tol, const bool count_twos_as_dirichlet)
static Teuchos::ArrayRCP< Scalar > GetMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal.