10#ifndef XPETRA_TPETRAMULTIVECTOR_DECL_HPP
11#define XPETRA_TPETRAMULTIVECTOR_DECL_HPP
19#include "Tpetra_MultiVector.hpp"
20#include "Tpetra_Vector.hpp"
26template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
27const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
toTpetra(
const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &);
29template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
30Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
toTpetra(MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &);
32#ifndef DOXYGEN_SHOULD_SKIP_THIS
34template <
class S,
class LO,
class GO,
class N>
39template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
40RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
42template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
43RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
toXpetra(RCP<
const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
45template <
class Scalar,
50 :
public virtual MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
76 void replaceGlobalValue(GlobalOrdinal globalRow,
size_t vectorIndex,
const Scalar &value);
79 void sumIntoGlobalValue(GlobalOrdinal globalRow,
size_t vectorIndex,
const Scalar &value);
82 void replaceLocalValue(LocalOrdinal myRow,
size_t vectorIndex,
const Scalar &value);
85 void sumIntoLocalValue(LocalOrdinal myRow,
size_t vectorIndex,
const Scalar &value);
143 void scale(
const Scalar &alpha);
155 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);
170 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);
205 void randomize(
bool bUseXpetraImplementation =
false);
207 void randomize(
const Scalar &minVal,
const Scalar &maxVal,
bool bUseXpetraImplementation =
false);
252 void setSeed(
unsigned int seed);
283template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
288 return Teuchos::null;
291template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
296 return Teuchos::null;
299template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
304 return Teuchos::null;
307template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
312 return Teuchos::null;
315template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
319 return *tX.getTpetra_MultiVector();
322template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
326 return *tX.getTpetra_MultiVector();
329template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
332 return Teuchos::rcp_dynamic_cast<Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(X,
true)->getTpetra_MultiVector();
335template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
338 return Teuchos::rcp_dynamic_cast<const Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(X,
true)->getTpetra_MultiVector();
343#define XPETRA_TPETRAMULTIVECTOR_SHORT
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
static const EVerbosityLevel verbLevel_default
Kokkos::DualView< impl_scalar_type **, Kokkos::LayoutStride, typename node_type::device_type, Kokkos::MemoryUnmanaged > dual_view_type
void replaceMap(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
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").
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
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.
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).
TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dual_view_type dual_view_type
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
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").
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void beginExport(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").
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 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]).
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst()
Nonconst persisting (1-D) view of this multivector's local values.
size_t getLocalLength() const
Local number of rows on the calling process.
size_t getNumVectors() const
Number of columns in the multivector.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
void setSeed(unsigned int seed)
Set seed for Random function.
void reduce()
Sum values of a locally replicated multivector across all processes.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
virtual dual_view_type::t_host_const_um getLocalViewHost(Access::ReadOnlyStruct) const
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
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 endExport(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").
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec_
The Tpetra::MultiVector which this class wraps.
std::string description() const
A simple one-line description of this object.
void scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
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.
virtual dual_view_type::t_dev_const_um getLocalViewDevice(Access::ReadOnlyStruct) const
global_size_t getGlobalLength() const
Global number of rows in the multivector.
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,...
virtual ~TpetraMultiVector()
Destructor (virtual for memory safety of derived classes).
void endImport(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 abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
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 sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
Teuchos::ArrayRCP< const Scalar > get1dView() const
Const persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra multivector.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const
Fill the given array with a copy of this multivector's local values.
void beginImport(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").
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst()
Return non-const persisting pointers to values.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const
Fill the given array with a copy of this multivector's local values.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const
Return const persisting pointers to values.
size_t global_size_t
Global size_t object.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
CombineMode
Xpetra::Combine Mode enumerable type.