42 #ifndef TPETRA_EXPERIMENTAL_BLOCKVECTOR_DEF_HPP 43 #define TPETRA_EXPERIMENTAL_BLOCKVECTOR_DEF_HPP 45 #include "Tpetra_Experimental_BlockVector_decl.hpp" 48 namespace Experimental {
50 template<
class Scalar,
class LO,
class GO,
class Node>
56 template<
class Scalar,
class LO,
class GO,
class Node>
61 base_type (meshMap, pointMap, blockSize, 1)
64 template<
class Scalar,
class LO,
class GO,
class Node>
71 TEUCHOS_TEST_FOR_EXCEPTION(
73 "Tpetra::Experimental::BlockVector: Input MultiVector has " 77 template<
class Scalar,
class LO,
class GO,
class Node>
82 const size_t offset) :
83 base_type (X, newMeshMap, newPointMap, offset)
86 template<
class Scalar,
class LO,
class GO,
class Node>
90 const size_t offset) :
94 template<
class Scalar,
class LO,
class GO,
class Node>
98 template<
class Scalar,
class LO,
class GO,
class Node>
102 vPtr->setCopyOrView (Teuchos::View);
106 template<
class Scalar,
class LO,
class GO,
class Node>
110 return ((
const base_type*)
this)->replaceLocalValues (localRowIndex, 0, vals);
113 template<
class Scalar,
class LO,
class GO,
class Node>
117 return ((
const base_type*)
this)->replaceGlobalValues (globalRowIndex, 0, vals);
121 template<
class Scalar,
class LO,
class GO,
class Node>
125 return ((
const base_type*)
this)->sumIntoLocalValues (localRowIndex, 0, vals);
128 template<
class Scalar,
class LO,
class GO,
class Node>
132 return ((
const base_type*)
this)->sumIntoLocalValues (globalRowIndex, 0, vals);
135 template<
class Scalar,
class LO,
class GO,
class Node>
139 return ((
const base_type*)
this)->getLocalRowView (localRowIndex, 0, vals);
142 template<
class Scalar,
class LO,
class GO,
class Node>
146 return ((
const base_type*)
this)->getGlobalRowView (globalRowIndex, 0, vals);
160 template<
class Scalar,
class LO,
class GO,
class Node>
165 if (! this->isValidLocalMeshIndex (localRowIndex)) {
168 const size_t blockSize = this->getBlockSize ();
169 const size_t offset = localRowIndex * blockSize;
182 #define TPETRA_EXPERIMENTAL_BLOCKVECTOR_INSTANT(S,LO,GO,NODE) \ 183 template class Experimental::BlockVector< S, LO, GO, NODE >; 185 #endif // TPETRA_EXPERIMENTAL_BLOCKVECTOR_DEF_HPP Namespace Tpetra contains the class and methods constituting the Tpetra library.
bool sumIntoGlobalValues(const GO globalRowIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a global index.
bool getGlobalRowView(const GO globalRowIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a global index.
MultiVector for multiple degrees of freedom per mesh point.
Vector for multiple degrees of freedom per mesh point.
Kokkos::View< impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_vec_type
"Block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector...
Kokkos::View< impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_vec_type
"Block view" of all degrees of freedom at a mesh point.
bool sumIntoLocalValues(const LO localRowIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a local index.
bool replaceLocalValues(const LO localRowIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a local index.
vec_type getVectorView()
Get a Tpetra::Vector that views this BlockVector's data.
bool getLocalRowView(const LO localRowIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a local index. ...
little_vec_type getLocalBlock(const LO localRowIndex) const
Get a view of the degrees of freedom at the given mesh point, using a local index.
size_t getNumVectors() const
Number of columns in the multivector.
bool replaceGlobalValues(const GO globalRowIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a global index.
A distributed dense vector.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node, classic > > getVectorNonConst(const size_t j)
Return a Vector which is a nonconst view of column j.
BlockVector()
Default constructor.