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>
85 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dot;
86 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::norm1;
87 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::norm2;
88 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::normInf;
89 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::meanValue;
90 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceGlobalValue;
91 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::sumIntoGlobalValue;
92 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceLocalValue;
93 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::sumIntoLocalValue;
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;
166 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
170 typename dual_view_type::t_host_um getHostLocalView()
const;
172 typename dual_view_type::t_dev_um getDeviceLocalView()
const;
184 template<
class TargetDeviceType>
185 typename Kokkos::Impl::if_c<
186 std::is_same<typename dual_view_type::t_dev_um::execution_space::memory_space, typename TargetDeviceType::memory_space>::value,
187 typename dual_view_type::t_dev_um,
188 typename dual_view_type::t_host_um>::type
189 getLocalView()
const;
198 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
199 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
204 return tX.getTpetra_Vector();
207 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
208 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
213 return tX.getTpetra_Vector();
216 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
217 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
218 toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec)
223 return Teuchos::null;
226 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
227 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
228 toXpetra(RCP<
const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec)
231 return toXpetra(Teuchos::rcp_const_cast<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(vec));
236 #define XPETRA_TPETRAVECTOR_SHORT
237 #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.
#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)
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
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.