Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_TpetraVector_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_TPETRAVECTOR_DECL_HPP
11 #define XPETRA_TPETRAVECTOR_DECL_HPP
12 
14 
15 #include "Xpetra_Vector.hpp"
18 
20 #include "Xpetra_Utils.hpp"
21 
22 #include "Tpetra_Vector.hpp"
23 
24 namespace Xpetra {
25 
26 // TODO: move that elsewhere
27 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
28 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> toTpetra(Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>&);
29 
30 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> toTpetra(const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>&);
32 
33 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
34 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> toXpetra(RCP<const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
35 
36 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
37 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
38 
39 //
40 //
41 
42 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
43 class TpetraVector
44  : public virtual Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
45  public TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
46  public:
47  using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dot; // overloading, not hiding
48  using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::norm1; // overloading, not hiding
49  using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::norm2; // overloading, not hiding
50  using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::normInf; // overloading, not hiding
51  using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::meanValue; // overloading, not hiding
52  using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceGlobalValue; // overloading, not hiding
53  using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::sumIntoGlobalValue; // overloading, not hiding
54  using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceLocalValue; // overloading, not hiding
55  using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::sumIntoLocalValue; // overloading, not hiding
56 
58 
59 
61  TpetraVector(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>>& map, bool zeroOut = true);
62 
64  TpetraVector(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>>& map, const Teuchos::ArrayView<const Scalar>& A);
65 
67  virtual ~TpetraVector();
68 
70 
72 
73 
75  void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar& value);
76 
78  void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar& value);
79 
81  void replaceLocalValue(LocalOrdinal myRow, const Scalar& value);
82 
84  void sumIntoLocalValue(LocalOrdinal myRow, const Scalar& value);
85 
87 
89 
90 
92  typename Teuchos::ScalarTraits<Scalar>::magnitudeType norm1() const;
93 
95  typename Teuchos::ScalarTraits<Scalar>::magnitudeType norm2() const;
96 
98  typename Teuchos::ScalarTraits<Scalar>::magnitudeType normInf() const;
99 
101  Scalar meanValue() const;
102 
104 
106 
107 
109  std::string description() const;
110 
112  void describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default) const;
113 
115 
117  Scalar dot(const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& a) const;
118 
120 
121 
123  TpetraVector(const Teuchos::RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& vec);
124 
126  RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> getTpetra_Vector() const;
127 
129 
130 }; // TpetraVector class
131 
132 // TODO: move that elsewhere
133 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
134 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
137  XPETRA_DYNAMIC_CAST(TpetraVectorClass, x, tX, "toTpetra");
138  return tX.getTpetra_Vector();
139 }
140 
141 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
142 RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
145  XPETRA_DYNAMIC_CAST(const TpetraVectorClass, x, tX, "toTpetra");
146  return tX.getTpetra_Vector();
147 }
148 
149 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
150 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
151 toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec) {
152  if (!vec.is_null())
154 
155  return Teuchos::null;
156 }
157 
158 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
159 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
160 toXpetra(RCP<const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec) {
161  // We cast away the const to wrap the Tpetra vector into an Xpetra object. But it's OK because the Xpetra vector is returned as const.
162  return toXpetra(Teuchos::rcp_const_cast<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(vec));
163 }
164 
165 } // namespace Xpetra
166 
167 #define XPETRA_TPETRAVECTOR_SHORT
168 #endif // XPETRA_TPETRAVECTOR_DECL_HPP
virtual ~TpetraVector()
Destructor.
Scalar dot(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &a) const
Computes dot product of this Vector against input Vector x.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of this Vector.
void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Adds specified value to existing value at the specified location.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const
Return 1-norm of this Vector.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
TpetraVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &map, bool zeroOut=true)
Sets all vector entries to zero.
void replaceLocalValue(LocalOrdinal myRow, const Scalar &value)
Replace current value at the specified location with specified values.
Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const
Compute Inf-norm of this 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, const Scalar &value)
Adds specified value to existing value at the specified location.
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_Vector() const
Get the underlying Tpetra multivector.
std::string description() const
Return a simple one-line description of this object.
void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Replace current value at the specified location with specified value.
Scalar meanValue() const
Compute mean (average) value of this Vector.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)