Xpetra Version of the Day
Loading...
Searching...
No Matches
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
50
51#include "Xpetra_Vector.hpp"
54
56#include "Xpetra_Utils.hpp"
57
58#include "Tpetra_Vector.hpp"
59
60namespace Xpetra {
61
62// TODO: move that elsewhere
63template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
64RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> toTpetra(Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>&);
65
66template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
67RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> toTpetra(const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>&);
68
69template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
70RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> toXpetra(RCP<const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
71
72template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
73RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec);
74
75//
76//
77
78template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node = KokkosClassic::DefaultNode::DefaultNodeType>
80 : public virtual Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>
81 , public TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>
82{
83
84 public:
85 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dot; // overloading, not hiding
86 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::norm1; // overloading, not hiding
87 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::norm2; // overloading, not hiding
88 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::normInf; // overloading, not hiding
89 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::meanValue; // overloading, not hiding
90 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceGlobalValue; // overloading, not hiding
91 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::sumIntoGlobalValue; // overloading, not hiding
92 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceLocalValue; // overloading, not hiding
93 using TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::sumIntoLocalValue; // overloading, not hiding
94
96
97
99 TpetraVector(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>>& map, bool zeroOut = true);
100
103
105 virtual ~TpetraVector();
106
108
110
111
113 void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar& value);
114
116 void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar& value);
117
119 void replaceLocalValue(LocalOrdinal myRow, const Scalar& value);
120
122 void sumIntoLocalValue(LocalOrdinal myRow, const Scalar& value);
123
125
127
128
131
134
137
139 Scalar meanValue() const;
140
142
144
145
147 std::string description() const;
148
151
153
156
158
159
161 TpetraVector(const Teuchos::RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& vec);
162
165
167
168}; // TpetraVector class
169
170
171// TODO: move that elsewhere
172template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
175{
177 XPETRA_DYNAMIC_CAST(TpetraVectorClass, x, tX, "toTpetra");
178 return tX.getTpetra_Vector();
179}
180
181template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
182RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
184{
186 XPETRA_DYNAMIC_CAST(const TpetraVectorClass, x, tX, "toTpetra");
187 return tX.getTpetra_Vector();
188}
189
190template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
191RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
192toXpetra(RCP<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec)
193{
194 if(!vec.is_null())
196
197 return Teuchos::null;
198}
199
200template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
201RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
202toXpetra(RCP<const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> vec)
203{
204 // 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.
205 return toXpetra(Teuchos::rcp_const_cast<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(vec));
206}
207
208} // namespace Xpetra
209
210#define XPETRA_TPETRAVECTOR_SHORT
211#endif // XPETRA_TPETRAVECTOR_DECL_HPP
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
static const EVerbosityLevel verbLevel_default
void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value)
Adds specified value to existing value at the specified location.
Scalar dot(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &a) const
Computes dot product of this Vector against input Vector x.
Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const
Compute Inf-norm of this Vector.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of this Vector.
void replaceLocalValue(LocalOrdinal myRow, const Scalar &value)
Replace current value at the specified location with specified values.
void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Replace current value at the specified location with specified value.
std::string description() const
Return a simple one-line description of this object.
Scalar meanValue() const
Compute mean (average) value of this Vector.
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_Vector() const
Get the underlying Tpetra multivector.
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::ScalarTraits< Scalar >::magnitudeType norm1() const
Return 1-norm of this Vector.
virtual ~TpetraVector()
Destructor.
TpetraVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, bool zeroOut=true)
Sets all vector entries to zero.
void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value)
Adds specified value to existing value at the specified location.
Xpetra namespace
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)