42 #ifndef TPETRA_VECTOR_DEF_HPP
43 #define TPETRA_VECTOR_DEF_HPP
53 #include "Tpetra_MultiVector.hpp"
55 #include "KokkosCompat_View.hpp"
56 #include "KokkosBlas1_nrm2w_squared.hpp"
57 #include "Teuchos_CommHelpers.hpp"
61 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
67 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
69 Vector (
const Teuchos::RCP<const map_type>& map,
74 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
77 const Teuchos::DataAccess copyOrView)
81 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
84 const Teuchos::RCP<const map_type>& map,
89 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
91 Vector (
const Teuchos::RCP<const map_type>& map,
92 const Teuchos::ArrayView<const Scalar>& values)
93 :
base_type (map, values, values.size (), 1)
96 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
98 Vector (
const Teuchos::RCP<const map_type>& map,
103 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
105 Vector (
const Teuchos::RCP<const map_type>& map,
111 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
118 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
122 this->base_type::replaceGlobalValue (globalRow, 0, value);
125 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
130 const bool atomic)
const
132 this->base_type::sumIntoGlobalValue (globalRow, 0, value, atomic);
135 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
139 this->base_type::replaceLocalValue (myRow, 0, value);
142 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
147 const bool atomic)
const
149 this->base_type::sumIntoLocalValue (globalRow, 0, value, atomic);
152 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
156 const size_t lda = this->getLocalLength ();
157 this->get1dCopy (A, lda);
160 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
166 this->dot (y, Teuchos::arrayView (&result, 1));
170 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
176 this->meanValue (Teuchos::arrayView (&mean, 1));
180 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
186 this->norm1 (Teuchos::arrayView (&norm, 1));
190 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
196 this->norm2 (Teuchos::arrayView (&norm, 1));
200 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
206 this->normInf (Teuchos::arrayView (&norm, 1));
211 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
215 return this->descriptionImpl (
"Tpetra::Vector");
218 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
221 const Teuchos::EVerbosityLevel verbLevel)
const
223 this->describeImpl (out,
"Tpetra::Vector", verbLevel);
226 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
227 Teuchos::RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
229 offsetView (
const Teuchos::RCP<const map_type>& subMap,
230 const size_t offset)
const
233 using Kokkos::subview;
237 const size_t newNumRows = subMap->getNodeNumElements ();
238 const bool tooManyElts = newNumRows + offset > this->getOrigNumLocalRows ();
240 const int myRank = this->getMap ()->getComm ()->getRank ();
241 TEUCHOS_TEST_FOR_EXCEPTION(
242 newNumRows + offset > this->getLocalLength (), std::runtime_error,
243 "Tpetra::Vector::offsetView(NonConst): Invalid input Map. The input "
244 "Map owns " << newNumRows <<
" entries on process " << myRank <<
". "
245 "offset = " << offset <<
". Yet, the Vector contains only "
246 << this->getOrigNumLocalRows () <<
" rows on this process.");
249 const std::pair<size_t, size_t> offsetPair (offset, offset + newNumRows);
251 return rcp (
new V (subMap,
252 subview (this->view_, offsetPair, ALL ()),
256 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
257 Teuchos::RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
258 Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
259 offsetViewNonConst (
const Teuchos::RCP<const map_type>& subMap,
262 typedef Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> V;
263 return Teuchos::rcp_const_cast<V> (this->offsetView (subMap, offset));
266 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
267 Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>
273 vec_type dst (src, Teuchos::Copy);
291 #define TPETRA_VECTOR_INSTANT(SCALAR,LO,GO,NODE) \
292 template class Vector< SCALAR , LO , GO , NODE >; \
293 template Vector< SCALAR , LO , GO , NODE > createCopy (const Vector< SCALAR , LO , GO , NODE >& src);
295 #endif // TPETRA_VECTOR_DEF_HPP
void sumIntoGlobalValue(const GlobalOrdinal globalRow, const Scalar &value, const bool atomic=base_type::useAtomicUpdatesByDefault) const
Add value to existing value, using global (row) index.
Vector()
Default constructor: makes a Vector with no rows or columns.
mag_type norm2() const
Return the two-norm of this Vector.
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.
base_type::mag_type mag_type
Type of a norm result.
void replaceGlobalValue(const GlobalOrdinal globalRow, const Scalar &value) const
Replace current value at the specified location with specified value.
Kokkos::DualView< impl_scalar_type **, Kokkos::LayoutLeft, execution_space > dual_view_type
Kokkos::DualView specialization used by this class.
Declaration of a function that prints strings from each process.
void get1dCopy(const Teuchos::ArrayView< Scalar > &A) const
Return multi-vector values in user-provided two-dimensional array (using Teuchos memory management cl...
One or more distributed dense vectors.
virtual std::string description() const
Return a one-line description of this object.
MultiVector< ST, LO, GO, NT > createCopy(const MultiVector< ST, LO, GO, NT > &src)
Return a deep copy of the given MultiVector.
dot_type dot(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &y) const
Return the dot product of this Vector and the input Vector x.
void sumIntoLocalValue(const LocalOrdinal myRow, const Scalar &value, const bool atomic=base_type::useAtomicUpdatesByDefault) const
Add value to existing value, using local (row) index.
typename Kokkos::Details::InnerProductSpaceTraits< impl_scalar_type >::dot_type dot_type
Type of an inner ("dot") product result.
typename Kokkos::ArithTraits< impl_scalar_type >::mag_type mag_type
Type of a norm result.
A distributed dense vector.
typename map_type::local_ordinal_type local_ordinal_type
The type of local indices that this class uses.
Scalar meanValue() const
Compute mean (average) value of this Vector.
mag_type normInf() const
Return the infinity-norm of this Vector.
base_type::dot_type dot_type
Type of an inner ("dot") product result.
Base class for distributed Tpetra objects that support data redistribution.
mag_type norm1() const
Return the one-norm of this Vector.
void replaceLocalValue(const LocalOrdinal myRow, const Scalar &value) const
Replace current value at the specified location with specified values.