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> 75 template<
class GlobalOrdinal>
76 const Epetra_MultiVector &
toEpetra(
const MultiVector<double,int,GlobalOrdinal> &);
77 template<
class GlobalOrdinal>
78 Epetra_MultiVector &
toEpetra(MultiVector<double, int,GlobalOrdinal> &);
79 template<
class GlobalOrdinal>
80 RCP<MultiVector<double, int, GlobalOrdinal> >
toXpetra(RCP<Epetra_MultiVector> vec);
87 template<
class EpetraGlobalOrdinal>
89 :
public virtual MultiVector<double, int, EpetraGlobalOrdinal>
121 void replaceGlobalValue(GlobalOrdinal globalRow,
size_t vectorIndex,
const Scalar &value) {
XPETRA_MONITOR(
"EpetraMultiVectorT::replaceGlobalValue");
vec_->ReplaceGlobalValue(globalRow, Teuchos::as<int>(vectorIndex), value); }
124 void sumIntoGlobalValue(GlobalOrdinal globalRow,
size_t vectorIndex,
const Scalar &value) {
XPETRA_MONITOR(
"EpetraMultiVectorT::sumIntoGlobalValue");
vec_->SumIntoGlobalValue(globalRow, Teuchos::as<int>(vectorIndex), value); }
127 void replaceLocalValue(LocalOrdinal myRow,
size_t vectorIndex,
const Scalar &value) {
XPETRA_MONITOR(
"EpetraMultiVectorT::replaceLocalValue");
vec_->ReplaceMyValue(myRow, Teuchos::as<int>(vectorIndex), value); }
130 void sumIntoLocalValue(LocalOrdinal myRow,
size_t vectorIndex,
const Scalar &value) {
XPETRA_MONITOR(
"EpetraMultiVectorT::sumIntoLocalValue");
vec_->SumIntoMyValue(myRow, Teuchos::as<int>(vectorIndex), value); }
175 for (
size_t j = 0; j < numVecs; ++j) {
176 vec_->Scale (alpha[j]);
184 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(A), beta,
toEpetra(B), gamma); }
202 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); }
205 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(A),
toEpetra(B), scalarThis); }
238 if (bUseXpetraImplementation)
284 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 298 template<
class TargetDeviceType>
299 typename Kokkos::Impl::if_c<
300 Kokkos::Impl::is_same<
301 typename dual_view_type::t_dev_um::execution_space::memory_space,
302 typename TargetDeviceType::memory_space>::value,
303 typename dual_view_type::t_dev_um,
304 typename dual_view_type::t_host_um>::type
305 getLocalView ()
const {
309 typename dual_view_type::t_host_um getHostLocalView ()
const {
310 typedef Kokkos::View<
typename dual_view_type::t_host::data_type ,
312 typename dual_view_type::t_host::device_type ,
313 Kokkos::MemoryUnmanaged> epetra_view_type;
318 vec_->ExtractView(&data, &myLDA);
319 int localLength =
vec_->MyLength();
323 epetra_view_type test = epetra_view_type(data, localLength, numVectors);
324 typename dual_view_type::t_host_um ret = subview(test, Kokkos::ALL(), Kokkos::ALL());
329 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
330 throw std::runtime_error(
"Epetra does not support device views!");
331 typename dual_view_type::t_dev_um ret;
351 #ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES 355 #ifndef XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES 364 template<
class EpetraGlobalOrdinal>
368 const std::string tfecfFuncName(
"MultiVector(ArrayOfPtrs)");
370 ": ArrayOfPtrs.size() must be strictly positive and as large as ArrayOfPtrs.");
372 #ifdef HAVE_XPETRA_DEBUG 375 size_t localLength = map->getNodeNumElements();
376 for(
int j=0; j<ArrayOfPtrs.size(); j++) {
378 ": ArrayOfPtrs[" << j <<
"].size() (== " << ArrayOfPtrs[j].size() <<
379 ") is not equal to getLocalLength() (== " << localLength);
387 for(
int i=0; i<ArrayOfPtrs.size(); i++) {
388 arrayOfRawPtrs[i] = ArrayOfPtrs[i].
getRawPtr();
390 double** rawArrayOfRawPtrs =
const_cast<double**
>(arrayOfRawPtrs.getRawPtr());
396 template<
class EpetraGlobalOrdinal>
402 template<
class EpetraGlobalOrdinal>
408 template<
class EpetraGlobalOrdinal>
412 double ** arrayOfPointers;
414 vec_->ExtractView(&arrayOfPointers);
416 double * data = arrayOfPointers[j];
417 int localLength =
vec_->MyLength();
422 template<
class EpetraGlobalOrdinal>
426 double ** arrayOfPointers;
428 vec_->ExtractView(&arrayOfPointers);
430 double * data = arrayOfPointers[j];
431 int localLength =
vec_->MyLength();
436 template<
class EpetraGlobalOrdinal>
444 template<
class EpetraGlobalOrdinal>
447 template<
class EpetraGlobalOrdinal>
450 template<
class EpetraGlobalOrdinal>
453 template<
class EpetraGlobalOrdinal>
458 vec_->NormWeighted(*eWeights.getEpetra_MultiVector(), norms.getRawPtr());
461 template<
class EpetraGlobalOrdinal>
464 template<
class EpetraGlobalOrdinal>
471 template<
class EpetraGlobalOrdinal>
477 template<
class EpetraGlobalOrdinal>
489 template<
class EpetraGlobalOrdinal>
501 template<
class EpetraGlobalOrdinal>
513 template<
class EpetraGlobalOrdinal>
525 template<
class EpetraGlobalOrdinal>
528 if (!map.is_null()) {
533 Epetra_SerialComm SComm;
534 Epetra_Map NewMap((EpetraGlobalOrdinal)
vec_->MyLength(), (EpetraGlobalOrdinal)
vec_->Map().IndexBase64(), SComm);
540 template<
class EpetraGlobalOrdinal>
545 const this_type* rhsPtr =
dynamic_cast<const this_type*
> (&rhs);
547 rhsPtr == NULL, std::invalid_argument,
"Xpetra::MultiVector::operator=: " 548 "The left-hand side (LHS) of the assignment has a different type than " 549 "the right-hand side (RHS). The LHS has type Xpetra::EpetraMultiVectorT " 550 "(which means it wraps an Epetra_MultiVector), but the RHS has some " 551 "other type. This probably means that the RHS wraps a Tpetra::Multi" 552 "Vector. Xpetra::MultiVector does not currently implement assignment " 553 "from a Tpetra object to an Epetra object, though this could be added " 554 "with sufficient interest.");
560 rhsImpl.
is_null (), std::logic_error,
"Xpetra::MultiVector::operator= " 561 "(in Xpetra::EpetraMultiVectorT::assign): *this (the right-hand side of " 562 "the assignment) has a null RCP<Epetra_MultiVector> inside. Please " 563 "report this bug to the Xpetra developers.");
565 lhsImpl.
is_null (), std::logic_error,
"Xpetra::MultiVector::operator= " 566 "(in Xpetra::EpetraMultiVectorT::assign): The left-hand side of the " 567 "assignment has a null RCP<Epetra_MultiVector> inside. Please report " 568 "this bug to the Xpetra developers.");
579 #endif // XPETRA_EPETRAMULTIVECTOR_HPP
size_t getLocalLength() const
Local number of rows on the calling process.
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...
const Epetra_CrsGraph & toEpetra(const RCP< const CrsGraph< int, GlobalOrdinal > > &graph)
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.
EpetraMultiVectorT< int > EpetraMultiVector
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).
#define TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(throw_exception_test, Exception, msg)
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.
RCP< Epetra_MultiVector > vec_
The Epetra_MultiVector which this class wraps.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void Xpetra_randomize()
Set multi-vector values to random numbers. XPetra implementation.
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
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< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
RCP< const CrsGraph< int, GlobalOrdinal > > toXpetra(const Epetra_CrsGraph &g)
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
MultiVector< double, int, GlobalOrdinal >::node_type Node
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update: this = beta*this + alpha*A.
void replaceMap(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
Replace the underlying Map in place.
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 doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
Exception throws when you call an unimplemented method of Xpetra.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
static void seedrandom(unsigned int s)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t global_size_t
Global size_t object.
static const EVerbosityLevel verbLevel_default
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
std::string description() const
A simple one-line description of this object.
EpetraMultiVectorT(const RCP< Epetra_MultiVector > &vec)
EpetraMultiVectorT constructor to wrap a Epetra_MultiVector object.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
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.
TypeTo as(const TypeFrom &t)
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
Teuchos::RCP< Vector< double, int, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
CombineMode
Xpetra::Combine Mode enumerable type.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
#define XPETRA_MONITOR(funcName)
EpetraMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
EpetraMultiVectorT< long long > EpetraMultiVector64
void scale(const Scalar &alpha)
Scale in place: this = alpha*this.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic MultiVector constuctor.
size_t getNumVectors() const
Number of columns in the multivector.
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)
Multiply a Vector A elementwise by a MultiVector B.
EpetraGlobalOrdinal GlobalOrdinal
RCP< Epetra_MultiVector > getEpetra_MultiVector() const
Get the underlying Epetra multivector.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void normWeighted(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &weights, const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
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.
virtual ~EpetraMultiVectorT()
MultiVector destructor.