46 #ifndef BELOS_MUELU_ADAPTER_HPP 
   47 #define BELOS_MUELU_ADAPTER_HPP 
   51 #ifdef HAVE_XPETRA_EPETRA 
   52 #include <Epetra_config.h> 
   53 #include <BelosOperator.hpp> 
   60 #ifdef HAVE_MUELU_AMGX 
   64 #include <BelosOperatorT.hpp> 
   67 #include "MueLu_Hierarchy.hpp" 
   71   using Teuchos::rcpFromRef;
 
   99             class Node          = KokkosClassic::DefaultNode::DefaultNodeType>
 
  101       public OperatorT<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
 
  102 #ifdef HAVE_XPETRA_TPETRA 
  103     , 
public OperatorT<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
 
  113 #ifdef HAVE_MUELU_AMGX 
  128     void Apply(
const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x, Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& y, ETrans trans = NOTRANS )
 const {
 
  131                          "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
 
  136 #ifdef HAVE_MUELU_AMGX 
  138         Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tX = Xpetra::toTpetra(x);
 
  139         Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tY = Xpetra::toTpetra(y);
 
  141         AMGX_->apply(tX, tY);
 
  150 #ifdef HAVE_XPETRA_TPETRA 
  157     void Apply(
const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x, Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& y, ETrans trans = NOTRANS )
 const {
 
  160                          "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
 
  165 #ifdef HAVE_MUELU_AMGX 
  171         Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> & temp_x = 
const_cast<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
>(x);
 
  173         const Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tX(rcpFromRef(temp_x));
 
  174         Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tY(rcpFromRef(y));
 
  182 #ifdef HAVE_MUELU_AMGX 
  187 #ifdef HAVE_XPETRA_EPETRA 
  188 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES 
  202     public OperatorT<Xpetra::MultiVector<double, int, int, Xpetra::EpetraNode> >
 
  203 #ifdef HAVE_XPETRA_TPETRA 
  205 #if ((defined(EPETRA_HAVE_OMP) && (defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_INT))) || \ 
  206     (!defined(EPETRA_HAVE_OMP) && (defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT)))) 
  207     , 
public OperatorT<Tpetra::MultiVector<double, int, int, Xpetra::EpetraNode> >
 
  210 #ifdef HAVE_XPETRA_EPETRA 
  211     , 
public OperatorT<Epetra_MultiVector>
 
  212     , 
public Belos::Operator<double>
 
  223 #ifdef HAVE_MUELU_AMGX 
  228     void Apply(
const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x, Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& y, ETrans trans = NOTRANS )
 const {
 
  231                          "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
 
  236 #ifdef HAVE_MUELU_AMGX 
  238         Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tX = Xpetra::toTpetra(x);
 
  239         Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tY = Xpetra::toTpetra(y);
 
  241         AMGX_->apply(tX, tY);
 
  248 #ifdef HAVE_XPETRA_TPETRA 
  249 #if ((defined(EPETRA_HAVE_OMP) && (defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_INT))) || \ 
  250     (!defined(EPETRA_HAVE_OMP) && (defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT)))) 
  251     void Apply ( 
const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x, Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& y, ETrans trans=NOTRANS )
 const {
 
  253                          "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
 
  258 #ifdef HAVE_MUELU_AMGX 
  264         Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> & temp_x = 
const_cast<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
>(x);
 
  266         const Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>  tX(rcpFromRef(temp_x));
 
  267         Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>        tY(rcpFromRef(y));
 
  277 #ifdef HAVE_XPETRA_EPETRA 
  286                          "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
 
  290       const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> tX(rcpFromRef(temp_x));
 
  291       Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node>       tY(rcpFromRef(y));
 
  304     void Apply(
const Belos::MultiVec<double>& x, Belos::MultiVec<double>& y, ETrans trans = NOTRANS )
 const {
 
  309                                  "Belos::MueLuOp::Apply, x and/or y cannot be dynamic cast to an Epetra_MultiVector.");
 
  311       Apply(*vec_x, *vec_y, trans);
 
  316     RCP<MueLu::Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node> >    
Hierarchy_;
 
  317 #ifdef HAVE_MUELU_AMGX 
  318     RCP<MueLu::AMGXOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node> > 
AMGX_;
 
  321 #endif // !EPETRA_NO_32BIT_GLOBAL_INDICES 
  322 #endif // HAVE_XPETRA_EPETRA 
  325 #ifdef HAVE_XPETRA_EPETRA 
  326 #ifndef  EPETRA_NO_64BIT_GLOBAL_INDICES 
  339   class MueLuOp<double, int, long long, Xpetra::
EpetraNode> :
 
  340     public OperatorT<Xpetra::MultiVector<double, int, long long, Xpetra::EpetraNode> >
 
  341 #ifdef HAVE_XPETRA_TPETRA 
  343 #if ((defined(EPETRA_HAVE_OMP) && (defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \ 
  344     (!defined(EPETRA_HAVE_OMP) && (defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_LONG_LONG)))) 
  345     , 
public OperatorT<Tpetra::MultiVector<double, int, long long, Xpetra::EpetraNode> >
 
  348 #ifdef HAVE_XPETRA_EPETRA 
  349     , 
public OperatorT<Epetra_MultiVector>
 
  350     , 
public Belos::Operator<double>
 
  361 #ifdef HAVE_MUELU_AMGX 
  366     void Apply(
const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x, Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& y, ETrans trans = NOTRANS )
 const {
 
  369                          "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
 
  374 #ifdef HAVE_MUELU_AMGX 
  376         Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tX = Xpetra::toTpetra(x);
 
  377         Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tY = Xpetra::toTpetra(y);
 
  379         AMGX_->apply(tX, tY);
 
  386 #ifdef HAVE_XPETRA_TPETRA 
  387 #if ((defined(EPETRA_HAVE_OMP) && (defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \ 
  388     (!defined(EPETRA_HAVE_OMP) && (defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_LONG_LONG)))) 
  389     void Apply ( 
const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x, Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& y, ETrans trans=NOTRANS )
 const {
 
  391                          "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
 
  396 #ifdef HAVE_MUELU_AMGX 
  402         Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> & temp_x = 
const_cast<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
>(x);
 
  404         const Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>  tX(rcpFromRef(temp_x));
 
  405         Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>        tY(rcpFromRef(y));
 
  415 #ifdef HAVE_XPETRA_EPETRA 
  424                          "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
 
  428       const Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node> tX(rcpFromRef(temp_x));
 
  429       Xpetra::EpetraMultiVectorT<GlobalOrdinal,Node>       tY(rcpFromRef(y));
 
  442     void Apply(
const Belos::MultiVec<double>& x, Belos::MultiVec<double>& y, ETrans trans = NOTRANS )
 const {
 
  447                                  "Belos::MueLuOp::Apply, x and/or y cannot be dynamic cast to an Epetra_MultiVector.");
 
  449       Apply(*vec_x, *vec_y, trans);
 
  454     RCP<MueLu::Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node> >    
Hierarchy_;
 
  455 #ifdef HAVE_MUELU_AMGX 
  456     RCP<MueLu::AMGXOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node> > 
AMGX_;
 
  459 #endif // !EPETRA_NO_64BIT_GLOBAL_INDICES 
  460 #endif // HAVE_XPETRA_EPETRA 
  463 #endif // BELOS_MUELU_ADAPTER_HPP 
MueLu::DefaultLocalOrdinal LocalOrdinal
 
MueLuOpFailure is thrown when a return value from an MueLu call on an Xpetra::Operator or MueLu::Hier...
 
MueLuOp derives from Belos::OperatorT and administrates a MueLu::Hierarchy. It implements the apply c...
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
void Apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &y, ETrans trans=NOTRANS) const 
This routine takes the Xpetra::MultiVector x and applies the operator to it resulting in the Xpetra::...
 
RCP< MueLu::AMGXOperator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > AMGX_
 
MueLuOpFailure(const std::string &what_arg)
 
MueLu::DefaultScalar Scalar
 
MueLu::DefaultGlobalOrdinal GlobalOrdinal
 
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Hierarchy_
 
MueLuOp(const RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &H)
Default constructor. 
 
virtual ~MueLuOp()
Destructor. 
 
MueLuOp(const RCP< MueLu::AMGXOperator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A)
 
Kokkos::Compat::KokkosSerialWrapperNode EpetraNode
 
Provides methods to build a multigrid hierarchy and apply multigrid cycles. 
 
Adapter for AmgX library from Nvidia.