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" 
   24template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   29template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   31    Vector(
const Teuchos::RCP<const map_type>& 
map,
 
 
   35template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   41template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   44           const Teuchos::RCP<const map_type>& 
map,
 
 
   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)
 
 
   54template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   56    Vector(
const Teuchos::RCP<const map_type>& 
map,
 
 
   60template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   62    Vector(
const Teuchos::RCP<const map_type>& 
map,
 
 
   67template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   69    Vector(
const Teuchos::RCP<const map_type>& 
map,
 
 
   73template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   79template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   82  this->base_type::replaceGlobalValue(
globalRow, 0, value);
 
 
   85template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   90  this->base_type::sumIntoGlobalValue(
globalRow, 0, value, atomic);
 
 
   93template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   96  this->base_type::replaceLocalValue(
myRow, 0, value);
 
 
   99template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  104  this->base_type::sumIntoLocalValue(
globalRow, 0, value, atomic);
 
 
  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);
 
 
  114template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  119  this->dot(
y, Teuchos::arrayView(&
result, 1));
 
 
  123template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  128  this->meanValue(Teuchos::arrayView(&
mean, 1));
 
 
  132template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  137  this->norm1(Teuchos::arrayView(&
norm, 1));
 
 
  141template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  146  this->norm2(Teuchos::arrayView(&
norm, 1));
 
 
  150template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  155  this->normInf(Teuchos::arrayView(&
norm, 1));
 
 
  159template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  162  return this->descriptionImpl(
"Tpetra::Vector");
 
 
  165template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  168             const Teuchos::EVerbosityLevel 
verbLevel)
 const {
 
 
  172template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  173Teuchos::RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
 
  176               const size_t offset)
 const {
 
  178  using Kokkos::subview;
 
  185    const int myRank = this->getMap()->getComm()->getRank();
 
  188        "Tpetra::Vector::offsetView(NonConst): Invalid input Map.  The input " 
  192            << 
offset << 
".  Yet, the Vector contains only " 
  193            << this->getOrigNumLocalRows() << 
" rows on this process.");
 
  197  return rcp(
new V(subMap, wrapped_dual_view_type(this->view_, Kokkos::pair<int, int>(offset, offset + newNumRows), ALL())));
 
  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));
 
  209template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  210Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>
 
  215  vec_type dst(src, Teuchos::Copy);
 
 
  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); 
 
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.