46 #ifndef XPETRA_EPETRAINTMULTIVECTOR_HPP
47 #define XPETRA_EPETRAINTMULTIVECTOR_HPP
52 #include "Xpetra_MultiVector.hpp"
57 #include "Epetra_IntMultiVector.h"
62 template<
class GlobalOrdinal,
class Node>
63 Epetra_IntMultiVector &
toEpetra(MultiVector<int, int, GlobalOrdinal, Node> &);
65 template<
class GlobalOrdinal,
class Node>
66 const Epetra_IntMultiVector &
toEpetra(
const MultiVector<int, int, GlobalOrdinal, Node> &);
70 template<
class EpetraGlobalOrdinal,
class Node>
72 :
public MultiVector<int,int,EpetraGlobalOrdinal, Node>
86 "Xpetra::EpetraIntMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
92 "Xpetra::EpetraIntMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
98 "Xpetra::EpetraIntMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
120 "Xpetra::EpetraIntMultiVectorT::setSeed(): Functionnality not available in Epetra"); }
129 Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > >
getVector(
size_t j)
const {
130 return Teuchos::null;
134 Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > >
getVectorNonConst(
size_t j) {
135 return Teuchos::null;
140 Teuchos::ArrayRCP<const int>
getData(
size_t j)
const {
141 return Teuchos::ArrayRCP<const int>();
147 return Teuchos::ArrayRCP<int>();
156 const Teuchos::ArrayView<int> &dots)
const {
158 "This function is not implemented in Epetra_IntMultiVector");
167 "This function is not implemented in Epetra_IntMultiVector");
174 void scale (Teuchos::ArrayView< const int > alpha) {
177 "Xpetra::EpetraIntMultiVectorT::scale(): Functionnality not available in Epetra");
184 "Xpetra::EpetraIntMultiVectorT::update(): Functionnality not available in Epetra");
188 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) {
191 "Xpetra::EpetraIntMultiVectorT::update(): Functionnality not available in Epetra");
195 void norm1(
const Teuchos::ArrayView<Teuchos::ScalarTraits<int>::magnitudeType> &norms)
const {
198 "Xpetra::EpetraIntMultiVectorT::norm1(): Functionnality not available in Epetra");
202 void norm2(
const Teuchos::ArrayView<Teuchos::ScalarTraits<int>::magnitudeType> &norms)
const {
205 "Xpetra::EpetraIntMultiVectorT::norm2(): Functionnality not available in Epetra"); }
208 void normInf(
const Teuchos::ArrayView<Teuchos::ScalarTraits<int>::magnitudeType> &norms)
const {
211 "Xpetra::EpetraIntMultiVectorT::normInf(): Functionnality not available in Epetra"); }
214 void meanValue(
const Teuchos::ArrayView<int> &means)
const {
217 "Xpetra::EpetraIntMultiVectorT::meanValue(): Functionnality not available in Epetra");
221 void maxValue(
const Teuchos::ArrayView<int> &maxs)
const {
224 "Xpetra::EpetraIntMultiVectorT::maxValue(): Functionnality not available in Epetra");
228 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) {
231 "Xpetra::EpetraIntMultiVectorT::multiply(): Functionnality not available in Epetra");
238 "Xpetra_EpetraIntMultiVector: elementWiseMultiply not implemented because Epetra_IntMultiVector does not support this operation");
295 return std::string(
"");
299 void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
const { }
310 Teuchos::RCP<const Map<int, GlobalOrdinal, Node> >
getMap ()
const {
311 return Teuchos::null;
333 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
336 typename dual_view_type::t_host_um getHostLocalView ()
const {
337 throw std::runtime_error(
"EpetraIntVector does not support device views! Must be implemented extra...");
338 #ifndef __NVCC__ //prevent nvcc warning
339 typename dual_view_type::t_host_um ret;
341 TEUCHOS_UNREACHABLE_RETURN(ret);
344 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
345 throw std::runtime_error(
"Epetra does not support device views!");
346 #ifndef __NVCC__ //prevent nvcc warning
347 typename dual_view_type::t_dev_um ret;
349 TEUCHOS_UNREACHABLE_RETURN(ret);
352 template<
class TargetDeviceType>
353 typename Kokkos::Impl::if_c<
355 typename dual_view_type::t_dev_um::execution_space::memory_space,
356 typename TargetDeviceType::memory_space>::value,
357 typename dual_view_type::t_dev_um,
358 typename dual_view_type::t_host_um>::type
359 getLocalView ()
const {
360 return this->MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::template getLocalView<TargetDeviceType>();
381 #ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES
384 :
public virtual MultiVector<int,int,int,EpetraNode>
398 vec_ = rcp(
new Epetra_IntMultiVector(toEpetra<GlobalOrdinal,Node>(map), NumVectors, zeroOut));
403 vec_ = rcp(
new Epetra_IntMultiVector(toEpetra<GlobalOrdinal,Node>(source)));
409 "Xpetra::EpetraIntMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
423 ierr = vec_->PutScalar(value);
442 Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > >
getVector(
size_t )
const {
447 Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > >
getVectorNonConst(
size_t ) {
453 Teuchos::ArrayRCP<const int>
getData(
size_t j)
const {
456 int ** arrayOfPointers;
457 vec_->ExtractView(&arrayOfPointers);
458 int * data = arrayOfPointers[j];
459 int localLength = vec_->MyLength();
461 return ArrayRCP<int>(data, 0, localLength,
false);
469 int ** arrayOfPointers;
470 vec_->ExtractView(&arrayOfPointers);
471 int * data = arrayOfPointers[j];
472 int localLength = vec_->MyLength();
474 return ArrayRCP<int>(data, 0, localLength,
false);
483 const Teuchos::ArrayView<int> &)
const {
488 "This function is not implemented in Epetra_IntMultiVector");
495 "This function is not available in Epetra_IntMultiVector");
513 void scale (Teuchos::ArrayView< const int > ) {
527 void update(
const int &,
const MultiVector<int,int,GlobalOrdinal,Node> &,
const int &,
const MultiVector<int,int,GlobalOrdinal,Node> &,
const int &) {
551 void multiply(Teuchos::ETransp , Teuchos::ETransp ,
const int &,
const MultiVector<int,int,GlobalOrdinal,Node> &,
const MultiVector<int,int,GlobalOrdinal,Node> &,
const int &) {
XPETRA_MONITOR(
"EpetraIntMultiVectorT::multiply"); TEUCHOS_TEST_FOR_EXCEPTION(1,
Xpetra::Exceptions::NotImplemented,
"Not available in Epetra"); }
556 TEUCHOS_TEST_FOR_EXCEPTION(1,
Xpetra::Exceptions::NotImplemented,
"Xpetra_EpetraIntMultiVector: elementWiseMultiply not implemented because Epetra_IntMultiVector does not support this operation");
566 vec_->ReplaceGlobalValue(globalRow, vectorIndex, value);
571 vec_->SumIntoGlobalValue(globalRow, vectorIndex, value);
576 vec_->ReplaceMyValue(myRow, vectorIndex, value);
581 vec_->SumIntoMyValue(myRow, vectorIndex, value);
591 return vec_->NumVectors();
597 return vec_->MyLength();
606 auto vv = toEpetra<GlobalOrdinal,Node>(vec);
607 return ( (
getLocalLength() == Teuchos::as<size_t>(vv.MyLength())) &&
621 std::ostringstream oss;
622 oss << Teuchos::Describable::description();
629 void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
const {
635 using Teuchos::VERB_DEFAULT;
636 using Teuchos::VERB_NONE;
637 using Teuchos::VERB_LOW;
638 using Teuchos::VERB_MEDIUM;
639 using Teuchos::VERB_HIGH;
640 using Teuchos::VERB_EXTREME;
642 if (verbLevel > Teuchos::VERB_NONE)
655 Teuchos::RCP<const Map<int, GlobalOrdinal, Node> >
getMap ()
const {
656 RCP<const Epetra_BlockMap> map = rcp(
new Epetra_BlockMap(vec_->Map()));
667 const Epetra_IntMultiVector & v = *tSource.getEpetra_IntMultiVector();
668 int err = vec_->Import(v, *tImporter.getEpetra_Import(),
toEpetra(CM));
669 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
"Catch error code returned by Epetra.");
679 const Epetra_IntMultiVector & v = *tDest.getEpetra_IntMultiVector();
680 int err = vec_->Import(v, *tImporter.getEpetra_Import(),
toEpetra(CM));
681 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
"Catch error code returned by Epetra.");
691 const Epetra_IntMultiVector & v = *tSource.getEpetra_IntMultiVector();
692 int err = vec_->Import(v, *tExporter.getEpetra_Export(),
toEpetra(CM));
693 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
"Catch error code returned by Epetra.");
703 const Epetra_IntMultiVector & v = *tDest.getEpetra_IntMultiVector();
704 int err = vec_->Export(v, *tExporter.getEpetra_Export(),
toEpetra(CM));
705 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
"Catch error code returned by Epetra.");
711 if (!map.is_null()) {
716 Epetra_SerialComm SComm;
720 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
"Catch error code returned by Epetra.");
726 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
729 typename dual_view_type::t_host_um getHostLocalView ()
const {
730 typedef Kokkos::View<
typename dual_view_type::t_host::data_type ,
732 typename dual_view_type::t_host::device_type ,
733 Kokkos::MemoryUnmanaged> epetra_view_type;
738 vec_->ExtractView(&data, &myLDA);
739 int localLength = vec_->MyLength();
740 int numVectors = vec_->NumVectors();
743 epetra_view_type test = epetra_view_type(data, localLength, numVectors);
744 typename dual_view_type::t_host_um ret = subview(test, Kokkos::ALL(), Kokkos::ALL());
749 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
750 throw std::runtime_error(
"Epetra does not support device views!");
751 #ifndef __NVCC__ //prevent nvcc warning
752 typename dual_view_type::t_dev_um ret;
754 TEUCHOS_UNREACHABLE_RETURN(ret);
767 template<
class TargetDeviceType>
768 typename Kokkos::Impl::if_c<
770 typename dual_view_type::t_dev_um::execution_space::memory_space,
771 typename TargetDeviceType::memory_space>::value,
772 typename dual_view_type::t_dev_um,
773 typename dual_view_type::t_host_um>::type
774 getLocalView ()
const {
775 return this->MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::template getLocalView<TargetDeviceType>();
788 const this_type* rhsPtr =
dynamic_cast<const this_type*
> (&rhs);
789 TEUCHOS_TEST_FOR_EXCEPTION(
790 rhsPtr == NULL, std::invalid_argument,
"Xpetra::MultiVector::operator=: "
791 "The left-hand side (LHS) of the assignment has a different type than "
792 "the right-hand side (RHS). The LHS has type Xpetra::EpetraIntMultiVectorT "
793 "(which means it wraps an Epetra_IntMultiVector), but the RHS has some "
794 "other type. This probably means that the RHS wraps either an "
795 "Tpetra::MultiVector, or an Epetra_MultiVector. Xpetra::MultiVector "
796 "does not currently implement assignment from a Tpetra object to an "
797 "Epetra object, though this could be added with sufficient interest.");
799 RCP<const Epetra_IntMultiVector> rhsImpl = rhsPtr->getEpetra_IntMultiVector ();
802 TEUCHOS_TEST_FOR_EXCEPTION(
803 rhsImpl.is_null (), std::logic_error,
"Xpetra::MultiVector::operator= "
804 "(in Xpetra::EpetraIntMultiVectorT::assign): *this (the right-hand side of "
805 "the assignment) has a null RCP<Epetra_IntMultiVector> inside. Please "
806 "report this bug to the Xpetra developers.");
807 TEUCHOS_TEST_FOR_EXCEPTION(
808 lhsImpl.is_null (), std::logic_error,
"Xpetra::MultiVector::operator= "
809 "(in Xpetra::EpetraIntMultiVectorT::assign): The left-hand side of the "
810 "assignment has a null RCP<Epetra_IntMultiVector> inside. Please report "
811 "this bug to the Xpetra developers.");
820 RCP< Epetra_IntMultiVector >
vec_;
825 #ifndef XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES
828 :
public virtual MultiVector<int,int,long long,EpetraNode>
842 vec_ = rcp(
new Epetra_IntMultiVector(toEpetra<GlobalOrdinal,Node>(map), NumVectors, zeroOut));
847 vec_ = rcp(
new Epetra_IntMultiVector(toEpetra<GlobalOrdinal,Node>(source)));
853 "Xpetra::EpetraIntMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
867 ierr = vec_->PutScalar(value);
886 Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > >
getVector(
size_t )
const {
891 Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > >
getVectorNonConst(
size_t ) {
897 Teuchos::ArrayRCP<const int>
getData(
size_t j)
const {
900 int ** arrayOfPointers;
901 vec_->ExtractView(&arrayOfPointers);
902 int * data = arrayOfPointers[j];
903 int localLength = vec_->MyLength();
905 return ArrayRCP<int>(data, 0, localLength,
false);
913 int ** arrayOfPointers;
914 vec_->ExtractView(&arrayOfPointers);
915 int * data = arrayOfPointers[j];
916 int localLength = vec_->MyLength();
918 return ArrayRCP<int>(data, 0, localLength,
false);
927 const Teuchos::ArrayView<int> &)
const {
938 "This function is not available in Epetra_IntMultiVector");
945 "This function is not implemented in Epetra_IntMultiVector");
955 void scale (Teuchos::ArrayView< const int > ) {
969 void update(
const int &,
const MultiVector<int,int,GlobalOrdinal,Node> &,
const int &,
const MultiVector<int,int,GlobalOrdinal,Node> &,
const int &) {
993 void multiply(Teuchos::ETransp , Teuchos::ETransp ,
const int &,
const MultiVector<int,int,GlobalOrdinal,Node> &,
const MultiVector<int,int,GlobalOrdinal,Node> &,
const int &) {
XPETRA_MONITOR(
"EpetraIntMultiVectorT::multiply"); TEUCHOS_TEST_FOR_EXCEPTION(1,
Xpetra::Exceptions::NotImplemented,
"Not available in Epetra"); }
998 TEUCHOS_TEST_FOR_EXCEPTION(1,
Xpetra::Exceptions::NotImplemented,
"Xpetra_EpetraIntMultiVector: elementWiseMultiply not implemented because Epetra_IntMultiVector does not support this operation");
1008 vec_->ReplaceGlobalValue(globalRow, vectorIndex, value);
1013 vec_->SumIntoGlobalValue(globalRow, vectorIndex, value);
1018 vec_->ReplaceMyValue(myRow, vectorIndex, value);
1023 vec_->SumIntoMyValue(myRow, vectorIndex, value);
1033 return vec_->NumVectors();
1047 auto vv = toEpetra<GlobalOrdinal, Node>(vec);
1048 return ( (
getLocalLength() == Teuchos::as<size_t>(vv.MyLength())) &&
1061 std::ostringstream oss;
1062 oss << Teuchos::Describable::description();
1069 void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
const {
1075 using Teuchos::VERB_DEFAULT;
1076 using Teuchos::VERB_NONE;
1077 using Teuchos::VERB_LOW;
1078 using Teuchos::VERB_MEDIUM;
1079 using Teuchos::VERB_HIGH;
1080 using Teuchos::VERB_EXTREME;
1082 if (verbLevel > Teuchos::VERB_NONE)
1095 Teuchos::RCP<const Map<int, GlobalOrdinal, Node> >
getMap ()
const {
1096 RCP<const Epetra_BlockMap> map = rcp(
new Epetra_BlockMap(vec_->Map()));
1107 const Epetra_IntMultiVector & v = *tSource.getEpetra_IntMultiVector();
1108 int err = vec_->Import(v, *tImporter.getEpetra_Import(),
toEpetra(CM));
1109 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
"Catch error code returned by Epetra.");
1119 const Epetra_IntMultiVector & v = *tDest.getEpetra_IntMultiVector();
1120 int err = vec_->Import(v, *tImporter.getEpetra_Import(),
toEpetra(CM));
1121 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
"Catch error code returned by Epetra.");
1131 const Epetra_IntMultiVector & v = *tSource.getEpetra_IntMultiVector();
1132 int err = vec_->Import(v, *tExporter.getEpetra_Export(),
toEpetra(CM));
1133 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
"Catch error code returned by Epetra.");
1143 const Epetra_IntMultiVector & v = *tDest.getEpetra_IntMultiVector();
1144 int err = vec_->Export(v, *tExporter.getEpetra_Export(),
toEpetra(CM));
1145 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
"Catch error code returned by Epetra.");
1151 if (!map.is_null()) {
1156 Epetra_SerialComm SComm;
1160 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
"Catch error code returned by Epetra.");
1166 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
1169 typename dual_view_type::t_host_um getHostLocalView ()
const {
1170 typedef Kokkos::View<
typename dual_view_type::t_host::data_type ,
1172 typename dual_view_type::t_host::device_type ,
1173 Kokkos::MemoryUnmanaged> epetra_view_type;
1178 vec_->ExtractView(&data, &myLDA);
1179 int localLength = vec_->MyLength();
1180 int numVectors = vec_->NumVectors();
1183 epetra_view_type test = epetra_view_type(data, localLength, numVectors);
1184 typename dual_view_type::t_host_um ret = subview(test, Kokkos::ALL(), Kokkos::ALL());
1189 typename dual_view_type::t_dev_um getDeviceLocalView()
const {
1190 throw std::runtime_error(
"Epetra does not support device views!");
1191 #ifndef __NVCC__ //prevent nvcc warning
1192 typename dual_view_type::t_dev_um ret;
1194 TEUCHOS_UNREACHABLE_RETURN(ret);
1207 template<
class TargetDeviceType>
1208 typename Kokkos::Impl::if_c<
1210 typename dual_view_type::t_dev_um::execution_space::memory_space,
1211 typename TargetDeviceType::memory_space>::value,
1212 typename dual_view_type::t_dev_um,
1213 typename dual_view_type::t_host_um>::type
1214 getLocalView ()
const {
1215 return this->MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::template getLocalView<TargetDeviceType>();
1228 const this_type* rhsPtr =
dynamic_cast<const this_type*
> (&rhs);
1229 TEUCHOS_TEST_FOR_EXCEPTION(
1230 rhsPtr == NULL, std::invalid_argument,
"Xpetra::MultiVector::operator=: "
1231 "The left-hand side (LHS) of the assignment has a different type than "
1232 "the right-hand side (RHS). The LHS has type Xpetra::EpetraIntMultiVectorT "
1233 "(which means it wraps an Epetra_IntMultiVector), but the RHS has some "
1234 "other type. This probably means that the RHS wraps either an "
1235 "Tpetra::MultiVector, or an Epetra_MultiVector. Xpetra::MultiVector "
1236 "does not currently implement assignment from a Tpetra object to an "
1237 "Epetra object, though this could be added with sufficient interest.");
1239 RCP<const Epetra_IntMultiVector> rhsImpl = rhsPtr->getEpetra_IntMultiVector ();
1242 TEUCHOS_TEST_FOR_EXCEPTION(
1243 rhsImpl.is_null (), std::logic_error,
"Xpetra::MultiVector::operator= "
1244 "(in Xpetra::EpetraIntMultiVectorT::assign): *this (the right-hand side of "
1245 "the assignment) has a null RCP<Epetra_IntMultiVector> inside. Please "
1246 "report this bug to the Xpetra developers.");
1247 TEUCHOS_TEST_FOR_EXCEPTION(
1248 lhsImpl.is_null (), std::logic_error,
"Xpetra::MultiVector::operator= "
1249 "(in Xpetra::EpetraIntMultiVectorT::assign): The left-hand side of the "
1250 "assignment has a null RCP<Epetra_IntMultiVector> inside. Please report "
1251 "this bug to the Xpetra developers.");
1254 *lhsImpl = *rhsImpl;
1267 #endif // XPETRA_EPETRAINTMULTIVECTOR_HPP
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. (copy).
EpetraIntMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Sets all vector entries to zero.
Teuchos::ArrayRCP< int > getDataNonConst(size_t j)
EpetraIntMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Sets all vector entries to zero.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
RCP< Epetra_IntMultiVector > getEpetra_IntMultiVector() const
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
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.
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 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.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t) const
Return a Vector which is a const view of column j.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< int, GlobalOrdinal, Node > &importer, 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.
void replaceMap(const RCP< const Map< int, GlobalOrdinal, Node > > &map)
void putScalar(const int &value)
Initialize all values in a multi-vector with specified value.
void meanValue(const Teuchos::ArrayView< int > &) const
Compute mean (average) value of each vector in multi-vector.
void scale(const int &)
Scale the current values of a multi-vector, this = alpha*this.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t)
Return a Vector which is a nonconst view of column j.
void abs(const MultiVector< int, int, GlobalOrdinal, Node > &A)
Puts element-wise absolute values of input Multi-vector in target: A = abs(this)
EpetraGlobalOrdinal GlobalOrdinal
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]) ...
Teuchos::ArrayRCP< int > getDataNonConst(size_t j)
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.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
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,j).
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
Exception throws to report errors in the internal logical of the program.
const Epetra_CrsGraph & toEpetra(const RCP< const CrsGraph< int, GlobalOrdinal, Node > > &graph)
size_t getNumVectors() const
Returns the number of vectors in the multi-vector.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< int, GlobalOrdinal, Node > &importer, CombineMode CM)
void replaceMap(const RCP< const Map< int, GlobalOrdinal, Node > > &map)
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 > > &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 scale(const int &)
Scale the current values of a multi-vector, this = alpha*this.
void meanValue(const Teuchos::ArrayView< int > &) const
Compute mean (average) value of each vector in multi-vector.
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< int, GlobalOrdinal, Node > &exporter, CombineMode CM)
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.
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute 2-norm of each vector in multi-vector.
EpetraIntMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
size_t getNumVectors() const
Returns the number of vectors in the multi-vector.
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 sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
void scale(const int &alpha)
Scale the current values of a multi-vector, this = alpha*this.
RCP< Epetra_IntMultiVector > vec_
The Epetra_IntMultiVector which this class wraps.
const RCP< const Comm< int > > getComm() const
void putScalar(const int &value)
Initialize all values in a multi-vector with specified value.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute Inf-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].
std::string description() const
Return a simple one-line description of this object.
void doImport(const DistObject< int, int, GlobalOrdinal, Node > &source, const Import< int, GlobalOrdinal, Node > &importer, CombineMode CM)
Teuchos::ArrayRCP< const int > getData(size_t j) const
const RCP< const Comm< int > > getComm() const
size_t getLocalLength() const
Returns the local vector length on the calling processor of vectors in the multi-vector.
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.
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,j).
void randomize(bool=true)
Set multi-vector values to random numbers.
void setSeed(unsigned int)
Set seed for Random function.
void scale(Teuchos::ArrayView< const int > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
~EpetraIntMultiVectorT()
Destructor.
void putScalar(const int &value)
Initialize all values in a multi-vector with specified value.
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.
size_t getNumVectors() const
Returns the number of vectors in the multi-vector.
Teuchos::RCP< const Map< int, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
RCP< Epetra_IntMultiVector > getEpetra_IntMultiVector() const
void randomize(bool bUseXpetraImplementation=true)
Set multi-vector values to random numbers.
void abs(const MultiVector< int, int, GlobalOrdinal, Node > &)
Puts element-wise absolute values of input Multi-vector in target: A = abs(this)
~EpetraIntMultiVectorT()
Destructor.
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.
~EpetraIntMultiVectorT()
Destructor.
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute Inf-norm of each vector in multi-vector.
size_t getLocalLength() const
Returns the local vector length on the calling processor of vectors in the multi-vector.
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
Teuchos::ArrayRCP< int > getDataNonConst(size_t j)
global_size_t getGlobalLength() const
Returns the global vector length of vectors in the multi-vector.
void maxValue(const Teuchos::ArrayView< int > &) const
Compute max value of each vector in multi-vector.
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.
global_size_t getGlobalLength() const
Returns the global vector length of vectors in the multi-vector.
const RCP< const Comm< int > > getComm() const
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
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
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.
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) const
Return a Vector which is a const view of column j.
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,j).
Exception throws when you call an unimplemented method of Xpetra.
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.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute 1-norm of each vector in multi-vector.
RCP< Epetra_IntMultiVector > getEpetra_IntMultiVector() const
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 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.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
void scale(Teuchos::ArrayView< const int >)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
size_t global_size_t
Global size_t object.
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &norms) const
Compute 2-norm of each vector in multi-vector.
void maxValue(const Teuchos::ArrayView< int > &) const
Compute max value of each vector in multi-vector.
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. (copy).
void randomize(bool=true)
Set multi-vector values to random numbers.
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 norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute 2-norm of each vector in multi-vector.
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 doImport(const DistObject< int, int, GlobalOrdinal, Node > &source, const Import< int, GlobalOrdinal, Node > &importer, CombineMode CM)
void doImport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< int, GlobalOrdinal, Node > &exporter, CombineMode CM)
void replaceMap(const RCP< const Map< int, GlobalOrdinal, Node > > &map)
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &norms) const
Compute 1-norm 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.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< int >::magnitudeType > &) const
Compute 1-norm of each vector in multi-vector.
RCP< Epetra_IntMultiVector > vec_
The Epetra_IntMultiVector which this class wraps.
std::string description() const
Return a simple one-line description of this object.
EpetraIntMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< int, GlobalOrdinal, Node > &exporter, CombineMode CM)
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).
void maxValue(const Teuchos::ArrayView< int > &maxs) const
Compute max value of each vector in multi-vector.
CombineMode
Xpetra::Combine Mode enumerable type.
void setSeed(unsigned int)
Set seed for Random function.
#define XPETRA_MONITOR(funcName)
void meanValue(const Teuchos::ArrayView< int > &means) const
Compute mean (average) value of each vector in multi-vector.
std::string description() const
Return a simple one-line description of this object.
void doExport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< int, GlobalOrdinal, Node > &importer, CombineMode CM)
Teuchos::ArrayRCP< const int > getData(size_t j) const
EpetraIntMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Sets all vector entries to zero.
void doImport(const DistObject< int, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< int, GlobalOrdinal, Node > &exporter, CombineMode CM)
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
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.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t)
Return a Vector which is a nonconst view of column j.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
EpetraIntMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
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).