42#ifndef TPETRA_BLOCKVECTOR_DEF_HPP
43#define TPETRA_BLOCKVECTOR_DEF_HPP
47 template<
class Scalar,
class LO,
class GO,
class Node>
53 template<
class Scalar,
class LO,
class GO,
class Node>
60 template<
class Scalar,
class LO,
class GO,
class Node>
66 template<
class Scalar,
class LO,
class GO,
class Node>
74 template<
class Scalar,
class LO,
class GO,
class Node>
82 X_mv.getNumVectors () != 1, std::invalid_argument,
83 "Tpetra::BlockVector: Input MultiVector has "
84 <<
X_mv.getNumVectors () <<
" != 1 columns.");
87 template<
class Scalar,
class LO,
class GO,
class Node>
95 template<
class Scalar,
class LO,
class GO,
class Node>
104 template<
class Scalar,
class LO,
class GO,
class Node>
112 template<
class Scalar,
class LO,
class GO,
class Node>
115 Teuchos::RCP<vec_type>
vPtr = this->mv_.getVectorNonConst (0);
119 template<
class Scalar,
class LO,
class GO,
class Node>
126 template<
class Scalar,
class LO,
class GO,
class Node>
133 template<
class Scalar,
class LO,
class GO,
class Node>
140 template<
class Scalar,
class LO,
class GO,
class Node>
147#ifdef TPETRA_ENABLE_DEPRECATED_CODE
148 template<
class Scalar,
class LO,
class GO,
class Node>
156 template<
class Scalar,
class LO,
class GO,
class Node>
159 BlockVector<Scalar, LO, GO, Node>::
160 getGlobalRowView (
const GO globalRowIndex, Scalar*& vals) {
161 return ((base_type*)
this)->getGlobalRowView (globalRowIndex, 0, vals);
164 template<
class Scalar,
class LO,
class GO,
class Node>
165 typename BlockVector<Scalar, LO, GO, Node>::little_host_vec_type
167 BlockVector<Scalar, LO, GO, Node>::
168 getLocalBlock (
const LO localRowIndex)
170 if (! this->isValidLocalMeshIndex (localRowIndex)) {
171 return little_host_vec_type ();
174 const size_t blockSize = this->getBlockSize ();
175 const size_t offset = localRowIndex * blockSize;
176 return little_host_vec_type (this->getRawPtr () + offset, blockSize);
181 template<
class Scalar,
class LO,
class GO,
class Node>
182 typename BlockVector<Scalar, LO, GO, Node>::const_little_host_vec_type
190 template<
class Scalar,
class LO,
class GO,
class Node>
191 typename BlockVector<Scalar, LO, GO, Node>::little_host_vec_type
195 return ((base_type*)
this)->getLocalBlockHost(
localRowIndex, 0,
199 template<
class Scalar,
class LO,
class GO,
class Node>
200 typename BlockVector<Scalar, LO, GO, Node>::little_host_vec_type
204 return ((base_type*)
this)->getLocalBlockHost(localRowIndex, 0,
205 Access::OverwriteAll);
215#define TPETRA_BLOCKVECTOR_INSTANT(S,LO,GO,NODE) \
216 template class BlockVector< S, LO, GO, NODE >;
Vector for multiple degrees of freedom per mesh point.
const_little_host_vec_type getLocalBlockHost(const LO localRowIndex, Access::ReadOnlyStruct) const
Get a view of the degrees of freedom at the given mesh point, using a local index.
bool sumIntoGlobalValues(const GO globalRowIndex, const Scalar vals[])
Sum into all values at the given mesh point, using a global index.
bool replaceGlobalValues(const GO globalRowIndex, const Scalar vals[])
Replace all values at the given mesh point, using a global index.
BlockVector()
Default constructor.
bool sumIntoLocalValues(const LO localRowIndex, const Scalar vals[])
Sum into 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 replaceLocalValues(const LO localRowIndex, const Scalar vals[])
Replace all values at the given mesh point, using a local index.
Struct that holds views of the contents of a CrsMatrix.
A distributed dense vector.
Namespace Tpetra contains the class and methods constituting the Tpetra library.