42 #ifndef TPETRA_EXPERIMENTAL_BLOCKMULTIVECTOR_DECL_HPP 43 #define TPETRA_EXPERIMENTAL_BLOCKMULTIVECTOR_DECL_HPP 45 #include <Tpetra_MultiVector.hpp> 49 namespace Experimental {
51 #ifndef DOXYGEN_SHOULD_SKIP_THIS 53 template<
class S,
class LO,
class GO,
class N>
class BlockCrsMatrix;
54 #endif // DOXYGEN_SHOULD_SKIP_THIS 146 class GO = Details::DefaultTypes::global_ordinal_type,
153 typedef Teuchos::ScalarTraits<Scalar> STS;
156 typedef Scalar packet_type;
250 const map_type& pointMap,
267 const map_type& meshMap,
287 makePointMap (
const map_type& meshMap,
const LO blockSize);
322 void scale (
const Scalar& val);
345 bool replaceLocalValues (
const LO localRowIndex,
const LO colIndex,
const Scalar vals[])
const;
357 bool replaceGlobalValues (
const GO globalRowIndex,
const LO colIndex,
const Scalar vals[])
const;
369 bool sumIntoLocalValues (
const LO localRowIndex,
const LO colIndex,
const Scalar vals[])
const;
381 bool sumIntoGlobalValues (
const GO globalRowIndex,
const LO colIndex,
const Scalar vals[])
const;
393 bool getLocalRowView (
const LO localRowIndex,
const LO colIndex, Scalar*& vals)
const;
405 bool getGlobalRowView (
const GO globalRowIndex,
const LO colIndex, Scalar*& vals)
const;
419 getLocalBlock (
const LO localRowIndex,
const LO colIndex)
const;
436 const Teuchos::ArrayView<const LO>& permuteToLIDs,
437 const Teuchos::ArrayView<const LO>& permuteFromLIDs);
441 const Teuchos::ArrayView<const LO>& exportLIDs,
442 Teuchos::Array<impl_scalar_type>& exports,
443 const Teuchos::ArrayView<size_t>& numPacketsPerLID,
444 size_t& constantNumPackets,
449 const Teuchos::ArrayView<const impl_scalar_type> &imports,
450 const Teuchos::ArrayView<size_t> &numPacketsPerLID,
451 size_t constantNumPackets,
464 return static_cast<size_t> (1);
475 return meshLocalIndex != Teuchos::OrdinalTraits<LO>::invalid () &&
481 template<
class Scalar2,
class LO2,
class GO2,
class Node2>
507 impl_scalar_type* mvData_;
514 replaceLocalValuesImpl (
const LO localRowIndex,
516 const Scalar vals[])
const;
519 sumIntoLocalValuesImpl (
const LO localRowIndex,
521 const Scalar vals[])
const;
523 static Teuchos::RCP<const mv_type>
526 static Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
533 #endif // TPETRA_EXPERIMENTAL_BLOCKMULTIVECTOR_DECL_HPP Namespace Tpetra contains the class and methods constituting the Tpetra library.
LO local_ordinal_type
The type of local indices.
mv_type mv_
The Tpetra::MultiVector used to represent the data.
KokkosClassic::DefaultNode::DefaultNodeType node_type
Default value of Node template parameter.
void putScalar(const Scalar &val)
Fill all entries with the given value val.
mv_type getMultiVectorView()
Get a Tpetra::MultiVector that views this BlockMultiVector's data.
Tpetra::Map< LO, GO, Node > map_type
The specialization of Tpetra::Map that this class uses.
bool replaceGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a global index.
bool sumIntoGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a global index.
MultiVector for multiple degrees of freedom per mesh point.
static map_type makePointMap(const map_type &meshMap, const LO blockSize)
Create and return the point Map corresponding to the given mesh Map and block size.
Node node_type
The Kokkos Node type.
virtual void unpackAndCombine(const Teuchos::ArrayView< const LO > &importLIDs, const Teuchos::ArrayView< const impl_scalar_type > &imports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t constantNumPackets, Tpetra::Distributor &distor, Tpetra::CombineMode CM)
Perform any unpacking and combining after communication.
int local_ordinal_type
Default value of LocalOrdinal template parameter.
void scale(const Scalar &val)
Multiply all entries in place by the given value val.
mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the matrix.
little_vec_type getLocalBlock(const LO localRowIndex, const LO colIndex) const
Get a view of the degrees of freedom at the given mesh point.
bool isNodeLocalElement(LocalOrdinal localIndex) const
Whether the given local index is valid for this Map on the calling process.
LO getBlockSize() const
Get the number of degrees of freedom per mesh point.
impl_scalar_type * getRawPtr() const
Raw pointer to the MultiVector's data.
Declaration and definition of LittleBlock and LittleVector.
Nonowning view of a set of degrees of freedom corresponding to a mesh point in a block vector or mult...
LO getNumVectors() const
Get the number of columns (vectors) in the BlockMultiVector.
Sets up and executes a communication plan for a Tpetra DistObject.
CombineMode
Rule for combining data in an Import or Export.
bool getLocalRowView(const LO localRowIndex, const LO colIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a local index. ...
Abstract base class for objects that can be the source of an Import or Export operation.
double scalar_type
Default value of Scalar template parameter.
bool isValidLocalMeshIndex(const LO meshLocalIndex) const
True if and only if meshLocalIndex is a valid local index in the mesh Map.
size_t getStride() const
Stride between columns in the multivector.
Kokkos::Details::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
map_type getPointMap() const
Get this BlockMultiVector's (previously computed) point Map.
virtual bool checkSizes(const Tpetra::SrcDistObject &source)
Compare the source and target (this) objects for compatibility.
GO global_ordinal_type
The type of global indices.
LittleVector< impl_scalar_type, LO > little_vec_type
"Block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector...
Scalar scalar_type
The type of entries in the matrix.
size_t getStrideY() const
Stride between consecutive local entries in the same row.
LittleVector< const impl_scalar_type, LO > const_little_vec_type
"Const block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector...
Constant block CRS matrix class.
BlockMultiVector()
Default constructor.
virtual void copyAndPermute(const Tpetra::SrcDistObject &source, size_t numSameIDs, const Teuchos::ArrayView< const LO > &permuteToLIDs, const Teuchos::ArrayView< const LO > &permuteFromLIDs)
Perform copies and permutations that are local to this process.
size_t getNumVectors() const
Number of columns in the multivector.
bool replaceLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using local row and column indices.
size_t getStrideX() const
Stride between consecutive local entries in the same column.
virtual void packAndPrepare(const Tpetra::SrcDistObject &source, const Teuchos::ArrayView< const LO > &exportLIDs, Teuchos::Array< impl_scalar_type > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets, Tpetra::Distributor &distor)
Perform any packing or preparation required for communication.
Base class for distributed Tpetra objects that support data redistribution.
bool sumIntoLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a local index.
Tpetra::MultiVector< Scalar, LO, GO, Node > mv_type
The specialization of Tpetra::MultiVector that this class uses.
bool getGlobalRowView(const GO globalRowIndex, const LO colIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a global index.