10 #ifndef XPETRA_TPETRAVECTOR_DECL_HPP
11 #define XPETRA_TPETRAVECTOR_DECL_HPP
22 #include "Tpetra_Vector.hpp"
27 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
28 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
toTpetra(Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>&);
30 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
31 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
toTpetra(
const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>&);
33 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
34 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
toXpetra(RCP<
const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
36 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
37 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
42 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
44 :
public virtual Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
45 public TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
47 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dot;
48 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::norm1;
49 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::norm2;
50 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::normInf;
51 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::meanValue;
52 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceGlobalValue;
53 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::sumIntoGlobalValue;
54 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceLocalValue;
55 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::sumIntoLocalValue;
61 TpetraVector(
const Teuchos::RCP<
const Map<LocalOrdinal, GlobalOrdinal, Node>>& map,
bool zeroOut =
true);
64 TpetraVector(
const Teuchos::RCP<
const Map<LocalOrdinal, GlobalOrdinal, Node>>& map,
const Teuchos::ArrayView<const Scalar>& A);
92 typename Teuchos::ScalarTraits<Scalar>::magnitudeType
norm1()
const;
95 typename Teuchos::ScalarTraits<Scalar>::magnitudeType
norm2()
const;
98 typename Teuchos::ScalarTraits<Scalar>::magnitudeType
normInf()
const;
112 void describe(Teuchos::FancyOStream& out,
const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default)
const;
117 Scalar
dot(
const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& a)
const;
123 TpetraVector(
const Teuchos::RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& vec);
126 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
getTpetra_Vector()
const;
133 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
134 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
138 return tX.getTpetra_Vector();
141 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
142 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
146 return tX.getTpetra_Vector();
149 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
150 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
151 toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec) {
155 return Teuchos::null;
158 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
159 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
160 toXpetra(RCP<
const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec) {
162 return toXpetra(Teuchos::rcp_const_cast<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(vec));
167 #define XPETRA_TPETRAVECTOR_SHORT
168 #endif // XPETRA_TPETRAVECTOR_DECL_HPP
virtual ~TpetraVector()
Destructor.
Scalar dot(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &a) const
Computes dot product of this Vector against input Vector x.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of this Vector.
void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Adds specified value to existing value at the specified location.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const
Return 1-norm of this Vector.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
TpetraVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &map, bool zeroOut=true)
Sets all vector entries to zero.
void replaceLocalValue(LocalOrdinal myRow, const Scalar &value)
Replace current value at the specified location with specified values.
Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const
Compute Inf-norm of this Vector.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value)
Adds specified value to existing value at the specified location.
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_Vector() const
Get the underlying Tpetra multivector.
std::string description() const
Return a simple one-line description of this object.
void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Replace current value at the specified location with specified value.
Scalar meanValue() const
Compute mean (average) value of this Vector.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)