Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Vector_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef TPETRA_VECTOR_DEF_HPP
11#define TPETRA_VECTOR_DEF_HPP
12
15
16#include "Tpetra_MultiVector.hpp"
18#include "KokkosCompat_View.hpp"
19#include "KokkosBlas1_nrm2w_squared.hpp"
20#include "Teuchos_CommHelpers.hpp"
21
22namespace Tpetra {
23
24template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
28
29template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31 Vector(const Teuchos::RCP<const map_type>& map,
32 const bool zeroOut)
33 : base_type(map, 1, zeroOut) {}
34
35template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40
41template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
47
48template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
50 Vector(const Teuchos::RCP<const map_type>& map,
51 const Teuchos::ArrayView<const Scalar>& values)
52 : base_type(map, values, values.size(), 1) {}
53
54template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
56 Vector(const Teuchos::RCP<const map_type>& map,
57 const dual_view_type& view)
58 : base_type(map, view) {}
59
60template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
62 Vector(const Teuchos::RCP<const map_type>& map,
63 const dual_view_type& view,
65 : base_type(map, view, origView) {}
66
67template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
69 Vector(const Teuchos::RCP<const map_type>& map,
70 const wrapped_dual_view_type& view)
71 : base_type(map, view) {}
72
73template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
78
79template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
82 this->base_type::replaceGlobalValue(globalRow, 0, value);
83}
84
85template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
88 const Scalar& value,
89 const bool atomic) {
90 this->base_type::sumIntoGlobalValue(globalRow, 0, value, atomic);
91}
92
93template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
95 replaceLocalValue(const LocalOrdinal myRow, const Scalar& value) {
96 this->base_type::replaceLocalValue(myRow, 0, value);
97}
98
99template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
102 const Scalar& value,
103 const bool atomic) {
104 this->base_type::sumIntoLocalValue(globalRow, 0, value, atomic);
105}
106
107template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
109 get1dCopy(const Teuchos::ArrayView<Scalar>& A) const {
110 const size_t lda = this->getLocalLength();
111 this->get1dCopy(A, lda);
112}
113
114template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
119 this->dot(y, Teuchos::arrayView(&result, 1));
120 return result;
121}
122
123template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
124Scalar
126 meanValue() const {
127 Scalar mean;
128 this->meanValue(Teuchos::arrayView(&mean, 1));
129 return mean;
130}
131
132template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
135 norm1() const {
137 this->norm1(Teuchos::arrayView(&norm, 1));
138 return norm;
139}
140
141template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
144 norm2() const {
146 this->norm2(Teuchos::arrayView(&norm, 1));
147 return norm;
148}
149
150template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
153 normInf() const {
155 this->normInf(Teuchos::arrayView(&norm, 1));
156 return norm;
157}
158
159template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
161 description() const {
162 return this->descriptionImpl("Tpetra::Vector");
163}
164
165template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
167 describe(Teuchos::FancyOStream& out,
168 const Teuchos::EVerbosityLevel verbLevel) const {
169 this->describeImpl(out, "Tpetra::Vector", verbLevel);
170}
171
172template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
173Teuchos::RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
175 offsetView(const Teuchos::RCP<const map_type>& subMap,
176 const size_t offset) const {
177 using Kokkos::ALL;
178 using Kokkos::subview;
179 using Teuchos::rcp;
181
182 const size_t newNumRows = subMap->getLocalNumElements();
183 const bool tooManyElts = newNumRows + offset > this->getOrigNumLocalRows();
184 if (tooManyElts) {
185 const int myRank = this->getMap()->getComm()->getRank();
187 newNumRows + offset > this->getLocalLength(), std::runtime_error,
188 "Tpetra::Vector::offsetView(NonConst): Invalid input Map. The input "
189 "Map owns "
190 << newNumRows << " entries on process " << myRank << ". "
191 "offset = "
192 << offset << ". Yet, the Vector contains only "
193 << this->getOrigNumLocalRows() << " rows on this process.");
194 }
195
196 // Need 'this->' to get view_ from parent class.
197 return rcp(new V(subMap, wrapped_dual_view_type(this->view_, Kokkos::pair<int, int>(offset, offset + newNumRows), ALL())));
198}
199
200template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
201Teuchos::RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
202Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
203 offsetViewNonConst(const Teuchos::RCP<const map_type>& subMap,
204 const size_t offset) {
205 typedef Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> V;
206 return Teuchos::rcp_const_cast<V>(this->offsetView(subMap, offset));
207}
208
209template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
210Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>
213 // The 2-argument copy constructor with second argument =
214 // Teuchos::Copy does a deep copy of its input.
215 vec_type dst(src, Teuchos::Copy);
216
217 // The Kokkos refactor version of Vector has view semantics, so
218 // returning the Vector directly, rather than through RCP, only
219 // does a shallow copy.
220 return dst;
221}
222
223} // namespace Tpetra
224
233#define TPETRA_VECTOR_INSTANT(SCALAR, LO, GO, NODE) \
234 template class Vector<SCALAR, LO, GO, NODE>; \
235 template Vector<SCALAR, LO, GO, NODE> createCopy(const Vector<SCALAR, LO, GO, NODE>& src);
236
237#endif // TPETRA_VECTOR_DEF_HPP
Declaration of a function that prints strings from each process.
Struct that holds views of the contents of a CrsMatrix.
One or more distributed dense vectors.
A distributed dense vector.
void replaceGlobalValue(const GlobalOrdinal globalRow, const Scalar &value)
Replace current value at the specified location with specified value.
void sumIntoLocalValue(const LocalOrdinal myRow, const Scalar &value, const bool atomic=base_type::useAtomicUpdatesByDefault)
Add value to existing value, using local (row) index.
mag_type normInf() const
Return the infinity-norm of this Vector.
Scalar meanValue() const
Compute mean (average) value of this Vector.
base_type::mag_type mag_type
Type of a norm result.
void replaceLocalValue(const LocalOrdinal myRow, const Scalar &value)
Replace current value at the specified location with specified values.
mag_type norm2() const
Return the two-norm of this Vector.
base_type::dot_type dot_type
Type of an inner ("dot") product result.
LocalOrdinal local_ordinal_type
This class' second template parameter; the type of local indices.
void sumIntoGlobalValue(const GlobalOrdinal globalRow, const Scalar &value, const bool atomic=base_type::useAtomicUpdatesByDefault)
Add value to existing value, using global (row) index.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Describe this object in a human-readable way to the given output stream.
Vector()
Default constructor: makes a Vector with no rows or columns.
virtual std::string description() const
Return a one-line description of this object.
void get1dCopy(const Teuchos::ArrayView< Scalar > &A) const
Return multi-vector values in user-provided two-dimensional array (using Teuchos memory management cl...
base_type::dual_view_type dual_view_type
Kokkos::DualView specialization used by this class.
dot_type dot(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &y) const
Return the dot product of this Vector and the input Vector x.
mag_type norm1() const
Return the one-norm of this Vector.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
MultiVector< ST, LO, GO, NT > createCopy(const MultiVector< ST, LO, GO, NT > &src)
Return a deep copy of the given MultiVector.