47 #ifndef XPETRA_BLOCKEDVECTOR_HPP
48 #define XPETRA_BLOCKEDVECTOR_HPP
60 #ifndef DOXYGEN_SHOULD_SKIP_THIS
62 template<
class S,
class LO,
class GO,
class N>
class Vector;
63 template<
class S,
class LO,
class GO,
class N>
class MapExtractor;
64 template<
class S,
class LO,
class GO,
class N>
class VectorFactory;
67 template <
class Scalar = double,
68 class LocalOrdinal = Map<>::local_ordinal_type,
69 class GlobalOrdinal =
typename Map<LocalOrdinal>::global_ordinal_type,
70 class Node =
typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
72 :
public virtual Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >,
93 #undef XPETRA_BLOCKEDVECTOR_SHORT
295 virtual void scale(
const Scalar &alpha) {
415 throw Xpetra::Exceptions::RuntimeError(
"BlockedVector::getLocalLength: routine not implemented. It has no value as one must iterate on the partial vectors.");
421 return this->
getBlockedMap()->getFullMap()->getGlobalNumElements();
437 return std::string(
"BlockedVector");
443 for(
size_t r = 0; r < this->
getBlockedMap()->getNumMaps(); r++)
478 for(
size_t r = 0; r < this->
getBlockedMap()->getNumMaps(); ++r) {
484 virtual void randomize(
bool bUseXpetraImplementation =
false) {
485 for(
size_t r = 0; r < this->
getBlockedMap()->getNumMaps(); ++r) {
496 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
497 typedef typename Kokkos::Details::ArithTraits<Scalar>::val_type impl_scalar_type;
498 typedef Kokkos::DualView<impl_scalar_type**, Kokkos::LayoutStride,
499 typename node_type::execution_space,
500 Kokkos::MemoryUnmanaged> dual_view_type;
501 typedef typename dual_view_type::host_mirror_space host_execution_space;
502 typedef typename dual_view_type::t_dev::execution_space dev_execution_space;
509 template<
class TargetDeviceType>
510 typename Kokkos::Impl::if_c<
511 Kokkos::Impl::is_same<
512 typename dev_execution_space::memory_space,
513 typename TargetDeviceType::memory_space>::value,
514 typename dual_view_type::t_dev_um,
515 typename dual_view_type::t_host_um>::type
516 getLocalView ()
const {
517 if(Kokkos::Impl::is_same<
518 typename host_execution_space::memory_space,
519 typename TargetDeviceType::memory_space
521 return getHostLocalView();
523 return getDeviceLocalView();
527 virtual typename dual_view_type::t_host_um getHostLocalView ()
const {
528 typename dual_view_type::t_host_um test;
531 virtual typename dual_view_type::t_dev_um getDeviceLocalView()
const {
532 typename dual_view_type::t_dev_um test;
587 #define XPETRA_BLOCKEDVECTOR_SHORT
588 #endif // XPETRA_BLOCKEDVECTOR_HPP
virtual void replaceGlobalValue(GlobalOrdinal, size_t, const Scalar &)
Replace value, using global (row) index.
virtual void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
virtual void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Export &, CombineMode)
Export (using an Importer).
Teuchos::RCP< MultiVector > getMultiVector(size_t r) const
return partial multivector associated with block row r
virtual void reciprocal(const MultiVector &)
Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
virtual void putScalar(const Scalar &value)
Set all values in the vector with the given value.
virtual void dot(const MultiVector &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
virtual void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Import &, CombineMode)
Export.
virtual void sumIntoGlobalValue(GlobalOrdinal, size_t, const Scalar &)
Add value to existing value, using global (row) index.
virtual Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this vector.
virtual void replaceLocalValue(LocalOrdinal, size_t, const Scalar &)
Replace value, using local (row) index.
virtual void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Replace value, using global (row) index.
virtual void abs(const MultiVector &A)
Put element-wise absolute values of input vector in target: A = abs(this).
void setMultiVector(size_t r, Teuchos::RCP< const Vector > v, bool bThyraMode)
set partial Vector associated with block row r
GlobalOrdinal global_ordinal_type
virtual void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value)
Add value to existing value, using local (row) index.
BlockedVector(Teuchos::RCP< const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > > bmap, Teuchos::RCP< Vector > v)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void Xpetra_randomize()
Set multi-vector values to random numbers. XPetra implementation.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of vector.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const
Compute 1-norm of vector.
Teuchos::RCP< const Map > getMap() const
Access function for the underlying Map this DistObject was constructed with.
virtual void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a vector, this[j] = alpha[j]*this[j].
Exception throws to report errors in the internal logical of the program.
virtual void replaceLocalValue(LocalOrdinal myRow, const Scalar &value)
Replace value, using local (row) index.
virtual void Xpetra_randomize()
Set vector values to random numbers. XPetra implementation.
virtual void update(const Scalar &alpha, const MultiVector &A, const Scalar &beta)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
virtual void scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
virtual void multiply(Teuchos::ETransp, Teuchos::ETransp, const Scalar &, const Vector &, const Vector &, const Scalar &)
Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
virtual void replaceMap(const RCP< const Map > &map)
virtual Teuchos::RCP< const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
virtual void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
virtual void update(const Scalar &alpha, const MultiVector &A, const Scalar &beta)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
virtual Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
virtual void abs(const MultiVector &)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
virtual void assign(const MultiVector &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
virtual void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Export &, CombineMode)
Import (using an Exporter).
virtual void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Add value to existing value, using global (row) index.
virtual void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Import &, CombineMode)
Import.
BlockedVector(Teuchos::RCP< const Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mapExtractor, Teuchos::RCP< Vector > v)
BlockedVector(Teuchos::RCP< const Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mapExtractor, Teuchos::RCP< const Vector > v)
virtual void setSeed(unsigned int seed)
Set seed for Random function.
virtual void multiply(Teuchos::ETransp, Teuchos::ETransp, const Scalar &, const MultiVector &, const MultiVector &, const Scalar &)
Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
virtual void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
virtual global_size_t getGlobalLength() const
Global number of rows in the Vector.
virtual void elementWiseMultiply(Scalar, const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const MultiVector &, Scalar)
Element-wise multiply of a Vector A with a MultiVector B.
virtual void replaceMap(const RCP< const Map > &map)
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const =0
The Map describing the parallel distribution of this object.
virtual Teuchos::RCP< const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
virtual size_t getLocalLength() const
Local number of rows on the calling process.
LocalOrdinal local_ordinal_type
virtual void assign(const MultiVector &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void setMultiVector(size_t r, Teuchos::RCP< const MultiVector > v, bool bThyraMode)
set partial multivector associated with block row r
Teuchos::RCP< MultiVector > Merge() const
merge BlockedVector blocks to a single Vector
BlockedVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & operator=(const MultiVector &rhs)
Assignment operator: Does a deep copy.
virtual void sumIntoLocalValue(LocalOrdinal, size_t, const Scalar &)
Add value to existing value, using local (row) index.
virtual void meanValue(const Teuchos::ArrayView< Scalar > &) const
Compute mean (average) value of each vector in vector. The outcome of this routine is undefined for n...
Teuchos::RCP< MultiVector > getMultiVector(size_t r, bool bThyraMode) const
return partial Vector associated with block row r
virtual void scale(const Scalar &alpha)
Scale the current values of a vector, this = alpha*this.
size_t global_size_t
Global size_t object.
virtual Scalar meanValue() const
Compute mean (average) value of this Vector.
static const EVerbosityLevel verbLevel_default
virtual void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
virtual Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const
Compute Inf-norm in vector.
BlockedVector(const Teuchos::RCP< const BlockedMap > &map, bool zeroOut=true)
Constructor.
Teuchos::RCP< MultiVector > Merge() const
merge BlockedMultiVector blocks to a single MultiVector
virtual Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this vector.
virtual void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
virtual void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
virtual void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
virtual size_t getNumVectors() const
Number of columns in the Vector.
virtual void update(const Scalar &alpha, const MultiVector &A, const Scalar &beta, const MultiVector &B, const Scalar &gamma)
Update vector with scaled values of A and B, this = gamma*this + alpha*A + beta*B.
#define XPETRA_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual bool isSameSize(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &) const
Local number of rows on the calling process.
virtual void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
virtual void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
virtual Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
virtual void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
CombineMode
Xpetra::Combine Mode enumerable type.
virtual ~BlockedVector()
Destructor.
#define XPETRA_MONITOR(funcName)
Teuchos::RCP< MultiVector > getMultiVector(size_t r) const
return partial Vector associated with block row r
virtual std::string description() const
A simple one-line description of this object.
virtual Scalar dot(const Vector &A) const
Computes dot product of this Vector against input Vector x.
virtual void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
virtual void dot(const MultiVector &, const Teuchos::ArrayView< Scalar > &) const
Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
BlockedVector(Teuchos::RCP< const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > > bmap, Teuchos::RCP< const Vector > v)
Teuchos::RCP< const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > > getBlockedMap() const
Access function for the underlying Map this DistObject was constructed with.
virtual Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
virtual void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
virtual void elementWiseMultiply(Scalar, const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Vector &B, Scalar)
Element-wise multiply of a Vector A with a Vector B.
virtual void reciprocal(const MultiVector &A)
Put element-wise reciprocal values of input vector in target, this(i,j) = 1/A(i,j).