46 #ifndef XPETRA_TPETRAVECTOR_DECL_HPP 47 #define XPETRA_TPETRAVECTOR_DECL_HPP 58 #include "Tpetra_Vector.hpp" 63 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
64 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
toTpetra(Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>&);
66 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
67 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
toTpetra(
const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>&);
69 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
70 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
toXpetra(RCP<
const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
72 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
73 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
78 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node = KokkosClassic::DefaultNode::DefaultNodeType>
80 :
public virtual Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>
81 ,
public TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>
99 TpetraVector(
const Teuchos::RCP<
const Map<LocalOrdinal, GlobalOrdinal, Node>>& map,
bool zeroOut =
true);
102 TpetraVector(
const Teuchos::RCP<
const Map<LocalOrdinal, GlobalOrdinal, Node>>& map,
const Teuchos::ArrayView<const Scalar>& A);
130 typename Teuchos::ScalarTraits<Scalar>::magnitudeType
norm1()
const;
133 typename Teuchos::ScalarTraits<Scalar>::magnitudeType
norm2()
const;
136 typename Teuchos::ScalarTraits<Scalar>::magnitudeType
normInf()
const;
150 void describe(Teuchos::FancyOStream& out,
const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default)
const;
155 Scalar
dot(
const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& a)
const;
161 TpetraVector(
const Teuchos::RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& vec);
164 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
getTpetra_Vector()
const;
172 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
173 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
178 return tX.getTpetra_Vector();
181 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
182 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
187 return tX.getTpetra_Vector();
190 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
191 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
192 toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec)
197 return Teuchos::null;
200 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
201 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
202 toXpetra(RCP<
const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec)
205 return toXpetra(Teuchos::rcp_const_cast<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(vec));
210 #define XPETRA_TPETRAVECTOR_SHORT 211 #endif // XPETRA_TPETRAVECTOR_DECL_HPP
virtual ~TpetraVector()
Destructor.
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_Vector() const
Get the underlying Tpetra multivector.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of this Vector.
std::string description() const
Return a simple one-line description of this object.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const
Return 1-norm of this Vector.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Adds specified value to existing value at the specified location.
Scalar meanValue() const
Compute mean (average) value of this Vector.
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const
Compute Inf-norm of this Vector.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
void meanValue(const Teuchos::ArrayView< Scalar > &means) const
Compute mean (average) value of each vector in multi-vector. The outcome of this routine is undefined...
Scalar dot(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &a) const
Computes dot product of this Vector against input Vector x.
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.
#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.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value)
Adds specified value to existing value at the specified location.
void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Replace current value at the specified location with specified value.