46#ifndef MUELU_UTILITIES_KOKKOS_DECL_HPP
47#define MUELU_UTILITIES_KOKKOS_DECL_HPP
53#include "Teuchos_ScalarTraits.hpp"
54#include "Teuchos_ParameterList.hpp"
56#include "Xpetra_BlockedCrsMatrix_fwd.hpp"
57#include "Xpetra_CrsMatrix_fwd.hpp"
58#include "Xpetra_CrsMatrixWrap_fwd.hpp"
59#include "Xpetra_ExportFactory.hpp"
60#include "Xpetra_ImportFactory_fwd.hpp"
61#include "Xpetra_MapFactory_fwd.hpp"
62#include "Xpetra_Map_fwd.hpp"
63#include "Xpetra_MatrixFactory_fwd.hpp"
64#include "Xpetra_Matrix_fwd.hpp"
65#include "Xpetra_MultiVectorFactory_fwd.hpp"
66#include "Xpetra_MultiVector_fwd.hpp"
67#include "Xpetra_Operator_fwd.hpp"
68#include "Xpetra_VectorFactory_fwd.hpp"
69#include "Xpetra_Vector_fwd.hpp"
71#include "Xpetra_IO.hpp"
73#include "Kokkos_ArithTraits.hpp"
75#ifdef HAVE_MUELU_EPETRA
76#include "Epetra_MultiVector.h"
77#include "Epetra_CrsMatrix.h"
78#include "Xpetra_EpetraCrsMatrix_fwd.hpp"
79#include "Xpetra_EpetraMultiVector_fwd.hpp"
83#include "MueLu_Utilities.hpp"
84#include "MueLu_UtilitiesBase.hpp"
86#ifdef HAVE_MUELU_TPETRA
87#include "Tpetra_CrsMatrix.hpp"
88#include "Tpetra_Map.hpp"
89#include "Tpetra_MultiVector.hpp"
90#include "Xpetra_TpetraCrsMatrix_fwd.hpp"
91#include "Xpetra_TpetraMultiVector_fwd.hpp"
109#undef MUELU_UTILITIES_KOKKOS_SHORT
113 using TST = Teuchos::ScalarTraits<SC>;
118#ifdef HAVE_MUELU_EPETRA
137#ifdef HAVE_MUELU_TPETRA
225 LO niters = 10,
Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
230 LO niters = 10,
Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
234 static void MyOldScaleMatrix(
Matrix& Op,
const Teuchos::ArrayRCP<const SC>& scalingVector,
bool doInverse =
true,
235 bool doFillComplete =
true,
bool doOptimizeStorage =
true);
238 bool doFillComplete,
bool doOptimizeStorage);
241 bool doFillComplete,
bool doOptimizeStorage);
252 static Kokkos::View<bool*, typename NO::device_type>
DetectDirichletRows(
const Matrix& A,
const Magnitude& tol = Teuchos::ScalarTraits<
typename Teuchos::ScalarTraits<SC>::magnitudeType>::zero(),
const bool count_twos_as_dirichlet=
false);
261 static void FindNonZeros(
const typename Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dual_view_type::t_dev_const_um vals,
262 Kokkos::View<bool*, typename Node::device_type> nonzeros);
275 static Kokkos::View<bool*, typename NO::device_type>
DetectDirichletCols(
const Matrix& A,
const Kokkos::View<const bool*, typename NO::device_type>& dirichletRows);
286 const Kokkos::View<bool*, typename Node::device_type> & dirichletRows,
287 Kokkos::View<bool*, typename Node::device_type> dirichletCols,
288 Kokkos::View<bool*, typename Node::device_type> dirichletDomain);
291 static void ZeroDirichletRows(RCP<Matrix>& A,
const Kokkos::View<const bool*, typename NO::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
293 static void ZeroDirichletRows(RCP<MultiVector>& X,
const Kokkos::View<const bool*, typename NO::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
295 static void ZeroDirichletCols(RCP<Matrix>& A,
const Kokkos::View<const bool*, typename NO::device_type>& dirichletCols, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
298 const typename Teuchos::ScalarTraits<Scalar>::magnitudeType rowSumTol,
299 Kokkos::View<bool*, typename NO::device_type> & dirichletRows);
319 static RCP<Matrix>
Transpose(
Matrix& Op,
bool optimizeTranspose =
false,
const std::string & label = std::string()) {
330 static RCP<Xpetra::Vector<LocalOrdinal,LocalOrdinal,GlobalOrdinal,Node> >
CuthillMcKee(
const Matrix &Op);
336 static void ApplyOAZToMatrixRows(RCP<Matrix>& A,
const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows);
350 template <
class Node>
356 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType
Magnitude;
361#undef MUELU_UTILITIES_KOKKOS_SHORT
366#ifdef HAVE_MUELU_EPETRA
385#ifdef HAVE_MUELU_TPETRA
409 const auto rowMap = A.getRowMap();
410 auto diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,
true);
412 A.getLocalDiagCopy(*diag);
418 static RCP<Vector>
GetLumpedMatrixDiagonal(
Matrix const &A,
const bool doReciprocal=
false,
Magnitude tol = Teuchos::ScalarTraits<Scalar>::eps()*100,
Scalar tolReplacement = Teuchos::ScalarTraits<Scalar>::zero(),
const bool replaceSingleEntryRowWithZero =
false,
const bool useAverageAbsDiagVal =
false) {
424 static RCP<Vector>
GetInverse(RCP<const Vector> v,
Magnitude tol = Teuchos::ScalarTraits<SC>::eps()*100, SC tolReplacement = Teuchos::ScalarTraits<SC>::zero()) {
433 static RCP<Teuchos::FancyOStream>
MakeFancy(std::ostream& os) {
441 static Kokkos::View<bool*, typename Node::device_type>
DetectDirichletRows(
const Matrix& A,
const Magnitude& tol = Teuchos::ScalarTraits<SC>::zero(),
const bool count_twos_as_dirichlet=
false);
443 static Kokkos::View<bool*, typename Node::device_type>
DetectDirichletCols(
const Matrix& A,
const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows);
445 static void FindNonZeros(
const typename Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dual_view_type::t_dev_const_um vals,
446 Kokkos::View<bool*, typename Node::device_type> nonzeros);
449 const Kokkos::View<bool*, typename Node::device_type> & dirichletRows,
450 Kokkos::View<bool*, typename Node::device_type> dirichletCols,
451 Kokkos::View<bool*, typename Node::device_type> dirichletDomain);
453 static void ZeroDirichletRows(RCP<Matrix>& A,
const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
455 static void ZeroDirichletRows(RCP<MultiVector>& X,
const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
457 static void ZeroDirichletCols(RCP<Matrix>& A,
const Kokkos::View<const bool*, typename Node::device_type>& dirichletCols, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
460 const typename Teuchos::ScalarTraits<Scalar>::magnitudeType rowSumTol,
461 Kokkos::View<bool*, typename NO::device_type> & dirichletRows);
469 static Scalar PowerMethod(
const Matrix& A,
const Teuchos::RCP<Vector> &invDiag, LO niters = 10,
Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
474 bool doFillComplete =
true,
bool doOptimizeStorage =
true) {
475 SC one = Teuchos::ScalarTraits<SC>::one();
476 Teuchos::ArrayRCP<SC> sv(scalingVector.size());
478 for (
int i = 0; i < scalingVector.size(); ++i)
479 sv[i] = one / scalingVector[i];
481 for (
int i = 0; i < scalingVector.size(); ++i)
482 sv[i] = scalingVector[i];
485 switch (Op.getRowMap()->lib()) {
486 case Xpetra::UseTpetra:
490 case Xpetra::UseEpetra:
504 bool doFillComplete,
bool doOptimizeStorage) {
505 #ifdef HAVE_MUELU_TPETRA
506 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT
510 const RCP<const Tpetra::Map<LO,GO,NO> > rowMap = tpOp.getRowMap();
511 const RCP<const Tpetra::Map<LO,GO,NO> > domainMap = tpOp.getDomainMap();
512 const RCP<const Tpetra::Map<LO,GO,NO> > rangeMap = tpOp.getRangeMap();
514 size_t maxRowSize = tpOp.getLocalMaxNumRowEntries();
515 if (maxRowSize == Teuchos::as<size_t>(-1))
518 std::vector<SC> scaledVals(maxRowSize);
519 if (tpOp.isFillComplete())
522 if (Op.isLocallyIndexed() ==
true) {
523 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type cols;
524 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type vals;
526 for (
size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
527 tpOp.getLocalRowView(i, cols, vals);
528 size_t nnz = tpOp.getNumEntriesInLocalRow(i);
529 if (nnz > maxRowSize) {
531 scaledVals.resize(maxRowSize);
533 for (
size_t j = 0; j < nnz; ++j)
534 scaledVals[j] = vals[j]*scalingVector[i];
537 Teuchos::ArrayView<const LocalOrdinal> cols_view(cols.data(), nnz);
538 Teuchos::ArrayView<const SC> valview(&scaledVals[0], nnz);
539 tpOp.replaceLocalValues(i, cols_view, valview);
544 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::global_inds_host_view_type cols;
545 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type vals;
547 for (
size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
548 GO gid = rowMap->getGlobalElement(i);
549 tpOp.getGlobalRowView(gid, cols, vals);
550 size_t nnz = tpOp.getNumEntriesInGlobalRow(gid);
551 if (nnz > maxRowSize) {
553 scaledVals.resize(maxRowSize);
556 for (
size_t j = 0; j < nnz; ++j)
557 scaledVals[j] = vals[j]*scalingVector[i];
560 Teuchos::ArrayView<const GlobalOrdinal> cols_view(cols.data(), nnz);
561 Teuchos::ArrayView<const SC> valview(&scaledVals[0], nnz);
562 tpOp.replaceGlobalValues(gid, cols_view, valview);
567 if (doFillComplete) {
568 if (domainMap == Teuchos::null || rangeMap == Teuchos::null)
569 throw Exceptions::RuntimeError(
"In Utilities::Scaling: cannot fillComplete because the domain and/or range map hasn't been defined");
571 RCP<Teuchos::ParameterList> params = rcp(
new Teuchos::ParameterList());
572 params->set(
"Optimize Storage", doOptimizeStorage);
573 params->set(
"No Nonlocal Changes",
true);
574 Op.fillComplete(Op.getDomainMap(), Op.getRangeMap(), params);
580 throw Exceptions::RuntimeError(
"Matrix scaling is not possible because Tpetra has not been compiled with support for LO=GO=int.");
588#ifdef HAVE_MUELU_EPETRA
600 for (
int j = 0; j < nnz; ++j)
601 vals[j] *= scalingVector[i];
617 static RCP<Matrix>
Transpose(
Matrix& Op,
bool =
false,
const std::string & label = std::string(),
const Teuchos::RCP<Teuchos::ParameterList> ¶ms=Teuchos::null) {
618 switch (Op.getRowMap()->lib()) {
619 case Xpetra::UseTpetra:
621#ifdef HAVE_MUELU_TPETRA
622#ifdef HAVE_MUELU_TPETRA_INST_INT_INT
627 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A;
628 Tpetra::RowMatrixTransposer<SC,LO,GO,NO> transposer(rcpFromRef(tpetraOp),label);
630 using Teuchos::ParameterList;
632 RCP<ParameterList> transposeParams = params.is_null () ?
633 rcp (
new ParameterList) :
634 rcp (
new ParameterList (*params));
635 transposeParams->set (
"sort",
false);
636 A = transposer.createTranspose(transposeParams);
639 RCP<Xpetra::TpetraCrsMatrix<SC,LO,GO,NO> > AA = rcp(
new Xpetra::TpetraCrsMatrix<SC,LO,GO,NO>(A));
640 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
645 catch (std::exception& e) {
646 std::cout <<
"threw exception '" << e.what() <<
"'" << std::endl;
650 throw Exceptions::RuntimeError(
"Utilities::Transpose: Tpetra is not compiled with LO=GO=int. Add TPETRA_INST_INT_INT:BOOL=ON to your configuration!");
659 case Xpetra::UseEpetra:
661#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
662 Teuchos::TimeMonitor tm(*Teuchos::TimeMonitor::getNewTimer(
"ZZ Entire Transpose"));
669 RCP<Epetra_CrsMatrix> rcpA(A);
670 RCP<Xpetra::EpetraCrsMatrixT<GO,NO> > AA = rcp(
new Xpetra::EpetraCrsMatrixT<GO,NO> (rcpA));
671 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
673 AAAA->fillComplete(Op.getRangeMap(), Op.getDomainMap());
691 return Teuchos::null;
698 RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > coordinates = Teuchos::null;
701 if(paramList.isParameter (
"Coordinates") ==
false)
704#if defined(HAVE_MUELU_TPETRA)
705#if ( defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_INT)) || \
706 (!defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT))
711#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
712 typedef Tpetra::MultiVector<float, LO,GO,NO> tfMV;
713 RCP<tfMV> floatCoords = Teuchos::null;
719 typedef Tpetra::MultiVector<SC,LO,GO,NO> tdMV;
720 RCP<tdMV> doubleCoords = Teuchos::null;
721 if (paramList.isType<RCP<tdMV> >(
"Coordinates")) {
723 doubleCoords = paramList.get<RCP<tdMV> >(
"Coordinates");
724 paramList.remove(
"Coordinates");
726#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
727 else if (paramList.isType<RCP<tfMV> >(
"Coordinates")) {
729 floatCoords = paramList.get<RCP<tfMV> >(
"Coordinates");
730 paramList.remove(
"Coordinates");
731 doubleCoords = rcp(
new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
732 deep_copy(*doubleCoords, *floatCoords);
736 if(doubleCoords != Teuchos::null) {
737 coordinates = Teuchos::rcp(
new Xpetra::TpetraMultiVector<SC,LO,GO,NO>(doubleCoords));
738 TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
743#if defined(HAVE_MUELU_EPETRA)
744 RCP<Epetra_MultiVector> doubleEpCoords;
745 if (paramList.isType<RCP<Epetra_MultiVector> >(
"Coordinates")) {
746 doubleEpCoords = paramList.get<RCP<Epetra_MultiVector> >(
"Coordinates");
747 paramList.remove(
"Coordinates");
748 RCP<Xpetra::EpetraMultiVectorT<GO,NO> > epCoordinates = Teuchos::rcp(
new Xpetra::EpetraMultiVectorT<GO,NO>(doubleEpCoords));
749 coordinates = rcp_dynamic_cast<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> >(epCoordinates);
750 TEUCHOS_TEST_FOR_EXCEPT(doubleEpCoords->NumVectors() != Teuchos::as<int>(coordinates->getNumVectors()));
753 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(coordinates));
765 static RCP<Xpetra::Vector<LocalOrdinal,LocalOrdinal,GlobalOrdinal,Node> >
CuthillMcKee(
const Matrix &Op);
767 static void ApplyOAZToMatrixRows(RCP<Matrix>& A,
const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows);
774 template <
class View,
unsigned AppendValue >
779 template <
class MT,
unsigned T >
784 template <
unsigned U,
unsigned T>
786 typedef Kokkos::MemoryTraits<U|T>
type;
789 template <
class DataType,
unsigned T,
class... Pack >
792 using type = Kokkos::View< DataType, typename view_type::array_layout, typename view_type::device_type, typename CombineMemoryTraits<typename view_type::memory_traits,T>::type >;
797#define MUELU_UTILITIES_KOKKOS_SHORT
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
int NumMyElements() const
const Epetra_Map & RowMap() const
int ExtractMyRowView(int MyRow, int &NumEntries, double *&Values, int *&Indices) const
Exception throws to report errors in the internal logical of the program.
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MultiVector
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
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.
static Teuchos::RCP< Vector > GetLumpedMatrixDiagonal(Matrix const &A, const bool doReciprocal=false, Magnitude tol=Teuchos::ScalarTraits< Scalar >::magnitude(Teuchos::ScalarTraits< Scalar >::zero()), Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool replaceSingleEntryRowWithZero=false, const bool useAverageAbsDiagVal=false)
Extract Matrix Diagonal of lumped matrix.
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixWrap
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
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.
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
static RCP< const Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2TpetraCrs(RCP< const Matrix > Op)
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LO, GO, NO > > ExtractCoordinatesFromParameterList(ParameterList ¶mList)
Extract coordinates from parameter list and return them in a Xpetra::MultiVector.
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Matrix > Op)
static RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2NonConstTpetraCrs(RCP< Matrix > Op)
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< MultiVector > vec)
static Teuchos::Array< Magnitude > ResidualNorm(const Operator &Op, const MultiVector &X, const MultiVector &RHS)
static void MyOldScaleMatrix(Matrix &Op, const Teuchos::ArrayRCP< const SC > &scalingVector, bool doInverse=true, bool doFillComplete=true, bool doOptimizeStorage=true)
static const Epetra_CrsMatrix & Op2EpetraCrs(const Matrix &Op)
static Epetra_CrsMatrix & Op2NonConstEpetraCrs(Matrix &Op)
static RCP< Vector > GetInverse(RCP< const Vector > v, Magnitude tol=Teuchos::ScalarTraits< SC >::eps() *100, SC tolReplacement=Teuchos::ScalarTraits< SC >::zero())
static void ZeroDirichletRows(RCP< MultiVector > &X, const Kokkos::View< const bool *, typename Node::device_type > &dirichletRows, SC replaceWith=Teuchos::ScalarTraits< SC >::zero())
static const Tpetra::MultiVector< SC, LO, GO, NO > & MV2TpetraMV(const MultiVector &vec)
static const Epetra_MultiVector & MV2EpetraMV(const MultiVector &vec)
static RCP< Xpetra::Matrix< SC, LO, GO, NO > > Crs2Op(RCP< CrsMatrix > Op)
static RCP< Tpetra::MultiVector< SC, LO, GO, NO > > MV2NonConstTpetraMV2(MultiVector &vec)
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Matrix > Op)
Xpetra::MultiVector< CoordinateType, LocalOrdinal, GlobalOrdinal, Node > RealValuedMultiVector
static void ZeroDirichletRows(RCP< Matrix > &A, const Kokkos::View< const bool *, typename Node::device_type > &dirichletRows, SC replaceWith=Teuchos::ScalarTraits< SC >::zero())
static void MyOldScaleMatrix_Epetra(Matrix &Op, const Teuchos::ArrayRCP< Scalar > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
static Epetra_MultiVector & MV2NonConstEpetraMV(MultiVector &vec)
static void MyOldScaleMatrix_Tpetra(Matrix &Op, const Teuchos::ArrayRCP< Scalar > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
static RCP< const Tpetra::MultiVector< SC, LO, GO, NO > > MV2TpetraMV(RCP< MultiVector > const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
static RCP< Tpetra::MultiVector< SC, LO, GO, NO > > MV2NonConstTpetraMV(RCP< MultiVector > vec)
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
static RCP< Vector > GetMatrixDiagonal(const Matrix &A)
static RCP< Tpetra::RowMatrix< SC, LO, GO, NO > > Op2NonConstTpetraRow(RCP< Matrix > Op)
static RCP< Vector > GetLumpedMatrixDiagonal(Matrix const &A, const bool doReciprocal=false, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar tolReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool replaceSingleEntryRowWithZero=false, const bool useAverageAbsDiagVal=false)
static Scalar PowerMethod(const Matrix &A, const Teuchos::RCP< Vector > &invDiag, LO niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
static RCP< const Epetra_MultiVector > MV2EpetraMV(RCP< MultiVector > const vec)
Helper utility to pull out the underlying Epetra objects from an Xpetra object.
static RCP< Matrix > Transpose(Matrix &Op, bool=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Transpose a Xpetra::Matrix.
static RCP< MultiVector > Residual(const Operator &Op, const MultiVector &X, const MultiVector &RHS)
static Tpetra::MultiVector< SC, LO, GO, NO > & MV2NonConstTpetraMV(MultiVector &vec)
static Scalar PowerMethod(const Matrix &A, bool scaleByDiag=true, LO niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
static Tpetra::CrsMatrix< SC, LO, GO, NO > & Op2NonConstTpetraCrs(Matrix &Op)
typename Teuchos::ScalarTraits< Scalar >::coordinateType CoordinateType
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
static const RCP< const Tpetra::Map< LO, GO, NO > > Map2TpetraMap(const Map &map)
static RCP< const Tpetra::RowMatrix< SC, LO, GO, NO > > Op2TpetraRow(RCP< const Matrix > Op)
static const Tpetra::CrsMatrix< SC, LO, GO, NO > & Op2TpetraCrs(const Matrix &Op)
static const Epetra_Map & Map2EpetraMap(const Map &map)
static const Tpetra::CrsMatrix< SC, LO, GO, NO > & Op2TpetraCrs(const Matrix &Op)
static void ApplyOAZToMatrixRows(RCP< Matrix > &A, const Kokkos::View< const bool *, typename Node::device_type > &dirichletRows)
static RCP< Tpetra::RowMatrix< SC, LO, GO, NO > > Op2NonConstTpetraRow(RCP< Matrix > Op)
static RCP< const Tpetra::RowMatrix< SC, LO, GO, NO > > Op2TpetraRow(RCP< const Matrix > Op)
static Teuchos::Array< Magnitude > ResidualNorm(const Operator &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< const Tpetra::MultiVector< SC, LO, GO, NO > > MV2TpetraMV(RCP< MultiVector > const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
static Epetra_CrsMatrix & Op2NonConstEpetraCrs(Matrix &Op)
static void MyOldScaleMatrix_Tpetra(Matrix &Op, const Teuchos::ArrayRCP< SC > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Matrix > Op)
static RCP< Vector > GetMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal.
Teuchos::ScalarTraits< SC > TST
Xpetra::MultiVector< CoordinateType, LO, GO, NO > RealValuedMultiVector
static const Tpetra::MultiVector< SC, LO, GO, NO > & MV2TpetraMV(const MultiVector &vec)
static RCP< const Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2TpetraCrs(RCP< const Matrix > Op)
static RCP< Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > > ReverseCuthillMcKee(const Matrix &Op)
static const Epetra_CrsMatrix & Op2EpetraCrs(const Matrix &Op)
static SC PowerMethod(const Matrix &A, bool scaleByDiag=true, LO niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
Simple transpose for Tpetra::CrsMatrix types.
static Tpetra::MultiVector< SC, LO, GO, NO > & MV2NonConstTpetraMV(MultiVector &vec)
static RCP< MultiVector > RealValuedToScalarMultiVector(RCP< RealValuedMultiVector > X)
static RCP< Vector > GetMatrixDiagonalInverse(const Matrix &A, Magnitude tol=TST::eps() *100, const bool doLumped=false)
Extract Matrix Diagonal.
static Epetra_MultiVector & MV2NonConstEpetraMV(MultiVector &vec)
static void ZeroDirichletRows(RCP< Matrix > &A, const Kokkos::View< const bool *, typename NO::device_type > &dirichletRows, SC replaceWith=Teuchos::ScalarTraits< SC >::zero())
static void MyOldScaleMatrix_Epetra(Matrix &Op, const Teuchos::ArrayRCP< SC > &scalingVector, bool doFillComplete, bool doOptimizeStorage)
static RCP< MultiVector > Residual(const Operator &Op, const MultiVector &X, const MultiVector &RHS)
static void DetectDirichletColsAndDomains(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Kokkos::View< bool *, typename Node::device_type > &dirichletRows, Kokkos::View< bool *, typename Node::device_type > dirichletCols, Kokkos::View< bool *, typename Node::device_type > dirichletDomain)
Detects Dirichlet columns & domains from a list of Dirichlet rows.
static RCP< const Epetra_MultiVector > MV2EpetraMV(RCP< MultiVector > const vec)
Helper utility to pull out the underlying Epetra objects from an Xpetra object.
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
static RCP< Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > > CuthillMcKee(const Matrix &Op)
static SC PowerMethod(const Matrix &A, const Teuchos::RCP< Vector > &invDiag, LO niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
static RCP< Matrix > Transpose(Matrix &Op, bool optimizeTranspose=false, const std::string &label=std::string())
Transpose a Xpetra::Matrix.
static const RCP< const Tpetra::Map< LO, GO, NO > > Map2TpetraMap(const Map &map)
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 Kokkos::View< bool *, typename NO::device_type > DetectDirichletCols(const Matrix &A, const Kokkos::View< const bool *, typename NO::device_type > &dirichletRows)
Detect Dirichlet columns based on Dirichlet rows.
static RCP< Tpetra::MultiVector< SC, LO, GO, NO > > MV2NonConstTpetraMV2(MultiVector &vec)
static void ZeroDirichletRows(RCP< MultiVector > &X, const Kokkos::View< const bool *, typename NO::device_type > &dirichletRows, SC replaceWith=Teuchos::ScalarTraits< SC >::zero())
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< MultiVector > vec)
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
Extract Overlapped Matrix Diagonal.
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Matrix > Op)
static void FindNonZeros(const typename Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dual_view_type::t_dev_const_um vals, Kokkos::View< bool *, typename Node::device_type > nonzeros)
Find non-zero values in an ArrayRCP Compares the value to 2 * machine epsilon.
static const Epetra_Map & Map2EpetraMap(const Map &map)
static void ApplyRowSumCriterion(const Matrix &A, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType rowSumTol, Kokkos::View< bool *, typename NO::device_type > &dirichletRows)
static void MyOldScaleMatrix(Matrix &Op, const Teuchos::ArrayRCP< const SC > &scalingVector, bool doInverse=true, bool doFillComplete=true, bool doOptimizeStorage=true)
static Kokkos::View< bool *, typename NO::device_type > DetectDirichletRows(const Matrix &A, const Magnitude &tol=Teuchos::ScalarTraits< typename Teuchos::ScalarTraits< SC >::magnitudeType >::zero(), const bool count_twos_as_dirichlet=false)
Detect Dirichlet rows.
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList ¶mList)
static RCP< Xpetra::Matrix< SC, LO, GO, NO > > Crs2Op(RCP< CrsMatrix > Op)
static Tpetra::CrsMatrix< SC, LO, GO, NO > & Op2NonConstTpetraCrs(Matrix &Op)
typename TST::magnitudeType Magnitude
static void ZeroDirichletCols(RCP< Matrix > &A, const Kokkos::View< const bool *, typename NO::device_type > &dirichletCols, SC replaceWith=Teuchos::ScalarTraits< SC >::zero())
static RCP< Tpetra::MultiVector< SC, LO, GO, NO > > MV2NonConstTpetraMV(RCP< MultiVector > vec)
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
static const Epetra_MultiVector & MV2EpetraMV(const MultiVector &vec)
typename TST::coordinateType CoordinateType
static RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2NonConstTpetraCrs(RCP< Matrix > Op)
static RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
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 Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV2(Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec)
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 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 const Epetra_Map & Map2EpetraMap(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map)
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< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
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_CrsMatrix > Op2NonConstEpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraRow(RCP< const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList ¶mList)
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.
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Crs2Op(RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
Namespace for MueLu classes and methods.
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Kokkos::View< DataType, typename view_type::array_layout, typename view_type::device_type, typename CombineMemoryTraits< typename view_type::memory_traits, T >::type > type
Kokkos::View< DataType, Pack... > view_type
Kokkos::MemoryTraits< U|T > type