10#ifndef TPETRA_BLOCKMULTIVECTOR_DECL_HPP
11#define TPETRA_BLOCKMULTIVECTOR_DECL_HPP
15#include "Tpetra_MultiVector.hpp"
105template <
class Scalar,
112 using STS = Teuchos::ScalarTraits<Scalar>;
161 typedef typename little_vec_type::host_mirror_type little_host_vec_type;
168 typedef Kokkos::View<const impl_scalar_type*, device_type>
173 typedef Kokkos::View<const impl_scalar_type*, host_device_type>
174 const_little_host_vec_type;
300 static Teuchos::RCP<const map_type>
416 template <
class TargetMemorySpace>
490 const_little_host_vec_type getLocalBlockHost(
493 const Access::ReadOnlyStruct)
const;
495 little_host_vec_type getLocalBlockHost(
498 const Access::ReadWriteStruct);
503 little_host_vec_type getLocalBlockHost(
506 const Access::OverwriteAllStruct);
519 using dist_object_type
::
525 const size_t numSameIDs,
541 Kokkos::DualView<packet_type*,
543 Kokkos::DualView<
size_t*,
556 Kokkos::DualView<packet_type*,
559 Kokkos::DualView<
size_t*,
570 return static_cast<size_t>(1);
591 Teuchos::RCP<const map_type> pointMap_;
612 static Teuchos::RCP<const mv_type>
615 static Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
Forward declaration of Tpetra::BlockCrsMatrix.
Forward declaration of Tpetra::BlockMultiVector.
MultiVector for multiple degrees of freedom per mesh point.
virtual bool checkSizes(const Tpetra::SrcDistObject &source) override
Compare the source and target (this) objects for compatibility.
void putScalar(const Scalar &val)
Fill all entries with the given value val.
const mv_type & getMultiVectorView() const
Get a Tpetra::MultiVector that views this BlockMultiVector's data.
BlockMultiVector(BlockMultiVector< Scalar, LO, GO, Node > &&)=default
Move constructor (shallow move).
void blockWiseMultiply(const Scalar &alpha, const Kokkos::View< const impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &D, const BlockMultiVector< Scalar, LO, GO, Node > &X)
*this := alpha * D * X, where D is a block diagonal matrix.
BlockMultiVector< Scalar, LO, GO, Node > & operator=(BlockMultiVector< Scalar, LO, GO, Node > &&)=default
Move assigment (shallow move).
bool need_sync_host() const
Whether this object needs synchronization to the host.
bool sumIntoLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[])
Sum into all values at the given mesh point, using a local index.
bool need_sync() const
Whether this object needs synchronization to the given memory space.
GO global_ordinal_type
The type of global indices.
static Teuchos::RCP< const map_type > makePointMapRCP(const map_type &meshMap, const LO blockSize)
Create and return an owning RCP to the point Map corresponding to the given mesh Map and block size.
Kokkos::View< impl_scalar_type *, device_type > little_vec_type
"Block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector.
map_type meshMap_
Mesh Map given to constructor.
void blockJacobiUpdate(const Scalar &alpha, const Kokkos::View< const impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &D, const BlockMultiVector< Scalar, LO, GO, Node > &X, BlockMultiVector< Scalar, LO, GO, Node > &Z, const Scalar &beta)
Block Jacobi update .
BlockMultiVector(const BlockMultiVector< Scalar, LO, GO, Node > &)=default
Copy constructor (shallow copy).
BlockMultiVector< Scalar, LO, GO, Node > & operator=(const BlockMultiVector< Scalar, LO, GO, Node > &)=default
Copy assigment (shallow copy).
void update(const Scalar &alpha, const BlockMultiVector< Scalar, LO, GO, Node > &X, const Scalar &beta)
Update: this = beta*this + alpha*X.
typename mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the object.
bool sumIntoGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[])
Sum into all values at the given mesh point, using a global index.
BlockMultiVector()
Default constructor.
const map_type getPointMap() const
Get this BlockMultiVector's (previously computed) point Map.
bool replaceLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[])
Replace all values at the given mesh point, using local row and column indices.
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.
typename dist_object_type::buffer_device_type buffer_device_type
Kokkos::Device specialization used for communication buffers.
LO getNumVectors() const
Get the number of columns (vectors) in the BlockMultiVector.
Node node_type
The Kokkos Node type; a legacy thing that will go away at some point.
Tpetra::MultiVector< Scalar, LO, GO, Node > mv_type
The specialization of Tpetra::MultiVector that this class uses.
virtual void copyAndPermute(const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM) override
virtual void packAndPrepare(const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets) override
typename mv_type::device_type device_type
The Kokkos Device type.
size_t getStrideX() const
Stride between consecutive local entries in the same column.
void scale(const Scalar &val)
Multiply all entries in place by the given value val.
LO local_ordinal_type
The type of local indices.
Kokkos::View< const impl_scalar_type *, device_type > const_little_vec_type
"Const block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector.
virtual void unpackAndCombine(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode) override
LO getBlockSize() const
Get the number of degrees of freedom per mesh point.
bool need_sync_device() const
Whether this object needs synchronization to the device.
bool replaceGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[])
Replace all values at the given mesh point, using a global index.
size_t getStrideY() const
Stride between consecutive local entries in the same row.
bool isValidLocalMeshIndex(const LO meshLocalIndex) const
True if and only if meshLocalIndex is a valid local index in the mesh Map.
mv_type mv_
The Tpetra::MultiVector used to represent the data.
Struct that holds views of the contents of a CrsMatrix.
Base class for distributed Tpetra objects that support data redistribution.
virtual void packAndPrepare(const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets)
Pack data and metadata for communication (sends).
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
Kokkos::Device specialization for communication buffers.
typename ::Kokkos::ArithTraits< Scalar >::val_type packet_type
The type of each datum being sent or received in an Import or Export.
virtual void unpackAndCombine(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode)
Perform any unpacking and combining after communication.
size_t getStride() const
Stride between columns in the multivector.
bool need_sync_device() const
Whether this MultiVector needs synchronization to the device.
bool need_sync_host() const
Whether this MultiVector needs synchronization to the host.
typename map_type::device_type device_type
This class' preferred Kokkos device type.
size_t getNumVectors() const
Number of columns in the multivector.
typename Kokkos::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
Abstract base class for objects that can be the source of an Import or Export operation.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
CombineMode
Rule for combining data in an Import or Export.