10 #ifndef TPETRA_VECTOR_DEF_HPP
11 #define TPETRA_VECTOR_DEF_HPP
16 #include "Tpetra_MultiVector.hpp"
18 #include "KokkosCompat_View.hpp"
19 #include "KokkosBlas1_nrm2w_squared.hpp"
20 #include "Teuchos_CommHelpers.hpp"
24 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
30 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
32 Vector (
const Teuchos::RCP<const map_type>& map,
37 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
40 const Teuchos::DataAccess copyOrView)
44 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
47 const Teuchos::RCP<const map_type>& map,
52 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
54 Vector (
const Teuchos::RCP<const map_type>& map,
55 const Teuchos::ArrayView<const Scalar>& values)
56 :
base_type (map, values, values.size (), 1)
59 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
61 Vector (
const Teuchos::RCP<const map_type>& map,
66 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
68 Vector (
const Teuchos::RCP<const map_type>& map,
75 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
77 Vector (
const Teuchos::RCP<const map_type>& map,
82 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
89 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
93 this->base_type::replaceGlobalValue (globalRow, 0, value);
96 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
103 this->base_type::sumIntoGlobalValue (globalRow, 0, value, atomic);
106 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
110 this->base_type::replaceLocalValue (myRow, 0, value);
113 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
120 this->base_type::sumIntoLocalValue (globalRow, 0, value, atomic);
123 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
127 const size_t lda = this->getLocalLength ();
128 this->get1dCopy (A, lda);
131 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
137 this->dot (y, Teuchos::arrayView (&result, 1));
141 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
147 this->meanValue (Teuchos::arrayView (&mean, 1));
151 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
157 this->norm1 (Teuchos::arrayView (&norm, 1));
161 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
167 this->norm2 (Teuchos::arrayView (&norm, 1));
171 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
177 this->normInf (Teuchos::arrayView (&norm, 1));
182 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
186 return this->descriptionImpl (
"Tpetra::Vector");
189 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
192 const Teuchos::EVerbosityLevel verbLevel)
const
194 this->describeImpl (out,
"Tpetra::Vector", verbLevel);
197 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
198 Teuchos::RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
200 offsetView (
const Teuchos::RCP<const map_type>& subMap,
201 const size_t offset)
const
204 using Kokkos::subview;
208 const size_t newNumRows = subMap->getLocalNumElements ();
209 const bool tooManyElts = newNumRows + offset > this->getOrigNumLocalRows ();
211 const int myRank = this->getMap ()->getComm ()->getRank ();
212 TEUCHOS_TEST_FOR_EXCEPTION(
213 newNumRows + offset > this->getLocalLength (), std::runtime_error,
214 "Tpetra::Vector::offsetView(NonConst): Invalid input Map. The input "
215 "Map owns " << newNumRows <<
" entries on process " << myRank <<
". "
216 "offset = " << offset <<
". Yet, the Vector contains only "
217 << this->getOrigNumLocalRows () <<
" rows on this process.");
221 return rcp (
new V (subMap, wrapped_dual_view_type(this->view_,Kokkos::pair<int,int>(offset,offset+newNumRows),ALL())));
224 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
225 Teuchos::RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
226 Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
227 offsetViewNonConst (
const Teuchos::RCP<const map_type>& subMap,
230 typedef Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> V;
231 return Teuchos::rcp_const_cast<V> (this->offsetView (subMap, offset));
234 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
235 Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>
241 vec_type dst (src, Teuchos::Copy);
259 #define TPETRA_VECTOR_INSTANT(SCALAR,LO,GO,NODE) \
260 template class Vector< SCALAR , LO , GO , NODE >; \
261 template Vector< SCALAR , LO , GO , NODE > createCopy (const Vector< SCALAR , LO , GO , NODE >& src);
263 #endif // TPETRA_VECTOR_DEF_HPP
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.
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...
void sumIntoLocalValue(const LocalOrdinal myRow, const Scalar &value, const bool atomic=base_type::useAtomicUpdatesByDefault)
Add value to existing value, using local (row) index.
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 sumIntoGlobalValue(const GlobalOrdinal globalRow, const Scalar &value, const bool atomic=base_type::useAtomicUpdatesByDefault)
Add value to existing value, using global (row) index.
void replaceLocalValue(const LocalOrdinal myRow, const Scalar &value)
Replace current value at the specified location with specified values.
void replaceGlobalValue(const GlobalOrdinal globalRow, const Scalar &value)
Replace current value at the specified location with specified value.
LocalOrdinal local_ordinal_type
This class' second template parameter; the type of local indices.
A distributed dense vector.
base_type::dual_view_type dual_view_type
Kokkos::DualView specialization used by this class.
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.
mag_type norm1() const
Return the one-norm of this Vector.