Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_TpetraMultiVector_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Xpetra: A linear algebra interface package
4 //
5 // Copyright 2012 NTESS and the Xpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef XPETRA_TPETRAMULTIVECTOR_DECL_HPP
11 #define XPETRA_TPETRAMULTIVECTOR_DECL_HPP
12 
15 
19 #include "Tpetra_MultiVector.hpp"
20 #include "Tpetra_Vector.hpp"
21 #include "Xpetra_Utils.hpp"
22 
23 namespace Xpetra {
24 
25 // TODO: move that elsewhere
26 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
27 const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &toTpetra(const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &);
28 
29 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
30 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &toTpetra(MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &);
31 
32 #ifndef DOXYGEN_SHOULD_SKIP_THIS
33 // forward declaration of TpetraVector, needed to prevent circular inclusions
34 template <class S, class LO, class GO, class N>
35 class TpetraVector;
36 #endif
37 
38 // Because we aren't including the header...
39 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > vec);
41 
42 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
43 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > toXpetra(RCP<const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > vec);
44 
45 template <class Scalar,
46  class LocalOrdinal,
47  class GlobalOrdinal,
48  class Node>
50  : public virtual MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
51  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
53 
54  public:
56 
57 
59  TpetraMultiVector(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &map, size_t NumVectors, bool zeroOut = true);
60 
62  TpetraMultiVector(const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &source, const Teuchos::DataAccess copyOrView = Teuchos::Copy);
63 
65  TpetraMultiVector(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &map, const Teuchos::ArrayView<const Scalar> &A, size_t LDA, size_t NumVectors);
66 
68  TpetraMultiVector(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &map, const Teuchos::ArrayView<const Teuchos::ArrayView<const Scalar> > &ArrayOfPtrs, size_t NumVectors);
69 
70  virtual ~TpetraMultiVector();
71 
73 
74 
76  void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value);
77 
79  void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value);
80 
82  void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value);
83 
85  void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value);
86 
88  void putScalar(const Scalar &value);
89 
91  void reduce();
92 
94 
96 
97 
99  Teuchos::RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > getVector(size_t j) const;
100 
102  Teuchos::RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > getVectorNonConst(size_t j);
103 
105  Teuchos::ArrayRCP<const Scalar> getData(size_t j) const;
106 
108  Teuchos::ArrayRCP<Scalar> getDataNonConst(size_t j);
109 
111  void get1dCopy(Teuchos::ArrayView<Scalar> A, size_t LDA) const;
112 
114  void get2dCopy(Teuchos::ArrayView<const Teuchos::ArrayView<Scalar> > ArrayOfPtrs) const;
115 
117  Teuchos::ArrayRCP<const Scalar> get1dView() const;
118 
120  Teuchos::ArrayRCP<Teuchos::ArrayRCP<const Scalar> > get2dView() const;
121 
123  Teuchos::ArrayRCP<Scalar> get1dViewNonConst();
124 
126  Teuchos::ArrayRCP<Teuchos::ArrayRCP<Scalar> > get2dViewNonConst();
127 
129 
131 
132 
134  void dot(const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &A, const Teuchos::ArrayView<Scalar> &dots) const;
135 
138 
141 
143  void scale(const Scalar &alpha);
144 
146  void scale(Teuchos::ArrayView<const Scalar> alpha);
147 
149  void scale(const Scalar &alpha, const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &A);
150 
152  void update(const Scalar &alpha, const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &A, const Scalar &beta);
153 
155  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);
156 
158  void norm1(const Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType> &norms) const;
159 
161  void norm2(const Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType> &norms) const;
162 
164  void normInf(const Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType> &norms) const;
165 
167  void meanValue(const Teuchos::ArrayView<Scalar> &means) const;
168 
170  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);
171 
173 
175 
177  size_t getNumVectors() const;
178 
180  size_t getLocalLength() const;
181 
184 
185  // \! Checks to see if the local length, number of vectors and size of Scalar type match
187 
189 
191 
193  std::string description() const;
194 
196  void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default) const;
197 
199 
201  void elementWiseMultiply(Scalar scalarAB, const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &A, const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &B, Scalar scalarThis); // definition at the end of this file
202  // TODO: void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis){ vec_->elementWiseMultiply(scalarAB, toTpetra(A), toTpetra(B), scalarThis); }
203 
205  void randomize(bool bUseXpetraImplementation = false);
206 
207  void randomize(const Scalar &minVal, const Scalar &maxVal, bool bUseXpetraImplementation = false);
208 
209  //{@
210  // Implements DistObject interface
211 
212  Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > getMap() const;
213 
215 
217 
219 
221 
223 
225 
227 
229 
231 
233 
235 
237 
238  void replaceMap(const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &map);
239 
241 
243 
244 
246  TpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &vec);
247 
249  RCP<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > getTpetra_MultiVector() const;
250 
252  void setSeed(unsigned int seed);
253 
255 
256  virtual typename dual_view_type::t_host_const_um getHostLocalView(Access::ReadOnlyStruct) const;
257 
258  virtual typename dual_view_type::t_dev_const_um getDeviceLocalView(Access::ReadOnlyStruct) const;
259 
260  virtual typename dual_view_type::t_host_um getHostLocalView(Access::OverwriteAllStruct) const;
261 
262  virtual typename dual_view_type::t_dev_um getDeviceLocalView(Access::OverwriteAllStruct) const;
263 
264  virtual typename dual_view_type::t_host_um getHostLocalView(Access::ReadWriteStruct) const;
265 
266  virtual typename dual_view_type::t_dev_um getDeviceLocalView(Access::ReadWriteStruct) const;
267 
268  protected:
271  virtual void
273 
274  private:
276  RCP<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > vec_;
277 
278 }; // TpetraMultiVector class
279 
280 // Things we actually need
281 
282 // Things we actually need
283 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
284 RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > toXpetra(RCP<Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > vec) {
285  if (!vec.is_null())
287 
288  return Teuchos::null;
289 }
290 
291 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
292 RCP<const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > toXpetra(RCP<const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > vec) {
293  if (!vec.is_null())
295 
296  return Teuchos::null;
297 }
298 
299 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
300 const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &toTpetra(const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &x) {
302  XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, x, tX, "toTpetra");
303  return *tX.getTpetra_MultiVector();
304 }
305 
306 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
307 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &toTpetra(MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &x) {
309  XPETRA_DYNAMIC_CAST(TpetraMultiVectorClass, x, tX, "toTpetra");
310  return *tX.getTpetra_MultiVector();
311 }
312 
313 } // namespace Xpetra
314 
315 #define XPETRA_TPETRAMULTIVECTOR_SHORT
316 #endif // XPETRA_TPETRAMULTIVECTOR_HPP
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
Teuchos::ArrayRCP< const Scalar > get1dView() const
Const persisting (1-D) view of this multivector&#39;s local values.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const
Return const persisting pointers to values.
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst()
Return non-const persisting pointers to values.
std::string description() const
A simple one-line description of this object.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Element-wise multiply of a Vector A with a TpetraMultiVector B.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
size_t getNumVectors() const
Number of columns in the multivector.
void setSeed(unsigned int seed)
Set seed for Random function.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void reduce()
Sum values of a locally replicated multivector across all processes.
void beginImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import data into this object using an Import object (&quot;forward mode&quot;).
void beginExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export data into this object using an Import object (&quot;reverse mode&quot;).
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dual_view_type dual_view_type
TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import data into this object using an Import object (&quot;forward mode&quot;).
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 endImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import data into this object using an Import object (&quot;forward mode&quot;).
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 abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
size_t getLocalLength() const
Local number of rows on the calling process.
void endExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export data into this object using an Import object (&quot;reverse mode&quot;).
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const
Fill the given array with a copy of this multivector&#39;s local values.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst()
Nonconst persisting (1-D) view of this multivector&#39;s local values.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec_
The Tpetra::MultiVector which this class wraps.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
void scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t global_size_t
Global size_t object.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic constuctor.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
Kokkos::DualView< impl_scalar_type **, Kokkos::LayoutStride, typename node_type::device_type, Kokkos::MemoryUnmanaged > dual_view_type
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
virtual dual_view_type::t_dev_const_um getDeviceLocalView(Access::ReadOnlyStruct) const
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const
Fill the given array with a copy of this multivector&#39;s local values.
virtual ~TpetraMultiVector()
Destructor (virtual for memory safety of derived classes).
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export data into this object using an Import object (&quot;reverse mode&quot;).
CombineMode
Xpetra::Combine Mode enumerable type.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra 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 dual_view_type::t_host_const_um getHostLocalView(Access::ReadOnlyStruct) const
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution 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.
global_size_t getGlobalLength() const
Global number of rows in the multivector.
void replaceMap(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)