46 #ifndef XPETRA_EPETRAVECTOR_HPP
47 #define XPETRA_EPETRAVECTOR_HPP
61 #include <Epetra_Vector.h>
66 template<
class GlobalOrdinal,
class Node>
69 template<
class GlobalOrdinal,
class Node>
73 template<
class EpetraGlobalOrdinal,
class Node>
75 :
public virtual Vector<double,int,EpetraGlobalOrdinal, Node>,
public EpetraMultiVectorT<EpetraGlobalOrdinal, Node>
183 "Xpetra::EpetraVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
195 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
198 typename dual_view_type::t_host_um getHostLocalView ()
const {
202 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
203 throw std::runtime_error(
"Epetra does not support device views!");
204 #ifndef __NVCC__ //prevent nvcc warning
205 typename dual_view_type::t_dev ret;
220 template<
class TargetDeviceType>
221 typename Kokkos::Impl::if_c<
222 Kokkos::Impl::is_same<
223 typename dual_view_type::t_dev_um::execution_space::memory_space,
224 typename TargetDeviceType::memory_space>::value,
225 typename dual_view_type::t_dev_um,
226 typename dual_view_type::t_host_um>::type
227 getLocalView ()
const {
228 return this->MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::template getLocalView<TargetDeviceType>();
243 #ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES
344 std::ostringstream oss;
383 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
386 typename dual_view_type::t_host_um getHostLocalView ()
const {
390 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
391 throw std::runtime_error(
"Epetra does not support device views!");
392 #ifndef __NVCC__ //prevent nvcc warning
393 typename dual_view_type::t_dev ret;
408 template<
class TargetDeviceType>
409 typename Kokkos::Impl::if_c<
410 Kokkos::Impl::is_same<
411 typename dual_view_type::t_dev_um::execution_space::memory_space,
412 typename TargetDeviceType::memory_space>::value,
413 typename dual_view_type::t_dev_um,
414 typename dual_view_type::t_host_um>::type
415 getLocalView ()
const {
416 return this->MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::template getLocalView<TargetDeviceType>();
430 #ifndef XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES
529 std::ostringstream oss;
568 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
571 typename dual_view_type::t_host_um getHostLocalView ()
const {
575 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
576 throw std::runtime_error(
"Epetra does not support device views!");
577 #ifndef __NVCC__ //prevent nvcc warning
578 typename dual_view_type::t_dev ret;
593 template<
class TargetDeviceType>
594 typename Kokkos::Impl::if_c<
595 Kokkos::Impl::is_same<
596 typename dual_view_type::t_dev_um::execution_space::memory_space,
597 typename TargetDeviceType::memory_space>::value,
598 typename dual_view_type::t_dev_um,
599 typename dual_view_type::t_host_um>::type
600 getLocalView ()
const {
601 return this->MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::template getLocalView<TargetDeviceType>();
616 #endif // XPETRA_EPETRAVECTOR_HPP
void replaceLocalValue(LocalOrdinal myRow, const Scalar &value)
Replace current value at the specified location with specified values.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value)
Adds specified value to existing value at the specified location.
EpetraVectorT(const Teuchos::RCP< const Map< int, GlobalOrdinal, Node > > &map, bool zeroOut=true)
Sets all vector entries to zero.
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.
virtual void Print(std::ostream &os) const
Scalar dot(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &a) const
Computes dot product of this Vector against input Vector x.
Scalar meanValue() const
Compute Weighted 2-norm (RMS Norm) of this Vector.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value)
Adds specified value to existing value at the specified location.
Exception throws to report errors in the internal logical of the program.
const Epetra_CrsGraph & toEpetra(const RCP< const CrsGraph< int, GlobalOrdinal, Node > > &graph)
EpetraVectorT(const RCP< Epetra_MultiVector > &mv, size_t j)
global_size_t getGlobalLength() const
Global number of rows in the multivector.
EpetraVectorT(const Teuchos::RCP< const Map< int, GlobalOrdinal, Node > > &map, bool zeroOut=true)
Sets all vector entries to zero.
std::string description() const
Return a simple one-line description of this object.
std::string description() const
Return a simple one-line description of this object.
Epetra_Vector * getEpetra_Vector() const
Get the underlying Epetra vector.
void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Replace current value at the specified location with specified value.
virtual ~EpetraVectorT()
Vector copy constructor.
EpetraGlobalOrdinal GlobalOrdinal
Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of this Vector.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
const RCP< const Epetra_MultiVector > internalRefToBaseMV_
virtual ~EpetraVectorT()
Vector copy constructor.
Epetra_Vector * getEpetra_Vector() const
Get the underlying Epetra vector.
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.
virtual std::string description() const
EpetraGlobalOrdinal GlobalOrdinal
EpetraVectorT(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, bool zeroOut=true)
Sets all vector entries to zero.
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute the dot product of each corresponding pair of vectors (columns) in A and B.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
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)
Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const
Return 1-norm of this Vector.
EpetraVectorT(const Teuchos::RCP< Epetra_Vector > &vec)
EpetraMultiVectorT constructor to wrap a Epetra_Vector object.
static const EVerbosityLevel verbLevel_default
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.
static magnitudeType magnitude(T a)
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.
void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Adds specified value to existing value at the specified location.
Scalar dot(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &a) const
Computes dot product of this Vector against input Vector x.
EpetraVectorT(const Teuchos::RCP< Epetra_Vector > &vec)
EpetraMultiVectorT constructor to wrap a Epetra_Vector object.
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
const RCP< const Epetra_MultiVector > internalRefToBaseMV_
virtual ~EpetraVectorT()
Vector copy constructor.
EpetraVectorT(const Teuchos::RCP< Epetra_Vector > &vec)
EpetraMultiVectorT constructor to wrap a Epetra_Vector object.
void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value)
Adds specified value to existing value at the specified location.
Scalar meanValue() const
Compute Weighted 2-norm (RMS Norm) of this Vector.
Epetra_Vector * getEpetra_Vector() const
Get the underlying Epetra vector.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of this Vector.
const RCP< const Epetra_MultiVector > internalRefToBaseMV_
void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Replace current value at the specified location with specified value.
RCP< Epetra_MultiVector > getEpetra_MultiVector() const
Get the underlying Epetra multivector.
void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Replace current value at the specified location with specified value.
#define XPETRA_MONITOR(funcName)
std::string description() const
Return a simple one-line description of this object.
EpetraVectorT(const RCP< Epetra_MultiVector > &mv, size_t j)
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of this Vector.
Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const
Compute Inf-norm of this Vector.
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 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 describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
Scalar meanValue() const
Compute Weighted 2-norm (RMS Norm) of this Vector.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.