MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_XpetraOperator_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_XPETRAOPERATOR_DECL_HPP
11#define MUELU_XPETRAOPERATOR_DECL_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14
15#include <Xpetra_Operator_fwd.hpp>
16#include <Xpetra_MultiVector_fwd.hpp>
17#include "MueLu_Level.hpp"
19
20namespace MueLu {
21
24template <class Scalar = DefaultScalar,
27 class Node = DefaultNode>
28class XpetraOperator : public Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
29 protected:
31
32 public:
34
35
38
40 virtual ~XpetraOperator();
41
43
45 const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getDomainMap() const;
46
48 const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getRangeMap() const;
49
51
55 void apply(const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& X,
56 Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Y,
57 Teuchos::ETransp mode = Teuchos::NO_TRANS,
58 Scalar /* alpha */ = Teuchos::ScalarTraits<Scalar>::one(),
59 Scalar /* beta */ = Teuchos::ScalarTraits<Scalar>::one()) const;
60
62 bool hasTransposeApply() const;
63
65 void residual(const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& X,
66 const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& B,
67 Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& R) const;
68
70
71
73 RCP<MueLu::Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node> > GetHierarchy() const;
74
76
77 private:
78 RCP<MueLu::Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Hierarchy_;
79};
80
81} // namespace MueLu
82
83#endif // MUELU_XPETRAOPERATOR_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Wraps an existing MueLu::Hierarchy as a Xpetra::Operator.
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
virtual ~XpetraOperator()
Destructor.
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Hierarchy_
void apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar=Teuchos::ScalarTraits< Scalar >::one(), Scalar=Teuchos::ScalarTraits< Scalar >::one()) const
Returns in Y the result of a Xpetra::Operator applied to a Xpetra::MultiVector X.
void residual(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const
Compute a residual R = B - (*this) * X.
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetHierarchy() const
Direct access to the underlying MueLu::Hierarchy.
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar