46 #ifndef XPETRA_TPETRAVECTOR_HPP 47 #define XPETRA_TPETRAVECTOR_HPP 62 #include "Tpetra_Vector.hpp" 67 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
68 RCP<Tpetra::Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> >
toTpetra(Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
70 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
71 RCP<Tpetra::Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> >
toTpetra(
const Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
73 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
74 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > >
toXpetra(RCP<
const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
76 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
77 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > >
toXpetra(RCP<Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
82 template <class Scalar = Vector<>::scalar_type,
87 :
public virtual Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>,
88 public TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
90 #undef XPETRA_TPETRAMULTIVECTOR_SHORT 115 :
TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> (map,A,map->getNodeNumElements(),1) { }
183 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 187 typename dual_view_type::t_host_um getHostLocalView ()
const {
191 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
205 template<
class TargetDeviceType>
206 typename Kokkos::Impl::if_c<
207 Kokkos::Impl::is_same<
208 typename dual_view_type::t_dev_um::execution_space::memory_space,
209 typename TargetDeviceType::memory_space>::value,
210 typename dual_view_type::t_dev_um,
211 typename dual_view_type::t_host_um>::type
212 getLocalView ()
const {
222 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
226 return tX.getTpetra_Vector();
229 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
233 return tX.getTpetra_Vector();
236 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
241 return Teuchos::null;
244 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
247 return toXpetra(Teuchos::rcp_const_cast<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > (vec));
254 #define XPETRA_TPETRAVECTOR_SHORT 255 #endif // XPETRA_TPETRAVECTOR_HPP void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
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]).
TpetraVector(const Teuchos::RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &vec)
TpetraMultiVector constructor to wrap a Tpetra::MultiVector object.
LocalOrdinal local_ordinal_type
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra multivector.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
GlobalOrdinal global_ordinal_type
RCP< const CrsGraph< int, GlobalOrdinal > > toXpetra(const Epetra_CrsGraph &g)
Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of this Vector.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace 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...
void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Adds specified value to existing value at the specified location.
Teuchos::ScalarTraits< Scalar >::magnitudeType normWeighted(const Vector &weights) const
Compute Weighted 2-norm (RMS Norm) of this Vector.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Scalar dot(const Vector &a) const
Computes dot product of this Vector against input Vector x.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const
Return 1-norm of this Vector.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
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.
virtual ~TpetraVector()
Destructor.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
static const EVerbosityLevel verbLevel_default
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 norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_Vector() const
Get the underlying Tpetra multivector.
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.
TpetraVector(const Teuchos::RCP< const Map > &map, const Teuchos::ArrayView< const Scalar > &A)
Set multi-vector values from an array using Teuchos memory management classes. (copy) ...
std::string description() const
Return a simple one-line description of this object.
#define XPETRA_MONITOR(funcName)
TpetraVector(const Teuchos::RCP< const Map > &map, bool zeroOut=true)
Sets all vector entries to zero.
void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Replace current value at the specified location with specified value.
void normWeighted(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &weights, const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Scalar meanValue() const
Compute mean (average) value of this Vector.