46 #ifndef XPETRA_TPETRAMULTIVECTOR_DECL_HPP
47 #define XPETRA_TPETRAMULTIVECTOR_DECL_HPP
55 #include "Tpetra_MultiVector.hpp"
56 #include "Tpetra_Vector.hpp"
62 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
63 const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> &
toTpetra(
const MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
65 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
66 Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> &
toTpetra(MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
68 #ifndef DOXYGEN_SHOULD_SKIP_THIS
70 template<
class S,
class LO,
class GO,
class N>
class TpetraVector;
75 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
76 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>
79 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > >
toXpetra(RCP<
const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
83 template <
class Scalar,
88 :
public virtual MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >
118 void replaceGlobalValue(GlobalOrdinal globalRow,
size_t vectorIndex,
const Scalar &value);
122 void sumIntoGlobalValue(GlobalOrdinal globalRow,
size_t vectorIndex,
const Scalar &value);
125 void replaceLocalValue(LocalOrdinal myRow,
size_t vectorIndex,
const Scalar &value);
128 void sumIntoLocalValue(LocalOrdinal myRow,
size_t vectorIndex,
const Scalar &value);
142 Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > >
getVector(
size_t j)
const;
145 Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > >
getVectorNonConst(
size_t j);
148 Teuchos::ArrayRCP< const Scalar >
getData(
size_t j)
const;
154 void get1dCopy(Teuchos::ArrayView< Scalar > A,
size_t LDA)
const;
157 void get2dCopy(Teuchos::ArrayView<
const Teuchos::ArrayView< Scalar > > ArrayOfPtrs)
const;
160 Teuchos::ArrayRCP< const Scalar >
get1dView()
const;
163 Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > >
get2dView()
const;
187 void scale(
const Scalar &alpha);
190 void scale(Teuchos::ArrayView< const Scalar > alpha);
199 void update(
const Scalar &alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A,
const Scalar &beta,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B,
const Scalar &gamma);
202 void norm1(
const Teuchos::ArrayView<
typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms)
const;
205 void norm2(
const Teuchos::ArrayView<
typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms)
const;
208 void normInf(
const Teuchos::ArrayView<
typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms)
const;
211 void meanValue(
const Teuchos::ArrayView< Scalar > &means)
const;
214 void multiply(Teuchos::ETransp transA, Teuchos::ETransp transB,
const Scalar &alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B,
const Scalar &beta);
240 void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
const;
249 void randomize(
bool bUseXpetraImplementation =
false);
254 Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal,Node> >
getMap()
const;
273 TpetraMultiVector(
const Teuchos::RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > &vec);
276 RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> >
getTpetra_MultiVector()
const;
279 void setSeed(
unsigned int seed);
282 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
295 template<
class TargetDeviceType>
296 typename Kokkos::Impl::if_c<
298 typename dual_view_type::t_dev_um::execution_space::memory_space,
299 typename TargetDeviceType::memory_space>::value,
300 typename dual_view_type::t_dev_um,
301 typename dual_view_type::t_host_um>::type
302 getLocalView ()
const;
304 typename dual_view_type::t_host_um getHostLocalView ()
const;
306 typename dual_view_type::t_dev_um
307 getDeviceLocalView()
const;
320 RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> >
vec_;
328 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
329 RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node > >
toXpetra(RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec) {
333 return Teuchos::null;
336 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
337 RCP<const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node > >
toXpetra(RCP<
const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec) {
341 return Teuchos::null;
344 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
348 return *tX.getTpetra_MultiVector();
351 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
355 return *tX.getTpetra_MultiVector();
360 #define XPETRA_TPETRAMULTIVECTOR_SHORT
361 #endif // XPETRA_TPETRAMULTIVECTOR_HPP
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
Teuchos::ArrayRCP< const Scalar > get1dView() const
Const persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const
Return const persisting pointers to values.
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]).
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst()
Return non-const persisting pointers to values.
std::string description() const
A simple one-line description of this object.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Element-wise multiply of a Vector A with a TpetraMultiVector B.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
size_t getNumVectors() const
Number of columns in the multivector.
void setSeed(unsigned int seed)
Set seed for Random function.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void reduce()
Sum values of a locally replicated multivector across all processes.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import data into this object using an Import object ("forward mode").
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 multiply(Teuchos::ETransp transA, Teuchos::ETransp transB, const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &beta)
Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
size_t getLocalLength() const
Local number of rows on the calling process.
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const
Fill the given array with a copy of this multivector's local values.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst()
Nonconst persisting (1-D) view of this multivector's local values.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
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 scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t global_size_t
Global size_t object.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic constuctor.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const
Fill the given array with a copy of this multivector's local values.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
virtual ~TpetraMultiVector()
Destructor (virtual for memory safety of derived classes).
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec_
The Tpetra::MultiVector which this class wraps.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export data into this object using an Import object ("reverse mode").
CombineMode
Xpetra::Combine Mode enumerable type.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra multivector.
TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass
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.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
global_size_t getGlobalLength() const
Global number of rows in the multivector.
void replaceMap(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)