All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Xpetra_TpetraVector_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Xpetra: A linear algebra interface package
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef XPETRA_TPETRAVECTOR_DECL_HPP
47 #define XPETRA_TPETRAVECTOR_DECL_HPP
48 
49 /* this file is automatically generated - do not edit (see script/tpetra.py) */
50 
52 
53 #include "Xpetra_Vector.hpp"
54 #include "Xpetra_MultiVector.hpp"
55 #include "Xpetra_TpetraMultiVector.hpp"
56 
57 #include "Xpetra_TpetraMap.hpp" //TMP
58 #include "Xpetra_Utils.hpp"
59 #include "Xpetra_TpetraImport.hpp"
60 #include "Xpetra_TpetraExport.hpp"
61 
62 #include "Tpetra_Vector.hpp"
63 
64 namespace Xpetra {
65 
66  // TODO: move that elsewhere
67  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
68  RCP<Tpetra::Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> > toTpetra(Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
69 
70  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
71  RCP<Tpetra::Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> > toTpetra(const Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
72 
73  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
74  RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
75 
76  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
77  RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
78 
79  //
80  //
81 
82  template <class Scalar = Vector<>::scalar_type,
83  class LocalOrdinal = typename Vector<Scalar>::local_ordinal_type,
84  class GlobalOrdinal = typename Vector<Scalar, LocalOrdinal>::global_ordinal_type,
85  class Node = typename Vector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
86  class TpetraVector
87  : public virtual Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>,
88  public TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
89  {
90 #undef XPETRA_TPETRAMULTIVECTOR_SHORT
91 #include "Xpetra_UseShortNames.hpp"
92 
93  public:
94 
95  using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dot; // overloading, not hiding
96  using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::norm1; // overloading, not hiding
97  using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::norm2; // overloading, not hiding
98  using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::normInf; // overloading, not hiding
99  using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::meanValue; // overloading, not hiding
100  using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceGlobalValue; // overloading, not hiding
101  using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::sumIntoGlobalValue; // overloading, not hiding
102  using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceLocalValue; // overloading, not hiding
103  using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::sumIntoLocalValue; // overloading, not hiding
104 
106 
107 
109  TpetraVector(const Teuchos::RCP<const Map> &map, bool zeroOut=true);
110 
113 
115  virtual ~TpetraVector();
116 
118 
120 
121 
123  void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value);
124 
126  void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value);
127 
129  void replaceLocalValue(LocalOrdinal myRow, const Scalar &value);
130 
132  void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value);
133 
135 
137 
138 
141 
144 
147 
149  Scalar meanValue() const;
150 
152 
154 
155 
157  std::string description() const;
158 
161 
163 
165  Scalar dot(const Vector &a) const;
166 
168  typename Teuchos::ScalarTraits< Scalar >::magnitudeType normWeighted(const Vector &weights) const;
169 
170 
172 
173 
175  TpetraVector(const Teuchos::RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &vec);
176 
178  RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > getTpetra_Vector() const;
179 
180 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
181 
183 
184  typename dual_view_type::t_host_um getHostLocalView () const;
185 
186  typename dual_view_type::t_dev_um getDeviceLocalView() const;
187 
198  template<class TargetDeviceType>
199  typename Kokkos::Impl::if_c<
200  Kokkos::Impl::is_same<
201  typename dual_view_type::t_dev_um::execution_space::memory_space,
202  typename TargetDeviceType::memory_space>::value,
203  typename dual_view_type::t_dev_um,
204  typename dual_view_type::t_host_um>::type
205  getLocalView () const;
206 #endif
207 
209 
210  }; // TpetraVector class
211 
212 
213  // TODO: move that elsewhere
214  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
217  XPETRA_DYNAMIC_CAST( TpetraVectorClass, x, tX, "toTpetra");
218  return tX.getTpetra_Vector();
219  }
220 
221  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
224  XPETRA_DYNAMIC_CAST(const TpetraVectorClass, x, tX, "toTpetra");
225  return tX.getTpetra_Vector();
226  }
227 
228  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
229  RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec) {
230  if (!vec.is_null())
232 
233  return Teuchos::null;
234  }
235 
236  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
237  RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec) {
238  // 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.
239  return toXpetra(Teuchos::rcp_const_cast<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > (vec));
240  }
241 
242 } // Xpetra namespace
243 
244 #define XPETRA_TPETRAVECTOR_SHORT
245 #endif // XPETRA_TPETRAVECTOR_DECL_HPP
246 
virtual ~TpetraVector()
Destructor.
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 normWeighted(const Vector &weights) const
Compute Weighted 2-norm (RMS Norm) of this Vector.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Scalar dot(const Vector &a) const
Computes dot product of this Vector against input Vector x.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const
Return 1-norm of this Vector.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
static const EVerbosityLevel verbLevel_default
void replaceLocalValue(LocalOrdinal myRow, const Scalar &value)
Replace current value at the specified location with specified values.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_Vector() const
Get the underlying Tpetra multivector.
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.
std::string description() const
Return a simple one-line description of this object.
TpetraVector(const Teuchos::RCP< const Map > &map, bool zeroOut=true)
Sets all vector entries to zero.
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.