46 #ifndef BELOS_MUELU_ADAPTER_HPP 47 #define BELOS_MUELU_ADAPTER_HPP 49 #ifdef HAVE_MUELU_EPETRA 50 #include <BelosOperator.hpp> 53 #ifdef HAVE_MUELU_AMGX 54 #include "MueLu_AMGXOperator.hpp" 60 #include "MueLu_Hierarchy.hpp" 64 using Teuchos::rcpFromRef;
78 template <
class Scalar,
79 class LocalOrdinal = int,
80 class GlobalOrdinal = LocalOrdinal,
81 class Node = KokkosClassic::DefaultNode::DefaultNodeType>
83 public OperatorT<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
84 #ifdef HAVE_MUELU_TPETRA 85 ,
public OperatorT<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
95 #ifdef HAVE_MUELU_AMGX 110 void Apply(
const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x, Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& y, ETrans trans = NOTRANS )
const {
113 "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
118 #ifdef HAVE_MUELU_AMGX 119 if (!AMGX_.is_null()) {
120 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tX = Xpetra::toTpetra(x);
121 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tY = Xpetra::toTpetra(y);
123 AMGX_->apply(tX, tY);
127 if (!Hierarchy_.is_null())
128 Hierarchy_->Iterate(x, y, 1,
true);
131 #ifdef HAVE_MUELU_TPETRA 138 void Apply(
const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x, Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& y, ETrans trans = NOTRANS )
const {
141 "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
146 #ifdef HAVE_MUELU_AMGX 147 if (!AMGX_.is_null())
151 if (!Hierarchy_.is_null()) {
152 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> & temp_x =
const_cast<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
>(x);
154 const Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tX(rcpFromRef(temp_x));
155 Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tY(rcpFromRef(y));
156 Hierarchy_->Iterate(tX, tY, 1,
true);
162 RCP<MueLu::Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
Hierarchy_;
163 #ifdef HAVE_MUELU_AMGX 164 RCP<MueLu::AMGXOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
AMGX_;
170 public OperatorT<Xpetra::MultiVector<double, int, int> >
171 #ifdef HAVE_MUELU_TPETRA 172 ,
public OperatorT<Tpetra::MultiVector<double, int, int> >
174 #ifdef HAVE_MUELU_EPETRA 176 ,
public Belos::Operator<double>
182 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal>::node_type
Node;
187 #ifdef HAVE_MUELU_AMGX 192 void Apply(
const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x, Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& y, ETrans trans = NOTRANS )
const {
195 "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
200 #ifdef HAVE_MUELU_AMGX 201 if (!AMGX_.is_null()) {
202 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tX = Xpetra::toTpetra(x);
203 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tY = Xpetra::toTpetra(y);
205 AMGX_->apply(tX, tY);
208 if (!Hierarchy_.is_null())
209 Hierarchy_->Iterate(x, y, 1,
true);
212 #ifdef HAVE_MUELU_TPETRA 213 void Apply (
const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x, Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& y, ETrans trans=NOTRANS )
const {
215 "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
220 #ifdef HAVE_MUELU_AMGX 221 if (!AMGX_.is_null())
225 if (!Hierarchy_.is_null()) {
226 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> & temp_x =
const_cast<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
>(x);
228 const Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tX(rcpFromRef(temp_x));
229 Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tY(rcpFromRef(y));
233 Hierarchy_->Iterate(tX, tY, 1,
true);
238 #ifdef HAVE_MUELU_EPETRA 245 void Apply(
const Epetra_MultiVector& x, Epetra_MultiVector& y, ETrans trans = NOTRANS)
const {
247 "Belos::MueLuOp::Apply, transpose mode != NOTRANS not supported by MueLu preconditionners.");
249 Epetra_MultiVector& temp_x =
const_cast<Epetra_MultiVector&
>(x);
251 const Xpetra::EpetraMultiVector tX(rcpFromRef(temp_x));
252 Xpetra::EpetraMultiVector tY(rcpFromRef(y));
257 Hierarchy_->Iterate(tX, tY, 1,
true);
265 void Apply(
const Belos::MultiVec<double>& x, Belos::MultiVec<double>& y, ETrans trans = NOTRANS )
const {
266 const Epetra_MultiVector* vec_x =
dynamic_cast<const Epetra_MultiVector*
>(&x);
267 Epetra_MultiVector* vec_y =
dynamic_cast<Epetra_MultiVector*
>(&y);
269 TEUCHOS_TEST_FOR_EXCEPTION(vec_x==NULL || vec_y==NULL,
MueLuOpFailure,
270 "Belos::MueLuOp::Apply, x and/or y cannot be dynamic cast to an Epetra_MultiVector.");
272 Apply(*vec_x, *vec_y, trans);
277 RCP<MueLu::Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
Hierarchy_;
278 #ifdef HAVE_MUELU_AMGX 279 RCP<MueLu::AMGXOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
AMGX_;
285 #endif // BELOS_MUELU_ADAPTER_HPP
void Apply(const Belos::MultiVec< double > &x, Belos::MultiVec< double > &y, ETrans trans=NOTRANS) const
This routine takes the Belos::MultiVec x and applies the operator to it resulting in the Belos::Multi...
void Apply(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &y, ETrans trans=NOTRANS) const
This routine takes the Tpetra::MultiVector x and applies the operator to it resulting in the Tpetra::...
MueLuOpFailure is thrown when a return value from an MueLu call on an Xpetra::Operator or MueLu::Hier...
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal >::node_type Node
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)
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Hierarchy_
void Apply(const Epetra_MultiVector &x, Epetra_MultiVector &y, ETrans trans=NOTRANS) const
This routine takes the Tpetra::MultiVector x and applies the operator to it resulting in the Tpetra::...
MueLuOp(const RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &H)
Default constructor.
RCP< MueLu::AMGXOperator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > AMGX_
void Apply(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &y, ETrans trans=NOTRANS) const
void Apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &y, ETrans trans=NOTRANS) const
MueLuOp(const RCP< MueLu::AMGXOperator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A)
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Hierarchy_
MueLuOp(const RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &H)
MueLuOp(const RCP< MueLu::AMGXOperator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A)
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Adapter for AmgX library from Nvidia.