10#ifndef XPETRA_EPETRAINTMULTIVECTOR_HPP
11#define XPETRA_EPETRAINTMULTIVECTOR_HPP
16#include "Xpetra_MultiVector.hpp"
21#include "Epetra_IntMultiVector.h"
23#if defined(XPETRA_ENABLE_DEPRECATED_CODE)
25#if defined(Xpetra_SHOW_DEPRECATED_WARNINGS)
26#warning "The header file Trilinos/packages/xpetra/src/MultiVector/Xpetra_EpetraIntMultiVector.hpp is deprecated."
30#error "The header file Trilinos/packages/xpetra/src/MultiVector/Xpetra_EpetraIntMultiVector.hpp is deprecated."
36template <
class GlobalOrdinal,
class Node>
39template <
class GlobalOrdinal,
class Node>
44template <
class EpetraGlobalOrdinal,
class Node>
46 :
public MultiVector<int, int, EpetraGlobalOrdinal, Node> {
58 "Xpetra::EpetraIntMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
64 "Xpetra::EpetraIntMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
70 "Xpetra::EpetraIntMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
85 void randomize(
bool bUseXpetraImplementation =
true) {
101 "Xpetra::EpetraIntMultiVectorT::setSeed(): Functionnality not available in Epetra");
111 return Teuchos::null;
116 return Teuchos::null;
139 "This function is not implemented in Epetra_IntMultiVector");
148 "This function is not implemented in Epetra_IntMultiVector");
158 "Xpetra::EpetraIntMultiVectorT::scale(): Functionnality not available in Epetra");
165 "Xpetra::EpetraIntMultiVectorT::update(): Functionnality not available in Epetra");
169 void update(
const int &alpha,
const MultiVector<int, int, GlobalOrdinal, Node> &A,
const int &beta,
const MultiVector<int, int, GlobalOrdinal, Node> &B,
const int &gamma) {
172 "Xpetra::EpetraIntMultiVectorT::update(): Functionnality not available in Epetra");
179 "Xpetra::EpetraIntMultiVectorT::norm1(): Functionnality not available in Epetra");
186 "Xpetra::EpetraIntMultiVectorT::norm2(): Functionnality not available in Epetra");
193 "Xpetra::EpetraIntMultiVectorT::normInf(): Functionnality not available in Epetra");
200 "Xpetra::EpetraIntMultiVectorT::meanValue(): Functionnality not available in Epetra");
207 "Xpetra::EpetraIntMultiVectorT::maxValue(): Functionnality not available in Epetra");
211 void multiply(
Teuchos::ETransp transA,
Teuchos::ETransp transB,
const int &alpha,
const MultiVector<int, int, GlobalOrdinal, Node> &A,
const MultiVector<int, int, GlobalOrdinal, Node> &B,
const int &beta) {
214 "Xpetra::EpetraIntMultiVectorT::multiply(): Functionnality not available in Epetra");
221 "Xpetra_EpetraIntMultiVector: elementWiseMultiply not implemented because Epetra_IntMultiVector does not support this operation");
280 return std::string(
"");
296 return Teuchos::null;
328#ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES
331 :
public virtual MultiVector<int, int, int, EpetraNode> {
362 "Xpetra::EpetraIntMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
376 ierr = vec_->PutScalar(value);
410 typedef Kokkos::View<
typename dual_view_type::t_host::data_type,
412 typename dual_view_type::t_host::device_type,
413 Kokkos::MemoryUnmanaged>
418 vec_->ExtractView(&data, &myLDA);
419 int localLength = vec_->MyLength();
423 epetra_view_type test = epetra_view_type(data, localLength, numVectors);
424 typename dual_view_type::t_host_um ret = subview(test, Kokkos::ALL(), Kokkos::ALL());
451 int **arrayOfPointers;
452 vec_->ExtractView(&arrayOfPointers);
453 int *data = arrayOfPointers[j];
454 int localLength = vec_->MyLength();
464 int **arrayOfPointers;
465 vec_->ExtractView(&arrayOfPointers);
466 int *data = arrayOfPointers[j];
467 int localLength = vec_->MyLength();
483 "This function is not implemented in Epetra_IntMultiVector");
490 "This function is not available in Epetra_IntMultiVector");
522 void update(
const int & ,
const MultiVector<int, int, GlobalOrdinal, Node> & ,
const int & ,
const MultiVector<int, int, GlobalOrdinal, Node> & ,
const int & ) {
561 void multiply(
Teuchos::ETransp ,
Teuchos::ETransp ,
const int & ,
const MultiVector<int, int, GlobalOrdinal, Node> & ,
const MultiVector<int, int, GlobalOrdinal, Node> & ,
const int & ) {
579 vec_->ReplaceGlobalValue(globalRow, vectorIndex, value);
584 vec_->SumIntoGlobalValue(globalRow, vectorIndex, value);
589 vec_->ReplaceMyValue(myRow, vectorIndex, value);
594 vec_->SumIntoMyValue(myRow, vectorIndex, value);
604 return vec_->NumVectors();
609 return vec_->MyLength();
618 auto vv = toEpetra<GlobalOrdinal, Node>(vec);
619 return ((
getLocalLength() == Teuchos::as<size_t>(vv.MyLength())) &&
633 std::ostringstream oss;
680 int err = vec_->
Import(v, *tImporter.getEpetra_Import(),
toEpetra(CM));
692 int err = vec_->
Import(v, *tImporter.getEpetra_Import(),
toEpetra(CM));
704 int err = vec_->
Import(v, *tExporter.getEpetra_Export(),
toEpetra(CM));
716 int err = vec_->
Export(v, *tExporter.getEpetra_Export(),
toEpetra(CM));
723 if (!map.is_null()) {
741 const this_type *rhsPtr =
dynamic_cast<const this_type *
>(&rhs);
743 rhsPtr == NULL, std::invalid_argument,
744 "Xpetra::MultiVector::operator=: "
745 "The left-hand side (LHS) of the assignment has a different type than "
746 "the right-hand side (RHS). The LHS has type Xpetra::EpetraIntMultiVectorT "
747 "(which means it wraps an Epetra_IntMultiVector), but the RHS has some "
748 "other type. This probably means that the RHS wraps either an "
749 "Tpetra::MultiVector, or an Epetra_MultiVector. Xpetra::MultiVector "
750 "does not currently implement assignment from a Tpetra object to an "
751 "Epetra object, though this could be added with sufficient interest.");
757 rhsImpl.
is_null(), std::logic_error,
758 "Xpetra::MultiVector::operator= "
759 "(in Xpetra::EpetraIntMultiVectorT::assign): *this (the right-hand side of "
760 "the assignment) has a null RCP<Epetra_IntMultiVector> inside. Please "
761 "report this bug to the Xpetra developers.");
763 lhsImpl.
is_null(), std::logic_error,
764 "Xpetra::MultiVector::operator= "
765 "(in Xpetra::EpetraIntMultiVectorT::assign): The left-hand side of the "
766 "assignment has a null RCP<Epetra_IntMultiVector> inside. Please report "
767 "this bug to the Xpetra developers.");
780#ifndef XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES
783 :
public virtual MultiVector<int, int, long long, EpetraNode> {
806 "Xpetra::EpetraIntMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
820 ierr = vec_->PutScalar(value);
863 int **arrayOfPointers;
864 vec_->ExtractView(&arrayOfPointers);
865 int *data = arrayOfPointers[j];
866 int localLength = vec_->MyLength();
876 int **arrayOfPointers;
877 vec_->ExtractView(&arrayOfPointers);
878 int *data = arrayOfPointers[j];
879 int localLength = vec_->MyLength();
901 "This function is not available in Epetra_IntMultiVector");
908 "This function is not implemented in Epetra_IntMultiVector");
932 void update(
const int & ,
const MultiVector<int, int, GlobalOrdinal, Node> & ,
const int & ,
const MultiVector<int, int, GlobalOrdinal, Node> & ,
const int & ) {
971 void multiply(
Teuchos::ETransp ,
Teuchos::ETransp ,
const int & ,
const MultiVector<int, int, GlobalOrdinal, Node> & ,
const MultiVector<int, int, GlobalOrdinal, Node> & ,
const int & ) {
989 vec_->ReplaceGlobalValue(globalRow, vectorIndex, value);
994 vec_->SumIntoGlobalValue(globalRow, vectorIndex, value);
999 vec_->ReplaceMyValue(myRow, vectorIndex, value);
1004 vec_->SumIntoMyValue(myRow, vectorIndex, value);
1014 return vec_->NumVectors();
1026 auto vv = toEpetra<GlobalOrdinal, Node>(vec);
1027 return ((
getLocalLength() == Teuchos::as<size_t>(vv.MyLength())) &&
1040 std::ostringstream oss;
1087 int err = vec_->
Import(v, *tImporter.getEpetra_Import(),
toEpetra(CM));
1099 int err = vec_->
Import(v, *tImporter.getEpetra_Import(),
toEpetra(CM));
1111 int err = vec_->
Import(v, *tExporter.getEpetra_Export(),
toEpetra(CM));
1123 int err = vec_->
Export(v, *tExporter.getEpetra_Export(),
toEpetra(CM));
1130 if (!map.is_null()) {
1148 const this_type *rhsPtr =
dynamic_cast<const this_type *
>(&rhs);
1150 rhsPtr == NULL, std::invalid_argument,
1151 "Xpetra::MultiVector::operator=: "
1152 "The left-hand side (LHS) of the assignment has a different type than "
1153 "the right-hand side (RHS). The LHS has type Xpetra::EpetraIntMultiVectorT "
1154 "(which means it wraps an Epetra_IntMultiVector), but the RHS has some "
1155 "other type. This probably means that the RHS wraps either an "
1156 "Tpetra::MultiVector, or an Epetra_MultiVector. Xpetra::MultiVector "
1157 "does not currently implement assignment from a Tpetra object to an "
1158 "Epetra object, though this could be added with sufficient interest.");
1164 rhsImpl.
is_null(), std::logic_error,
1165 "Xpetra::MultiVector::operator= "
1166 "(in Xpetra::EpetraIntMultiVectorT::assign): *this (the right-hand side of "
1167 "the assignment) has a null RCP<Epetra_IntMultiVector> inside. Please "
1168 "report this bug to the Xpetra developers.");
1170 lhsImpl.
is_null(), std::logic_error,
1171 "Xpetra::MultiVector::operator= "
1172 "(in Xpetra::EpetraIntMultiVectorT::assign): The left-hand side of the "
1173 "assignment has a null RCP<Epetra_IntMultiVector> inside. Please report "
1174 "this bug to the Xpetra developers.");
1177 *lhsImpl = *rhsImpl;
#define XPETRA_MONITOR(funcName)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
int Import(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
int Export(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
static const EVerbosityLevel verbLevel_default
virtual std::string description() const
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t) const
Return a Vector which is a const view of column j.
void update(const int &, const MultiVector< int, int, GlobalOrdinal, Node > &, const int &, const MultiVector< int, int, GlobalOrdinal, Node > &, const int &)
Update multi-vector with scaled values of A and B, this = gamma*this + alpha*A + beta*B.
EpetraIntMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &, size_t)
Set multi-vector values from array of pointers using Teuchos memory management classes....
void reciprocal(const MultiVector< int, int, GlobalOrdinal, Node > &)
Puts element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,...
std::string description() const
Return a simple one-line description of this object.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute 1-norm of each vector in multi-vector.
Teuchos::ArrayRCP< int > getDataNonConst(size_t j)
Teuchos::RCP< const Map< int, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
dual_view_type::t_host_um getLocalViewHost(Access::OverwriteAllStruct) const override
~EpetraIntMultiVectorT()
Destructor.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t)
Return a Vector which is a nonconst view of column j.
dual_view_type::t_host_um getLocalViewHost(Access::ReadWriteStruct) const override
void randomize(bool=true)
Set multi-vector values to random numbers.
RCP< Epetra_IntMultiVector > vec_
The Epetra_IntMultiVector which this class wraps.
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< int, GlobalOrdinal, Node > &importer, CombineMode CM)
dual_view_type::t_dev_const_um getLocalViewDevice(Access::ReadOnlyStruct) const override
void doImport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< int, GlobalOrdinal, Node > &exporter, CombineMode CM)
void multiply(Teuchos::ETransp, Teuchos::ETransp, const int &, const MultiVector< int, int, GlobalOrdinal, Node > &, const MultiVector< int, int, GlobalOrdinal, Node > &, const int &)
Matrix-Matrix multiplication, this = beta*this + alpha*op(A)*op(B).
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< int, GlobalOrdinal, Node > &exporter, CombineMode CM)
dual_view_type::t_dev_um getLocalViewDevice(Access::OverwriteAllStruct) const override
void update(const int &, const MultiVector< int, int, GlobalOrdinal, Node > &, const int &)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute Inf-norm of each vector in multi-vector.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute 2-norm of each vector in multi-vector.
void scale(Teuchos::ArrayView< const int >)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
size_t getNumVectors() const
Returns the number of vectors in the multi-vector.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
void elementWiseMultiply(int, const Vector< int, int, GlobalOrdinal, Node > &, const MultiVector< int, int, GlobalOrdinal, Node > &, int)
Element-wise multiply of a Vector A with a EpetraMultiVector B.
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dual_view_type dual_view_type
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void randomize(const Scalar &minVal, const Scalar &maxVal, bool bUseXpetraImplementation=true)
Set multi-vector values to random numbers.
void setSeed(unsigned int)
Set seed for Random function.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
void doImport(const DistObject< int, int, GlobalOrdinal, Node > &source, const Import< int, GlobalOrdinal, Node > &importer, CombineMode CM)
EpetraIntMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Sets all vector entries to zero.
dual_view_type::t_dev_um getLocalViewDevice(Access::ReadWriteStruct) const override
RCP< Epetra_IntMultiVector > getEpetra_IntMultiVector() const
void scale(const int &)
Scale the current values of a multi-vector, this = alpha*this.
void replaceMap(const RCP< const Map< int, GlobalOrdinal, Node > > &map)
void maxValue(const Teuchos::ArrayView< int > &) const
Compute max value of each vector in multi-vector.
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.
dual_view_type::t_host_const_um getLocalViewHost(Access::ReadOnlyStruct) const override
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
EpetraIntMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Teuchos::DataAccess copyOrView=Teuchos::Copy)
MultiVector copy constructor.
const RCP< const Comm< int > > getComm() const
global_size_t getGlobalLength() const
Returns the global vector length of vectors in the multi-vector.
void meanValue(const Teuchos::ArrayView< int > &) const
Compute mean (average) value of each vector in multi-vector.
void dot(const MultiVector< int, int, GlobalOrdinal, Node > &, const Teuchos::ArrayView< int > &) const
Computes dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i])
void putScalar(const int &value)
Initialize all values in a multi-vector with specified value.
void abs(const MultiVector< int, int, GlobalOrdinal, Node > &)
Puts element-wise absolute values of input Multi-vector in target: A = abs(this)
Teuchos::ArrayRCP< const int > getData(size_t j) const
size_t getLocalLength() const
Returns the local vector length on the calling processor of vectors in the multi-vector.
void randomize(const Scalar &minVal, const Scalar &maxVal, bool bUseXpetraImplementation=true)
Set multi-vector values to random numbers.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
const RCP< const Comm< int > > getComm() const
void dot(const MultiVector< int, int, GlobalOrdinal, Node > &, const Teuchos::ArrayView< int > &) const
Computes dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i])
void maxValue(const Teuchos::ArrayView< int > &) const
Compute max value of each vector in multi-vector.
global_size_t getGlobalLength() const
Returns the global vector length of vectors in the multi-vector.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
void putScalar(const int &value)
Initialize all values in a multi-vector with specified value.
void reciprocal(const MultiVector< int, int, GlobalOrdinal, Node > &)
Puts element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,...
RCP< Epetra_IntMultiVector > getEpetra_IntMultiVector() const
void doImport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< int, GlobalOrdinal, Node > &exporter, CombineMode CM)
size_t getLocalLength() const
Returns the local vector length on the calling processor of vectors in the multi-vector.
Teuchos::ArrayRCP< const int > getData(size_t j) const
void scale(const int &)
Scale the current values of a multi-vector, this = alpha*this.
void elementWiseMultiply(int, const Vector< int, int, GlobalOrdinal, Node > &, const MultiVector< int, int, GlobalOrdinal, Node > &, int)
Element-wise multiply of a Vector A with a EpetraMultiVector B.
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 norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute 2-norm of each vector in multi-vector.
void abs(const MultiVector< int, int, GlobalOrdinal, Node > &)
Puts element-wise absolute values of input Multi-vector in target: A = abs(this)
void replaceMap(const RCP< const Map< int, GlobalOrdinal, Node > > &map)
void doImport(const DistObject< int, int, GlobalOrdinal, Node > &source, const Import< int, GlobalOrdinal, Node > &importer, CombineMode CM)
Teuchos::RCP< const Map< int, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
void scale(Teuchos::ArrayView< const int >)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
EpetraIntMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Sets all vector entries to zero.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
~EpetraIntMultiVectorT()
Destructor.
EpetraIntMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &, size_t)
Set multi-vector values from array of pointers using Teuchos memory management classes....
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
void meanValue(const Teuchos::ArrayView< int > &) const
Compute mean (average) value of each vector in multi-vector.
void setSeed(unsigned int)
Set seed for Random function.
Teuchos::ArrayRCP< int > getDataNonConst(size_t j)
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< int, GlobalOrdinal, Node > &importer, CombineMode CM)
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void update(const int &, const MultiVector< int, int, GlobalOrdinal, Node > &, const int &, const MultiVector< int, int, GlobalOrdinal, Node > &, const int &)
Update multi-vector with scaled values of A and B, this = gamma*this + alpha*A + beta*B.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute 1-norm of each vector in multi-vector.
std::string description() const
Return a simple one-line description of this object.
size_t getNumVectors() const
Returns the number of vectors in the multi-vector.
void randomize(bool=true)
Set multi-vector values to random numbers.
void update(const int &, const MultiVector< int, int, GlobalOrdinal, Node > &, const int &)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
RCP< Epetra_IntMultiVector > vec_
The Epetra_IntMultiVector which this class wraps.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t)
Return a Vector which is a nonconst view of column j.
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute Inf-norm of each vector in multi-vector.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
EpetraIntMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t) const
Return a Vector which is a const view of column j.
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< int, GlobalOrdinal, Node > &exporter, CombineMode CM)
void multiply(Teuchos::ETransp, Teuchos::ETransp, const int &, const MultiVector< int, int, GlobalOrdinal, Node > &, const MultiVector< int, int, GlobalOrdinal, Node > &, const int &)
Matrix-Matrix multiplication, this = beta*this + alpha*op(A)*op(B).
std::string description() const
Return a simple one-line description of this object.
void setSeed(unsigned int seed)
Set seed for Random function.
void update(const int &alpha, const MultiVector< int, int, GlobalOrdinal, Node > &A, const int &beta)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
size_t getLocalLength() const
Returns the local vector length on the calling processor of vectors in the multi-vector.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
void meanValue(const Teuchos::ArrayView< int > &means) const
Compute mean (average) value of each vector in multi-vector.
Teuchos::RCP< const Map< int, 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.
size_t getNumVectors() const
Returns the number of vectors in the multi-vector.
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< int, GlobalOrdinal, Node > &importer, CombineMode CM)
EpetraIntMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Teuchos::DataAccess copyOrView=Teuchos::Copy)
MultiVector copy constructor.
void abs(const MultiVector< int, int, GlobalOrdinal, Node > &A)
Puts element-wise absolute values of input Multi-vector in target: A = abs(this)
void randomize(const Scalar &minVal, const Scalar &maxVal, bool bUseXpetraImplementation=true)
Set multi-vector values to random numbers.
EpetraIntMultiVectorT(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....
Teuchos::ArrayRCP< int > getDataNonConst(size_t j)
void putScalar(const int &value)
Initialize all values in a multi-vector with specified value.
void maxValue(const Teuchos::ArrayView< int > &maxs) const
Compute max value of each vector in multi-vector.
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
EpetraIntMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Sets all vector entries to zero.
void replaceMap(const RCP< const Map< int, GlobalOrdinal, Node > > &map)
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &norms) const
Compute 2-norm of each vector in multi-vector.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
void reciprocal(const MultiVector< int, int, GlobalOrdinal, Node > &A)
Puts element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,...
void multiply(Teuchos::ETransp transA, Teuchos::ETransp transB, const int &alpha, const MultiVector< int, int, GlobalOrdinal, Node > &A, const MultiVector< int, int, GlobalOrdinal, Node > &B, const int &beta)
Matrix-Matrix multiplication, this = beta*this + alpha*op(A)*op(B).
void doImport(const DistObject< int, int, GlobalOrdinal, Node > &source, const Import< int, GlobalOrdinal, Node > &importer, CombineMode CM)
void elementWiseMultiply(int scalarAB, const Vector< int, int, GlobalOrdinal, Node > &A, const MultiVector< int, int, GlobalOrdinal, Node > &B, int scalarThis)
Element-wise multiply of a Vector A with a EpetraMultiVector B.
const RCP< const Comm< int > > getComm() const
~EpetraIntMultiVectorT()
Destructor.
void dot(const MultiVector< int, int, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< int > &dots) const
Computes dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i])
EpetraGlobalOrdinal GlobalOrdinal
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< int, GlobalOrdinal, Node > &exporter, CombineMode CM)
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
Teuchos::ArrayRCP< const int > getData(size_t j) const
global_size_t getGlobalLength() const
Returns the global vector length of vectors in the multi-vector.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void scale(const int &alpha)
Scale the current values of a multi-vector, this = alpha*this.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
void scale(Teuchos::ArrayView< const int > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
RCP< Epetra_IntMultiVector > getEpetra_IntMultiVector() const
void update(const int &alpha, const MultiVector< int, int, GlobalOrdinal, Node > &A, const int &beta, const MultiVector< int, int, GlobalOrdinal, Node > &B, const int &gamma)
Update multi-vector with scaled values of A and B, this = gamma*this + alpha*A + beta*B.
void randomize(bool bUseXpetraImplementation=true)
Set multi-vector values to random numbers.
void doImport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< int, GlobalOrdinal, Node > &exporter, CombineMode CM)
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.
Exception throws when you call an unimplemented method of Xpetra.
Exception throws to report errors in the internal logical of the program.
virtual dual_view_type::t_host_const_um getLocalViewHost(Access::ReadOnlyStruct) const
Kokkos::DualView< impl_scalar_type **, Kokkos::LayoutStride, typename node_type::device_type, Kokkos::MemoryUnmanaged > dual_view_type
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
constexpr struct ReadWriteStruct ReadWrite
Tpetra::KokkosCompat::KokkosSerialWrapperNode EpetraNode
size_t global_size_t
Global size_t object.
const Epetra_CrsGraph & toEpetra(const RCP< const CrsGraph< int, GlobalOrdinal, Node > > &graph)
CombineMode
Xpetra::Combine Mode enumerable type.