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>
81 template<
class Scalar,
class LO,
class GO,
class Node>
85 vPtr->setCopyOrView (Teuchos::View);
89 template<
class Scalar,
class LO,
class GO,
class Node>
93 return ((
const base_type*)
this)->replaceLocalValues (localRowIndex, 0, vals);
96 template<
class Scalar,
class LO,
class GO,
class Node>
100 return ((
const base_type*)
this)->replaceGlobalValues (globalRowIndex, 0, vals);
104 template<
class Scalar,
class LO,
class GO,
class Node>
108 return ((
const base_type*)
this)->sumIntoLocalValues (localRowIndex, 0, vals);
111 template<
class Scalar,
class LO,
class GO,
class Node>
115 return ((
const base_type*)
this)->sumIntoLocalValues (globalRowIndex, 0, vals);
118 template<
class Scalar,
class LO,
class GO,
class Node>
122 return ((
const base_type*)
this)->getLocalRowView (localRowIndex, 0, vals);
125 template<
class Scalar,
class LO,
class GO,
class Node>
129 return ((
const base_type*)
this)->getGlobalRowView (globalRowIndex, 0, vals);
143 template<
class Scalar,
class LO,
class GO,
class Node>
153 const size_t offset = localRowIndex * blockSize * strideX;
166 #define TPETRA_EXPERIMENTAL_BLOCKVECTOR_INSTANT(S,LO,GO,NODE) \ 167 template class Experimental::BlockVector< S, LO, GO, NODE >; 169 #endif // TPETRA_EXPERIMENTAL_BLOCKVECTOR_DEF_HPP bool getLocalRowView(const LO localRowIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a local index. ...
Namespace Tpetra contains the class and methods constituting the Tpetra library.
mv_type mv_
The Tpetra::MultiVector used to represent the data.
MultiVector for multiple degrees of freedom per mesh point.
LO getBlockSize() const
Get the number of degrees of freedom per mesh point.
bool sumIntoGlobalValues(const GO globalRowIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a global index.
impl_scalar_type * getRawPtr() const
Raw pointer to the MultiVector's data.
LittleVector< impl_scalar_type, LO > little_vec_type
"Block view" of all degrees of freedom at a mesh point.
Nonowning view of a set of degrees of freedom corresponding to a mesh point in a block vector or mult...
vec_type getVectorView()
Get a Tpetra::Vector that views this BlockVector's data.
bool getGlobalRowView(const GO globalRowIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a global index.
bool isValidLocalMeshIndex(const LO meshLocalIndex) const
True if and only if meshLocalIndex is a valid local index in the mesh Map.
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.
bool sumIntoLocalValues(const LO localRowIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a local index.
A distributed dense vector.
bool replaceGlobalValues(const GO globalRowIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a global index.
size_t getNumVectors() const
Number of columns in the multivector.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node, classic > > getVectorNonConst(const size_t j)
Return a Vector which is a nonconst view of column j.
size_t getStrideX() const
Stride between consecutive local entries in the same column.
BlockVector()
Default constructor.
bool replaceLocalValues(const LO localRowIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a local index.