46 #ifndef XPETRA_TPETRAMULTIVECTOR_DEF_HPP
47 #define XPETRA_TPETRAMULTIVECTOR_DEF_HPP
56 #include "Tpetra_MultiVector.hpp"
57 #include "Tpetra_Vector.hpp"
63 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
65 : vec_(Teuchos::
rcp(new Tpetra::
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(map), NumVectors, zeroOut))) {
67 TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument,
"Xpetra::TpetraMultiVector(map,numVecs,zeroOut): numVecs = " << NumVectors <<
" < 1.");
71 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
74 : vec_(Teuchos::
rcp(new Tpetra::
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(Tpetra::createCopy(
toTpetra(source))))) { }
77 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
80 : vec_(Teuchos::
rcp(new Tpetra::
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(map), A, LDA, NumVectors))) {
82 TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument,
"Xpetra::TpetraMultiVector(map,A,LDA,numVecs): numVecs = " << NumVectors <<
" < 1.");
86 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
89 : vec_(Teuchos::
rcp(new Tpetra::
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(map), ArrayOfPtrs, NumVectors))) {
91 TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument,
"Xpetra::TpetraMultiVector(map,ArrayOfPtrs,numVecs): numVecs = " << NumVectors <<
" < 1.");
96 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
101 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
103 replaceGlobalValue(GlobalOrdinal globalRow,
size_t vectorIndex,
const Scalar &value) {
XPETRA_MONITOR(
"TpetraMultiVector::replaceGlobalValue"); vec_->replaceGlobalValue(globalRow, vectorIndex, value); }
106 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
108 sumIntoGlobalValue(GlobalOrdinal globalRow,
size_t vectorIndex,
const Scalar &value) {
XPETRA_MONITOR(
"TpetraMultiVector::sumIntoGlobalValue"); vec_->sumIntoGlobalValue(globalRow, vectorIndex, value); }
111 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
113 replaceLocalValue(LocalOrdinal myRow,
size_t vectorIndex,
const Scalar &value) {
XPETRA_MONITOR(
"TpetraMultiVector::replaceLocalValue"); vec_->replaceLocalValue(myRow, vectorIndex, value); }
116 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
118 sumIntoLocalValue(LocalOrdinal myRow,
size_t vectorIndex,
const Scalar &value) {
XPETRA_MONITOR(
"TpetraMultiVector::sumIntoLocalValue"); vec_->sumIntoLocalValue(myRow, vectorIndex, value); }
121 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
126 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
131 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
137 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
143 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
149 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
155 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
161 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
167 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
173 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
179 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
185 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
191 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
196 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
201 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
206 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
211 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
216 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
221 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
226 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
228 update(
const Scalar &alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A,
const Scalar &beta,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B,
const Scalar &gamma) {
XPETRA_MONITOR(
"TpetraMultiVector::update"); vec_->update(alpha,
toTpetra(A), beta,
toTpetra(B), gamma); }
231 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
236 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
241 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
246 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
251 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
253 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) {
XPETRA_MONITOR(
"TpetraMultiVector::multiply"); vec_->multiply(transA, transB, alpha,
toTpetra(A),
toTpetra(B), beta); }
257 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
262 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
267 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
272 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
277 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
282 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
287 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
292 if(bUseXpetraImplementation)
298 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
303 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
310 this->getTpetra_MultiVector()->doImport(*v,
toTpetra(importer),
toTpetra(CM));
313 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
320 this->getTpetra_MultiVector()->doExport(*v,
toTpetra(importer),
toTpetra(CM));
324 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
331 this->getTpetra_MultiVector()->doImport(*v,
toTpetra(exporter),
toTpetra(CM));
335 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
342 this->getTpetra_MultiVector()->doExport(*v,
toTpetra(exporter),
toTpetra(CM));
346 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
350 this->getTpetra_MultiVector()->replaceMap(
toTpetra(map));
353 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
354 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
355 template<
class Node2>
366 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
371 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
377 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
382 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
394 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
395 template<
class TargetDeviceType>
396 typename Kokkos::Impl::if_c<
397 Kokkos::Impl::is_same<
398 typename dual_view_type::t_dev_um::execution_space::memory_space,
399 typename TargetDeviceType::memory_space>::value,
400 typename dual_view_type::t_dev_um,
401 typename dual_view_type::t_host_um>::type
408 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
409 typename TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dual_view_type::t_host_um
410 TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
411 getHostLocalView ()
const {
412 return subview(vec_->template getLocalView<
typename TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dual_view_type::host_mirror_space> (),
413 Kokkos::ALL(), Kokkos::ALL());
416 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
417 typename TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dual_view_type::t_dev_um
418 TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
419 getDeviceLocalView()
const {
420 return subview(vec_->template getLocalView<
typename TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dual_view_type::t_dev_um::execution_space> (),
421 Kokkos::ALL(), Kokkos::ALL());
428 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
433 const this_type* rhsPtr =
dynamic_cast<const this_type*
> (&rhs);
435 rhsPtr == NULL, std::invalid_argument,
"Xpetra::MultiVector::operator=:"
436 " The left-hand side (LHS) of the assignment has a different type than "
437 "the right-hand side (RHS). The LHS has type Xpetra::TpetraMultiVector"
438 " (which means it wraps a Tpetra::MultiVector), but the RHS has some "
439 "other type. This probably means that the RHS wraps an "
440 "Epetra_MultiVector. Xpetra::MultiVector does not currently implement "
441 "assignment from an Epetra object to a Tpetra object, though this could"
442 " be added with sufficient interest.");
444 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> TMV;
446 RCP<TMV> lhsImpl = this->getTpetra_MultiVector ();
449 rhsImpl.
is_null (), std::logic_error,
"Xpetra::MultiVector::operator= "
450 "(in Xpetra::TpetraMultiVector::assign): *this (the right-hand side of "
451 "the assignment) has a null RCP<Tpetra::MultiVector> inside. Please "
452 "report this bug to the Xpetra developers.");
454 lhsImpl.
is_null (), std::logic_error,
"Xpetra::MultiVector::operator= "
455 "(in Xpetra::TpetraMultiVector::assign): The left-hand side of the "
456 "assignment has a null RCP<Tpetra::MultiVector> inside. Please report "
457 "this bug to the Xpetra developers.");
459 Tpetra::deep_copy (*lhsImpl, *rhsImpl);
463 #ifdef HAVE_XPETRA_EPETRA
465 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
466 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
469 template <
class Scalar>
471 :
public virtual MultiVector< Scalar, int, int, EpetraNode >
599 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) { }
614 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) { }
650 void randomize(
bool bUseXpetraImplementation =
false) { }
667 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
668 template<
class Node2>
688 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
703 template<
class TargetDeviceType>
704 typename Kokkos::Impl::if_c<
705 Kokkos::Impl::is_same<
706 typename dual_view_type::t_dev_um::execution_space::memory_space,
707 typename TargetDeviceType::memory_space>::value,
708 typename dual_view_type::t_dev_um,
709 typename dual_view_type::t_host_um>::type
710 getLocalView ()
const {
711 typename Kokkos::Impl::if_c<
712 Kokkos::Impl::is_same<
713 typename dual_view_type::t_dev_um::execution_space::memory_space,
714 typename TargetDeviceType::memory_space>::value,
715 typename dual_view_type::t_dev_um,
716 typename dual_view_type::t_host_um>::type dummy;
720 typename dual_view_type::t_host_um getHostLocalView ()
const {
723 return typename dual_view_type::t_host_um();
726 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
729 return typename dual_view_type::t_dev_um();
745 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
746 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
749 template <
class Scalar>
751 :
public virtual MultiVector< Scalar, int, long long, EpetraNode >
874 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) { }
889 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) { }
925 void randomize(
bool bUseXpetraImplementation =
false) { }
942 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
943 template<
class Node2>
963 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
978 template<
class TargetDeviceType>
979 typename Kokkos::Impl::if_c<
980 Kokkos::Impl::is_same<
981 typename dual_view_type::t_dev_um::execution_space::memory_space,
982 typename TargetDeviceType::memory_space>::value,
983 typename dual_view_type::t_dev_um,
984 typename dual_view_type::t_host_um>::type
985 getLocalView ()
const {
986 typename Kokkos::Impl::if_c<
987 Kokkos::Impl::is_same<
988 typename dual_view_type::t_dev_um::execution_space::memory_space,
989 typename TargetDeviceType::memory_space>::value,
990 typename dual_view_type::t_dev_um,
991 typename dual_view_type::t_host_um>::type dummy;
995 typename dual_view_type::t_host_um getHostLocalView ()
const {
998 return typename dual_view_type::t_host_um();
1001 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
1004 return typename dual_view_type::t_dev_um();
1019 #endif // TpetraMultiVector class (specialization GO=long long, NO=EpetraNode)
1021 #endif // HAVE_XPETRA_EPETRA
1028 #include "Xpetra_TpetraVector.hpp"
1031 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1032 void TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::elementWiseMultiply(Scalar scalarAB,
const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A,
const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &B, Scalar scalarThis) {
1038 XPETRA_DYNAMIC_CAST(
const tpv, A, tA,
"Xpetra::TpetraMultiVectorMatrix->multiply() only accept Xpetra::TpetraMultiVector as input arguments.");
1040 vec_->elementWiseMultiply(scalarAB, *tA.getTpetra_Vector(), *tB.getTpetra_MultiVector(), scalarThis);
1045 #define XPETRA_TPETRAMULTIVECTOR_SHORT
1046 #endif // XPETRA_TPETRAMULTIVECTOR_HPP
void replaceMap(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
size_t getLocalLength() const
Local number of rows on the calling process.
RCP< Map< LocalOrdinal, GlobalOrdinal, Node2 > > clone(const Map< LocalOrdinal, GlobalOrdinal, Node1 > &map, const RCP< Node2 > &node2)
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 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.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra multivector.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const
Return const persisting pointers to values.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
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 scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic constuctor.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void setSeed(unsigned int seed)
Set seed for Random function.
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::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst()
Return non-const persisting pointers to values.
void reduce()
Sum values of a locally replicated multivector across all processes.
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 normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
std::string description() const
A simple one-line description of this object.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
Create multivector by copying array of views of local data.
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).
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
virtual ~TpetraMultiVector()
Destructor (virtual for memory safety of derived classes).
void setSeed(unsigned int seed)
Set seed for Random function.
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.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic constuctor.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
size_t getNumVectors() const
Number of columns in the multivector.
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
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.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void Xpetra_randomize()
Set multi-vector values to random numbers. XPetra implementation.
std::string description() const
A simple one-line description of this object.
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.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Replace multi-vector values with scaled values of A, this = alpha*A.
void replaceMap(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const
Return const persisting pointers to values.
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.
TpetraMultiVector()
Default constructor.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
#define XPETRA_TPETRA_ETI_EXCEPTION(cl, obj, go, node)
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const
Fill the given array with a copy of this multivector's local values.
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 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...
TpetraMultiVector(const Teuchos::RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &vec)
TpetraMultiVector constructor to wrap a Tpetra::MultiVector object.
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const
Fill the given array with a copy of this multivector's local values.
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 norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
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).
Teuchos::ArrayRCP< Scalar > get1dViewNonConst()
Nonconst persisting (1-D) view of this multivector's local values.
TpetraMultiVector(const Teuchos::RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &vec)
TpetraMultiVector constructor to wrap a Tpetra::MultiVector object.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst()
Return non-const persisting pointers to values.
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)
Update multi-vector with scaled values of A and B, this = gamma*this + alpha*A + beta*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 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 putScalar(const Scalar &value)
Set all values in the multivector with the given value.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const
Fill the given array with a copy of this multivector's local values.
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 norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
std::string description() const
A simple one-line description of this object.
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
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]).
size_t getNumVectors() const
Number of columns in the multivector.
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.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
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 reduce()
Sum values of a locally replicated multivector across all processes.
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
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 sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
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).
virtual ~TpetraMultiVector()
Destructor (virtual for memory safety of derived classes).
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
static void seedrandom(unsigned int s)
#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.
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)
Update multi-vector with scaled values of A and B, this = gamma*this + alpha*A + beta*B.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
static const EVerbosityLevel verbLevel_default
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
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 doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const
Fill the given array with a copy of this multivector's local values.
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const
Fill the given array with a copy of this multivector's local values.
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).
global_size_t getGlobalLength() const
Global number of rows in the multivector.
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
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
Create multivector by copying array of views of local data.
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.
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 putScalar(const Scalar &value)
Set all values in the multivector with the given value.
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
size_t getLocalLength() const
Local number of rows on the calling process.
void scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Replace multi-vector values with scaled values of A, this = alpha*A.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Scalar > &A, size_t LDA, size_t NumVectors)
Create multivector by copying two-dimensional array of local data.
TpetraMultiVector(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
Copy constructor (performs a deep copy).
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Scalar > &A, size_t LDA, size_t NumVectors)
Create multivector by copying two-dimensional array of local data.
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...
Teuchos::ArrayRCP< const Scalar > get1dView() const
Const persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP< const Scalar > get1dView() const
Const persisting (1-D) view of this multivector's local values.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
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.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra multivector.
TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass
#define XPETRA_MONITOR(funcName)
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...
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst()
Nonconst persisting (1-D) view of this multivector's local values.
global_size_t getGlobalLength() const
Global number of rows in the multivector.
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 abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst()
Return non-const persisting pointers to values.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
TpetraMultiVector(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
Copy constructor (performs a deep copy).
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const
Return const persisting pointers to values.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
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 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.
global_size_t getGlobalLength() const
Global number of rows in the multivector.
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 replaceMap(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra multivector.