1 #ifndef MUELU_CREATE_TPETRA_PRECONDITIONER_HPP
2 #define MUELU_CREATE_TPETRA_PRECONDITIONER_HPP
8 #include <Tpetra_Operator.hpp>
9 #include <Tpetra_RowMatrix.hpp>
10 #include <Xpetra_TpetraBlockCrsMatrix.hpp>
11 #include <Tpetra_BlockCrsMatrix.hpp>
12 #include <Xpetra_CrsMatrix.hpp>
13 #include <Xpetra_MultiVector.hpp>
14 #include <Xpetra_MultiVectorFactory.hpp>
19 #include <MueLu_Hierarchy.hpp>
21 #include <MueLu_MLParameterListInterpreter.hpp>
22 #include <MueLu_ParameterListInterpreter.hpp>
23 #include <MueLu_TpetraOperator.hpp>
25 #include <MueLu_Utilities.hpp>
26 #include <MueLu_HierarchyUtils.hpp>
29 #if defined(HAVE_MUELU_AMGX)
32 #include "cuda_runtime.h"
45 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
57 typedef Xpetra::MultiVector<SC,LO,GO,NO> MultiVector;
58 typedef Xpetra::Matrix<SC,LO,GO,NO> Matrix;
60 typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> crs_matrix_type;
61 typedef Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> block_crs_matrix_type;
63 #if defined(HAVE_MUELU_AMGX)
64 std::string externalMG =
"use external multigrid package";
65 if (inParamList.
isParameter(externalMG) && inParamList.
get<std::string>(externalMG) ==
"amgx"){
76 if (crsA != Teuchos::null)
77 A = TpetraCrs_To_XpetraMatrix<SC,LO,GO,NO>(crsA);
78 else if (bcrsA != Teuchos::null) {
81 A =
rcp(
new Xpetra::CrsMatrixWrap<SC,LO,GO,NO>(temp));
108 RCP<Hierarchy> H = MueLu::CreateXpetraPreconditioner<SC,LO,GO,NO>(A, inParamList);
122 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
125 const std::string& xmlFileName)
129 return CreateTpetraPreconditioner<Scalar, LocalOrdinal, GlobalOrdinal, Node>(inA, paramList);
141 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
146 return CreateTpetraPreconditioner<Scalar, LocalOrdinal, GlobalOrdinal, Node>(inA, paramList);
157 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
165 typedef Xpetra::Matrix<SC,LO,GO,NO> Matrix;
169 RCP<Matrix> A = TpetraCrs_To_XpetraMatrix<SC,LO,GO,NO>(inA);
171 MueLu::ReuseXpetraPreconditioner<SC,LO,GO,NO>(A, H);
178 #endif //ifndef MUELU_CREATE_TPETRA_PRECONDITIONER_HPP
Various adapters that will create a MueLu preconditioner that is an Xpetra::Matrix.
MueLu::DefaultLocalOrdinal LocalOrdinal
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
bool isParameter(const std::string &name) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetHierarchy() const
Direct access to the underlying MueLu::Hierarchy.
Wraps an existing MueLu::Hierarchy as a Tpetra::Operator.
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Exception throws to report errors in the internal logical of the program.
void ReuseTpetraPreconditioner(const Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &inA, MueLu::TpetraOperator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op)
Helper function to reuse an existing MueLu preconditioner.
Teuchos::RCP< MueLu::TpetraOperator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > CreateTpetraPreconditioner(const Teuchos::RCP< Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &inA, Teuchos::ParameterList &inParamList)
Helper function to create a MueLu or AMGX preconditioner that can be used by Tpetra.Given a Tpetra::Operator, this function returns a constructed MueLu preconditioner.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Adapter for AmgX library from Nvidia.