46 #ifndef XPETRA_EPETRAMULTIVECTOR_HPP
47 #define XPETRA_EPETRAMULTIVECTOR_HPP
51 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
52 #include <Kokkos_Core.hpp>
53 #include <Kokkos_DualView.hpp>
68 #include "Epetra_SerialComm.h"
70 #include <Epetra_MultiVector.h>
71 #include <Epetra_Vector.h>
76 template<
class GlobalOrdinal,
class Node>
78 template<
class GlobalOrdinal,
class Node>
80 template<
class GlobalOrdinal,
class Node>
81 RCP<MultiVector<double, int, GlobalOrdinal, Node> >
toXpetra(RCP<Epetra_MultiVector> vec);
84 #ifndef DOXYGEN_SHOULD_SKIP_THIS
85 template<
class GlobalOrdinal,
class Node>
class EpetraVectorT;
88 template<
class EpetraGlobalOrdinal,
class Node>
90 :
public virtual MultiVector<double, int, EpetraGlobalOrdinal, Node>
104 "Xpetra::EpetraMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
110 "Xpetra::EpetraMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
116 "Xpetra::EpetraMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
149 return Teuchos::null;
154 return Teuchos::null;
191 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) { }
206 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) { }
242 void randomize(
bool bUseXpetraImplementation =
false) { }
273 "Xpetra::EpetraMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
282 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
296 template<
class TargetDeviceType>
297 typename Kokkos::Impl::if_c<
298 Kokkos::Impl::is_same<
299 typename dual_view_type::t_dev_um::execution_space::memory_space,
300 typename TargetDeviceType::memory_space>::value,
301 typename dual_view_type::t_dev_um,
302 typename dual_view_type::t_host_um>::type
303 getLocalView ()
const {
304 typename Kokkos::Impl::if_c<
305 Kokkos::Impl::is_same<
306 typename dual_view_type::t_dev_um::execution_space::memory_space,
307 typename TargetDeviceType::memory_space>::value,
308 typename dual_view_type::t_dev_um,
309 typename dual_view_type::t_host_um>::type dummy;
313 typename dual_view_type::t_host_um getHostLocalView ()
const {
314 return typename dual_view_type::t_host_um();
317 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
318 throw std::runtime_error(
"Epetra does not support device views!");
319 #ifndef __NVCC__ //prevent nvcc warning
320 typename dual_view_type::t_dev_um ret;
338 #ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES
341 :
public virtual MultiVector<double, int, int, EpetraNode>
365 const std::string tfecfFuncName(
"MultiVector(ArrayOfPtrs)");
367 ": ArrayOfPtrs.size() must be strictly positive and as large as ArrayOfPtrs.");
369 #ifdef HAVE_XPETRA_DEBUG
372 size_t localLength = map->getNodeNumElements();
373 for(
int j=0; j<ArrayOfPtrs.size(); j++) {
375 ": ArrayOfPtrs[" << j <<
"].size() (== " << ArrayOfPtrs[j].size() <<
376 ") is not equal to getLocalLength() (== " << localLength);
384 for(
int i=0; i<ArrayOfPtrs.size(); i++) {
385 arrayOfRawPtrs[i] = ArrayOfPtrs[i].
getRawPtr();
387 double** rawArrayOfRawPtrs =
const_cast<double**
>(arrayOfRawPtrs.getRawPtr());
430 double ** arrayOfPointers;
432 vec_->ExtractView(&arrayOfPointers);
434 double * data = arrayOfPointers[j];
435 int localLength = vec_->MyLength();
444 double ** arrayOfPointers;
446 vec_->ExtractView(&arrayOfPointers);
448 double * data = arrayOfPointers[j];
449 int localLength = vec_->MyLength();
464 vec_->Dot(*eA.getEpetra_MultiVector(), dots.
getRawPtr());
482 for (
size_t j = 0; j < numVecs; ++j) {
492 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) {
XPETRA_MONITOR(
"EpetraMultiVectorT::update"); vec_->Update(alpha, toEpetra<GlobalOrdinal,Node>(A), beta, toEpetra<GlobalOrdinal,Node>(B), gamma); }
507 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) {
XPETRA_MONITOR(
"EpetraMultiVectorT::multiply"); vec_->Multiply(
toEpetra(transA),
toEpetra(transB), alpha,
toEpetra(A),
toEpetra(B), beta); }
510 void elementWiseMultiply(
Scalar scalarAB,
const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B,
Scalar scalarThis) {
XPETRA_MONITOR(
"EpetraMultiVectorT::elementWiseMultiply"); vec_->Multiply(scalarAB, toEpetra<GlobalOrdinal,Node>(A), toEpetra<GlobalOrdinal,Node>(B), scalarThis); }
529 auto vv = toEpetra<GlobalOrdinal,Node>(vec);
530 return ( (vec_->MyLength() == vv.MyLength()) && (vec_->NumVectors() == vv.NumVectors()));
557 if (bUseXpetraImplementation)
621 if (!map.is_null()) {
652 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
666 template<
class TargetDeviceType>
667 typename Kokkos::Impl::if_c<
668 Kokkos::Impl::is_same<
669 typename dual_view_type::t_dev_um::execution_space::memory_space,
670 typename TargetDeviceType::memory_space>::value,
671 typename dual_view_type::t_dev_um,
672 typename dual_view_type::t_host_um>::type
673 getLocalView ()
const {
677 typename dual_view_type::t_host_um getHostLocalView ()
const {
678 typedef Kokkos::View<
typename dual_view_type::t_host::data_type ,
680 typename dual_view_type::t_host::device_type ,
681 Kokkos::MemoryUnmanaged> epetra_view_type;
686 vec_->ExtractView(&data, &myLDA);
687 int localLength = vec_->MyLength();
691 epetra_view_type test = epetra_view_type(data, localLength, numVectors);
692 typename dual_view_type::t_host_um ret = subview(test, Kokkos::ALL(), Kokkos::ALL());
697 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
698 throw std::runtime_error(
"Epetra does not support device views!");
699 #ifndef __NVCC__ //prevent nvcc warning
700 typename dual_view_type::t_dev_um ret;
715 const this_type* rhsPtr =
dynamic_cast<const this_type*
> (&rhs);
717 rhsPtr == NULL, std::invalid_argument,
"Xpetra::MultiVector::operator=: "
718 "The left-hand side (LHS) of the assignment has a different type than "
719 "the right-hand side (RHS). The LHS has type Xpetra::EpetraMultiVectorT "
720 "(which means it wraps an Epetra_MultiVector), but the RHS has some "
721 "other type. This probably means that the RHS wraps a Tpetra::Multi"
722 "Vector. Xpetra::MultiVector does not currently implement assignment "
723 "from a Tpetra object to an Epetra object, though this could be added "
724 "with sufficient interest.");
730 rhsImpl.
is_null (), std::logic_error,
"Xpetra::MultiVector::operator= "
731 "(in Xpetra::EpetraMultiVectorT::assign): *this (the right-hand side of "
732 "the assignment) has a null RCP<Epetra_MultiVector> inside. Please "
733 "report this bug to the Xpetra developers.");
735 lhsImpl.
is_null (), std::logic_error,
"Xpetra::MultiVector::operator= "
736 "(in Xpetra::EpetraMultiVectorT::assign): The left-hand side of the "
737 "assignment has a null RCP<Epetra_MultiVector> inside. Please report "
738 "this bug to the Xpetra developers.");
752 #ifndef XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES
755 :
public virtual MultiVector<double, int, long long, EpetraNode>
779 const std::string tfecfFuncName(
"MultiVector(ArrayOfPtrs)");
781 ": ArrayOfPtrs.size() must be strictly positive and as large as ArrayOfPtrs.");
783 #ifdef HAVE_XPETRA_DEBUG
786 size_t localLength = map->getNodeNumElements();
787 for(
int j=0; j<ArrayOfPtrs.size(); j++) {
789 ": ArrayOfPtrs[" << j <<
"].size() (== " << ArrayOfPtrs[j].size() <<
790 ") is not equal to getLocalLength() (== " << localLength);
798 for(
int i=0; i<ArrayOfPtrs.size(); i++) {
799 arrayOfRawPtrs[i] = ArrayOfPtrs[i].
getRawPtr();
801 double** rawArrayOfRawPtrs =
const_cast<double**
>(arrayOfRawPtrs.getRawPtr());
844 double ** arrayOfPointers;
846 vec_->ExtractView(&arrayOfPointers);
848 double * data = arrayOfPointers[j];
849 int localLength = vec_->MyLength();
858 double ** arrayOfPointers;
860 vec_->ExtractView(&arrayOfPointers);
862 double * data = arrayOfPointers[j];
863 int localLength = vec_->MyLength();
878 vec_->Dot(*eA.getEpetra_MultiVector(), dots.
getRawPtr());
896 for (
size_t j = 0; j < numVecs; ++j) {
906 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) {
XPETRA_MONITOR(
"EpetraMultiVectorT::update"); vec_->Update(alpha, toEpetra<GlobalOrdinal,Node>(A), beta, toEpetra<GlobalOrdinal,Node>(B), gamma); }
921 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) {
XPETRA_MONITOR(
"EpetraMultiVectorT::multiply"); vec_->Multiply(
toEpetra(transA),
toEpetra(transB), alpha,
toEpetra(A),
toEpetra(B), beta); }
924 void elementWiseMultiply(
Scalar scalarAB,
const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B,
Scalar scalarThis) {
XPETRA_MONITOR(
"EpetraMultiVectorT::elementWiseMultiply"); vec_->Multiply(scalarAB, toEpetra<GlobalOrdinal,Node>(A), toEpetra<GlobalOrdinal,Node>(B), scalarThis); }
943 auto vv = toEpetra<GlobalOrdinal,Node>(vec);
944 return ( (vec_->MyLength() == vv.MyLength()) && (vec_->NumVectors() == vv.NumVectors()));
971 if (bUseXpetraImplementation)
1035 if (!map.is_null()) {
1063 vec_->SetSeed(seed);
1066 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
1080 template<
class TargetDeviceType>
1081 typename Kokkos::Impl::if_c<
1082 Kokkos::Impl::is_same<
1083 typename dual_view_type::t_dev_um::execution_space::memory_space,
1084 typename TargetDeviceType::memory_space>::value,
1085 typename dual_view_type::t_dev_um,
1086 typename dual_view_type::t_host_um>::type
1087 getLocalView ()
const {
1091 typename dual_view_type::t_host_um getHostLocalView ()
const {
1092 typedef Kokkos::View<
typename dual_view_type::t_host::data_type ,
1094 typename dual_view_type::t_host::device_type ,
1095 Kokkos::MemoryUnmanaged> epetra_view_type;
1098 double* data = NULL;
1100 vec_->ExtractView(&data, &myLDA);
1101 int localLength = vec_->MyLength();
1105 epetra_view_type test = epetra_view_type(data, localLength, numVectors);
1106 typename dual_view_type::t_host_um ret = subview(test, Kokkos::ALL(), Kokkos::ALL());
1111 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
1112 throw std::runtime_error(
"Epetra does not support device views!");
1113 #ifndef __NVCC__ //prevent nvcc warning
1114 typename dual_view_type::t_dev_um ret;
1129 const this_type* rhsPtr =
dynamic_cast<const this_type*
> (&rhs);
1131 rhsPtr == NULL, std::invalid_argument,
"Xpetra::MultiVector::operator=: "
1132 "The left-hand side (LHS) of the assignment has a different type than "
1133 "the right-hand side (RHS). The LHS has type Xpetra::EpetraMultiVectorT "
1134 "(which means it wraps an Epetra_MultiVector), but the RHS has some "
1135 "other type. This probably means that the RHS wraps a Tpetra::Multi"
1136 "Vector. Xpetra::MultiVector does not currently implement assignment "
1137 "from a Tpetra object to an Epetra object, though this could be added "
1138 "with sufficient interest.");
1144 rhsImpl.
is_null (), std::logic_error,
"Xpetra::MultiVector::operator= "
1145 "(in Xpetra::EpetraMultiVectorT::assign): *this (the right-hand side of "
1146 "the assignment) has a null RCP<Epetra_MultiVector> inside. Please "
1147 "report this bug to the Xpetra developers.");
1149 lhsImpl.
is_null (), std::logic_error,
"Xpetra::MultiVector::operator= "
1150 "(in Xpetra::EpetraMultiVectorT::assign): The left-hand side of the "
1151 "assignment has a null RCP<Epetra_MultiVector> inside. Please report "
1152 "this bug to the Xpetra developers.");
1155 *lhsImpl = *rhsImpl;
1169 #endif // XPETRA_EPETRAMULTIVECTOR_HPP
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 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: this = gamma*this + alpha*A + beta*B.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
Set multi-vector values from array of pointers using Teuchos memory management classes. (copy).
size_t getNumVectors() const
Number of columns in the multivector.
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: this = beta*this + alpha*A.
std::string description() const
A simple one-line description 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.
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 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: this = gamma*this + alpha*A + beta*B.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
Set multi-vector values from array of pointers using Teuchos memory management classes. (copy).
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: this = gamma*this + alpha*A + beta*B.
void scale(const Scalar &alpha)
Scale in place: this = alpha*this.
void scale(const Scalar &alpha)
Scale in place: this = alpha*this.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update: this = beta*this + alpha*A.
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
virtual ~EpetraMultiVectorT()
MultiVector destructor.
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
RCP< Epetra_MultiVector > vec_
The Epetra_MultiVector which this class wraps.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
size_t getLocalLength() const
Local number of rows on the calling process.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
#define TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(throw_exception_test, Exception, msg)
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...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void Xpetra_randomize()
Set multi-vector values to random numbers. XPetra implementation.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic MultiVector constuctor.
EpetraMultiVectorT(const RCP< Epetra_MultiVector > &vec)
EpetraMultiVectorT constructor to wrap a Epetra_MultiVector object.
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 doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
Exception throws to report errors in the internal logical of the program.
const Epetra_CrsGraph & toEpetra(const RCP< const CrsGraph< int, GlobalOrdinal, Node > > &graph)
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
Checks to see if the local length, number of vectors and size of Scalar type match.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic MultiVector constuctor.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
RCP< Epetra_MultiVector > getEpetra_MultiVector() const
Get the underlying Epetra multivector.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Multiply a Vector A elementwise by a MultiVector B.
global_size_t getGlobalLength() const
Global number of rows in the multivector.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
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 replaceMap(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
Replace the underlying Map in place.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
EpetraMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Multiply a Vector A elementwise by a MultiVector B.
EpetraGlobalOrdinal GlobalOrdinal
RCP< Epetra_MultiVector > getEpetra_MultiVector() const
Get the underlying Epetra multivector.
global_size_t getGlobalLength() const
Global number of rows in the multivector.
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
std::string description() const
A simple one-line description of this object.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
size_t getNumVectors() const
Number of columns in the multivector.
void setSeed(unsigned int seed)
Set seed for Random function.
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
EpetraMultiVectorT(const RCP< Epetra_MultiVector > &vec)
EpetraMultiVectorT constructor to wrap a Epetra_MultiVector object.
virtual ~EpetraMultiVectorT()
MultiVector destructor.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
Set multi-vector values from array of pointers using Teuchos memory management classes. (copy).
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic MultiVector constuctor.
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 doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
size_t getNumVectors() const
Number of columns in the multivector.
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
void replaceMap(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
Replace the underlying Map in place.
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).
global_size_t getGlobalLength() const
Global number of rows in the multivector.
Exception throws when you call an unimplemented method of Xpetra.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
static void seedrandom(unsigned int s)
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t global_size_t
Global size_t object.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update: this = beta*this + alpha*A.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
size_t getLocalLength() const
Local number of rows on the calling process.
static const EVerbosityLevel verbLevel_default
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
Checks to see if the local length, number of vectors and size of Scalar type match.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void setSeed(unsigned int seed)
Set seed for Random function.
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).
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
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 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).
EpetraMultiVectorT(const RCP< Epetra_MultiVector > &vec)
EpetraMultiVectorT constructor to wrap a Epetra_MultiVector object.
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
void setSeed(unsigned int seed)
Set seed for Random function.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void scale(const Scalar &alpha)
Scale in place: this = alpha*this.
TypeTo as(const TypeFrom &t)
RCP< Epetra_MultiVector > vec_
The Epetra_MultiVector which this class wraps.
virtual ~EpetraMultiVectorT()
MultiVector destructor.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
void replaceMap(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
Replace the underlying Map in place.
EpetraMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
RCP< Epetra_MultiVector > getEpetra_MultiVector() const
Get the underlying Epetra multivector.
CombineMode
Xpetra::Combine Mode enumerable type.
#define XPETRA_MONITOR(funcName)
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Multiply a Vector A elementwise by a MultiVector B.
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).
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
EpetraMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
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).
std::string description() const
A simple one-line description of this object.
Teuchos::RCP< const Vector< double, int, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
Teuchos::RCP< Vector< double, int, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
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 norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const