46 #ifndef MUELU_UTILITIES_KOKKOS_DECL_HPP 47 #define MUELU_UTILITIES_KOKKOS_DECL_HPP 50 #if defined(HAVE_MUELU_KOKKOS_REFACTOR) 74 #ifdef HAVE_MUELU_EPETRA 75 #include <Epetra_MultiVector.h> 76 #include <Epetra_CrsMatrix.h> 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> 103 template <
class Scalar,
104 class LocalOrdinal = int,
105 class GlobalOrdinal = LocalOrdinal,
108 #undef MUELU_UTILITIES_KOKKOS_SHORT 114 #ifdef HAVE_MUELU_EPETRA 117 static RCP<const Epetra_MultiVector> MV2EpetraMV(RCP<MultiVector>
const vec) {
return Utilities::MV2EpetraMV(vec); }
133 #ifdef HAVE_MUELU_TPETRA 136 static RCP<const Tpetra::MultiVector<SC,LO,GO,NO> > MV2TpetraMV(RCP<MultiVector>
const vec) {
return Utilities::MV2TpetraMV(vec); }
140 static const Tpetra::MultiVector<SC,LO,GO,NO>& MV2TpetraMV(
const MultiVector& vec) {
return Utilities::MV2TpetraMV(vec); }
143 static RCP<const Tpetra::CrsMatrix<SC,LO,GO,NO> > Op2TpetraCrs(RCP<const Matrix> Op) {
return Utilities::Op2TpetraCrs(Op); }
146 static const Tpetra::CrsMatrix<SC,LO,GO,NO>& Op2TpetraCrs(
const Matrix& Op) {
return Utilities::Op2TpetraCrs(Op); }
149 static RCP<const Tpetra::RowMatrix<SC,LO,GO,NO> > Op2TpetraRow(RCP<const Matrix> Op) {
return Utilities::Op2TpetraRow(Op); }
152 static const RCP<const Tpetra::Map<LO, GO, NO> > Map2TpetraMap(
const Map& map) {
return Utilities::Map2TpetraMap(map); }
155 static RCP<Xpetra::Matrix<SC,LO,GO,NO> > Crs2Op(RCP<CrsMatrix> Op) {
return Utilities::Crs2Op(Op); }
192 static RCP<Vector> GetMatrixOverlappedDiagonal(
const Matrix& A);
211 static RCP<MultiVector> Residual(
const Operator& Op,
const MultiVector& X,
const MultiVector& RHS) {
215 static void PauseForDebugger();
232 static SC PowerMethod(
const Matrix& A,
bool scaleByDiag =
true,
233 LO niters = 10, Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
238 bool doFillComplete =
true,
bool doOptimizeStorage =
true);
241 bool doFillComplete,
bool doOptimizeStorage);
244 bool doFillComplete,
bool doOptimizeStorage);
246 static RCP<Teuchos::FancyOStream> MakeFancy(std::ostream& os) {
return Utilities::MakeFancy(os); }
279 static RCP<Matrix> Transpose(Matrix& Op,
bool optimizeTranspose =
false,
const std::string & label = std::string()) {
283 static RCP<Xpetra::MultiVector<double,LocalOrdinal,GlobalOrdinal,Node> > ExtractCoordinatesFromParameterList(ParameterList& paramList) {
299 template <
class Node>
300 class Utilities_kokkos<double,int,int,Node> :
public UtilitiesBase<double,int,int,Node> {
302 typedef double Scalar;
303 typedef int LocalOrdinal;
304 typedef int GlobalOrdinal;
308 #undef MUELU_UTILITIES_KOKKOS_SHORT 313 #ifdef HAVE_MUELU_EPETRA 316 static RCP<const Epetra_MultiVector> MV2EpetraMV(RCP<MultiVector>
const vec) {
return Utilities::MV2EpetraMV(vec); }
332 #ifdef HAVE_MUELU_TPETRA 335 static RCP<const Tpetra::MultiVector<SC,LO,GO,NO> > MV2TpetraMV(RCP<MultiVector>
const vec) {
return Utilities::MV2TpetraMV(vec); }
339 static const Tpetra::MultiVector<SC,LO,GO,NO>& MV2TpetraMV(
const MultiVector& vec) {
return Utilities::MV2TpetraMV(vec); }
342 static RCP<const Tpetra::CrsMatrix<SC,LO,GO,NO> > Op2TpetraCrs(RCP<const Matrix> Op) {
return Utilities::Op2TpetraCrs(Op); }
345 static const Tpetra::CrsMatrix<SC,LO,GO,NO>& Op2TpetraCrs(
const Matrix& Op) {
return Utilities::Op2TpetraCrs(Op); }
348 static RCP<const Tpetra::RowMatrix<SC,LO,GO,NO> > Op2TpetraRow(RCP<const Matrix> Op) {
return Utilities::Op2TpetraRow(Op); }
351 static const RCP<const Tpetra::Map<LO, GO, NO> > Map2TpetraMap(
const Map& map) {
return Utilities::Map2TpetraMap(map); }
353 static RCP<Xpetra::Matrix<SC,LO,GO,NO> > Crs2Op(RCP<CrsMatrix> Op) {
return Utilities::Crs2Op(Op); }
355 static ArrayRCP<SC> GetMatrixDiagonal(
const Matrix& A) {
361 static ArrayRCP<SC> GetLumpedMatrixDiagonal(
const Matrix& A) {
364 static RCP<Vector> GetLumpedMatrixDiagonal(RCP<const Matrix > A) {
367 static RCP<Vector> GetMatrixOverlappedDiagonal(
const Matrix& A) {
376 static RCP<MultiVector> Residual(
const Operator& Op,
const MultiVector& X,
const MultiVector& RHS) {
379 static void PauseForDebugger() {
382 static RCP<Teuchos::FancyOStream> MakeFancy(std::ostream& os) {
395 static Scalar PowerMethod(
const Matrix& A,
bool scaleByDiag =
true,
LO niters = 10, Magnitude tolerance = 1e-2,
bool verbose =
false,
unsigned int seed = 123) {
400 bool doFillComplete =
true,
bool doOptimizeStorage =
true) {
404 for (
int i = 0; i < scalingVector.
size(); ++i)
405 sv[i] = one / scalingVector[i];
407 for (
int i = 0; i < scalingVector.
size(); ++i)
408 sv[i] = scalingVector[i];
411 switch (Op.getRowMap()->lib()) {
413 MyOldScaleMatrix_Tpetra(Op, sv, doFillComplete, doOptimizeStorage);
417 MyOldScaleMatrix_Epetra(Op, sv, doFillComplete, doOptimizeStorage);
421 throw Exceptions::RuntimeError(
"Only Epetra and Tpetra matrices can be scaled.");
428 bool doFillComplete,
bool doOptimizeStorage) {
429 #ifdef HAVE_MUELU_TPETRA 430 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT 432 Tpetra::CrsMatrix<SC,LO,GO,NO>& tpOp = Op2NonConstTpetraCrs(Op);
434 const RCP<const Tpetra::Map<LO,GO,NO> > rowMap = tpOp.getRowMap();
435 const RCP<const Tpetra::Map<LO,GO,NO> > domainMap = tpOp.getDomainMap();
436 const RCP<const Tpetra::Map<LO,GO,NO> > rangeMap = tpOp.getRangeMap();
438 size_t maxRowSize = tpOp.getNodeMaxNumRowEntries();
439 if (maxRowSize == Teuchos::as<size_t>(-1))
442 std::vector<SC> scaledVals(maxRowSize);
443 if (tpOp.isFillComplete())
446 if (Op.isLocallyIndexed() ==
true) {
450 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
451 tpOp.getLocalRowView(i, cols, vals);
452 size_t nnz = tpOp.getNumEntriesInLocalRow(i);
453 if (nnz > maxRowSize) {
455 scaledVals.resize(maxRowSize);
457 for (
size_t j = 0; j < nnz; ++j)
458 scaledVals[j] = vals[j]*scalingVector[i];
462 tpOp.replaceLocalValues(i, cols, valview);
470 for (
size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
471 GO gid = rowMap->getGlobalElement(i);
472 tpOp.getGlobalRowView(gid, cols, vals);
473 size_t nnz = tpOp.getNumEntriesInGlobalRow(gid);
474 if (nnz > maxRowSize) {
476 scaledVals.resize(maxRowSize);
479 for (
size_t j = 0; j < nnz; ++j)
480 scaledVals[j] = vals[j]*scalingVector[i];
484 tpOp.replaceGlobalValues(gid, cols, valview);
489 if (doFillComplete) {
490 if (domainMap == Teuchos::null || rangeMap == Teuchos::null)
491 throw Exceptions::RuntimeError(
"In Utilities::Scaling: cannot fillComplete because the domain and/or range map hasn't been defined");
494 params->
set(
"Optimize Storage", doOptimizeStorage);
495 params->
set(
"No Nonlocal Changes",
true);
496 Op.fillComplete(Op.getDomainMap(), Op.getRangeMap(), params);
499 throw Exceptions::RuntimeError(
"Only Tpetra::CrsMatrix types can be scaled (Err.1)");
502 throw Exceptions::RuntimeError(
"Matrix scaling is not possible because Tpetra has not been compiled with support for LO=GO=int.");
505 throw Exceptions::RuntimeError(
"Matrix scaling is not possible because Tpetra has not been enabled.");
509 static void MyOldScaleMatrix_Epetra(Matrix& Op,
const Teuchos::ArrayRCP<Scalar>& scalingVector,
bool doFillComplete,
bool doOptimizeStorage) {
510 #ifdef HAVE_MUELU_EPETRA 522 for (
int j = 0; j < nnz; ++j)
523 vals[j] *= scalingVector[i];
527 throw Exceptions::RuntimeError(
"Only Epetra_CrsMatrix types can be scaled");
530 throw Exceptions::RuntimeError(
"Matrix scaling is not possible because Epetra has not been enabled.");
531 #endif // HAVE_MUELU_EPETRA 539 static RCP<Matrix> Transpose(Matrix& Op,
bool optimizeTranspose =
false,
const std::string & label = std::string()) {
540 switch (Op.getRowMap()->lib()) {
543 #ifdef HAVE_MUELU_TPETRA 544 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT 549 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A;
550 Tpetra::RowMatrixTransposer<SC,LO,GO,NO> transposer(rcpFromRef(tpetraOp),label);
551 A = transposer.createTranspose();
553 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
554 RCP<CrsMatrixWrap> AAAA =
rcp(
new CrsMatrixWrap(AAA));
558 catch (std::exception& e) {
559 std::cout <<
"threw exception '" << e.what() <<
"'" << std::endl;
560 throw Exceptions::RuntimeError(
"Utilities::Transpose failed, perhaps because matrix is not a Crs matrix");
563 throw Exceptions::RuntimeError(
"Utilities::Transpose: Tpetra is not compiled with LO=GO=int. Add TPETRA_INST_INT_INT:BOOL=ON to your configuration!");
566 throw Exceptions::RuntimeError(
"Utilities::Transpose: Tpetra is not compiled!");
572 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT) 576 EpetraExt::RowMatrix_Transpose transposer;
578 transposer.ReleaseTranspose();
580 RCP<Epetra_CrsMatrix> rcpA(A);
582 RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
583 RCP<CrsMatrixWrap> AAAA =
rcp(
new CrsMatrixWrap(AAA));
584 AAAA->fillComplete(Op.getRangeMap(), Op.getDomainMap());
588 throw Exceptions::RuntimeError(
"Epetra (Err. 2)");
593 throw Exceptions::RuntimeError(
"Only Epetra and Tpetra matrices can be transposed.");
597 return Teuchos::null;
602 static RCP<Xpetra::MultiVector<double,LO,GO,NO> > ExtractCoordinatesFromParameterList(ParameterList& paramList) {
603 RCP<Xpetra::MultiVector<double,LO,GO,NO> > coordinates = Teuchos::null;
606 if(paramList.isParameter (
"Coordinates") ==
false)
609 #if defined(HAVE_MUELU_TPETRA) 610 #if ( defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_INT)) || \ 611 (!defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT)) 616 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT) 617 typedef Tpetra::MultiVector<float, LO,GO,NO> tfMV;
618 RCP<tfMV> floatCoords = Teuchos::null;
624 typedef Tpetra::MultiVector<SC,LO,GO,NO> tdMV;
625 RCP<tdMV> doubleCoords = Teuchos::null;
626 if (paramList.isType<RCP<tdMV> >(
"Coordinates")) {
628 doubleCoords = paramList.get<RCP<tdMV> >(
"Coordinates");
629 paramList.remove(
"Coordinates");
631 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT) 632 else if (paramList.isType<RCP<tfMV> >(
"Coordinates")) {
634 floatCoords = paramList.get<RCP<tfMV> >(
"Coordinates");
635 paramList.remove(
"Coordinates");
636 doubleCoords =
rcp(
new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
637 deep_copy(*doubleCoords, *floatCoords);
641 if(doubleCoords != Teuchos::null) {
645 #endif // Tpetra instantiated on GO=int and EpetraNode 646 #endif // endif HAVE_TPETRA 648 #if defined(HAVE_MUELU_EPETRA) 649 RCP<Epetra_MultiVector> doubleEpCoords;
650 if (paramList.isType<RCP<Epetra_MultiVector> >(
"Coordinates")) {
651 doubleEpCoords = paramList.get<RCP<Epetra_MultiVector> >(
"Coordinates");
652 paramList.remove(
"Coordinates");
666 template <
class View,
unsigned AppendValue >
671 template <
class MT,
unsigned T >
672 struct CombineMemoryTraits {
676 template <
unsigned U,
unsigned T>
677 struct CombineMemoryTraits<
Kokkos::MemoryTraits<U>, T> {
678 typedef Kokkos::MemoryTraits<U|T> type;
681 template <
class DataType,
unsigned T,
class... Pack >
682 struct AppendTrait<
Kokkos::
View< DataType, Pack... >, T> {
683 typedef Kokkos::View< DataType, Pack... > view_type;
684 using type = Kokkos::View< DataType, typename view_type::array_layout, typename view_type::device_type, typename CombineMemoryTraits<typename view_type::memory_traits,T>::type >;
689 #define MUELU_UTILITIES_KOKKOS_SHORT 691 #endif // #if defined(HAVE_MUELU_KOKKOS_REFACTOR) 693 #endif // MUELU_UTILITIES_KOKKOS_DECL_HPP const Epetra_Map & RowMap() const
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< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string())
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.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
bool is_null(const std::shared_ptr< T > &p)
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 Teuchos::ScalarTraits< Scalar >::magnitudeType Distance2(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &v, LocalOrdinal i0, LocalOrdinal i1)
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 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.
int ExtractMyRowView(int MyRow, int &NumEntries, double *&Values, int *&Indices) const
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< Time > getNewTimer(const std::string &name)
Kokkos::View< const bool *, typename NO::device_type > DetectDirichletRows(const Xpetra::Matrix< SC, LO, GO, NO > &A, const typename Teuchos::ScalarTraits< SC >::magnitudeType &tol)
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)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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)
int NumMyElements() const
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV2(Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec)
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)
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.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
static Teuchos::ArrayRCP< Scalar > GetMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal.