Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_TpetraOperator.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef XPETRA_TPETRAOPERATOR_HPP
11#define XPETRA_TPETRAOPERATOR_HPP
12
14
15#include <Tpetra_Operator.hpp>
17
18#include "Xpetra_Map.hpp"
19#include "Xpetra_TpetraMap.hpp"
20#include "Xpetra_MultiVector.hpp"
21#include "Xpetra_TpetraMultiVector.hpp"
22#include "Xpetra_Operator.hpp"
23
24#include "Xpetra_Utils.hpp"
25
26namespace Xpetra {
27
28template <class Scalar,
29 class LocalOrdinal,
30 class GlobalOrdinal,
31 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
32class TpetraOperator : public Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
33 public:
35
38 XPETRA_MONITOR("TpetraOperator::getDomainMap()");
39 return toXpetra(op_->getDomainMap());
40 }
41
44 XPETRA_MONITOR("TpetraOperator::getRangeMap()");
45 return toXpetra(op_->getRangeMap());
46 }
47
49
54 virtual void
62
64 virtual bool hasTransposeApply() const {
65 return op_->hasTransposeApply();
66 }
67
69
71
72
74 std::string description() const {
75 XPETRA_MONITOR("TpetraOperator::description");
76 return op_->description();
77 }
78
81 XPETRA_MONITOR("TpetraOperator::describe");
82 op_->describe(out, verbLevel);
83 }
84
86
88
89
91
95
98
101
108
110
111 private:
114
115}; // TpetraOperator class
116
117} // namespace Xpetra
118
119#define XPETRA_TPETRAOPERATOR_SHORT
120#endif // XPETRA_TPETRAOPERATOR_HPP
#define XPETRA_MONITOR(funcName)
static const EVerbosityLevel verbLevel_default
virtual const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
The Map associated with the domain of this operator, which must be compatible with X....
TpetraOperator(const Teuchos::RCP< Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &op)
TpetraOperator constructor to wrap a Tpetra::Operator object.
std::string description() const
A simple one-line description of this object.
RCP< Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > op_
The Tpetra::Operator which this class wraps.
virtual void apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Computes the operator-multivector application.
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.
virtual bool hasTransposeApply() const
Whether this operator supports applying the transpose or conjugate transpose.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
virtual RCP< Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getOperator()
Gets the operator out.
virtual RCP< const Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getOperatorConst() const
Gets the operator out.
virtual const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
The Map associated with the range of this operator, which must be compatible with Y....
void residual(const Operator< SC, LO, GO, NO > &Aop, const MultiVector< SC, LO, GO, NO > &X_in, const MultiVector< SC, LO, GO, NO > &B_in, MultiVector< SC, LO, GO, NO > &R_in)
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > graph)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)