42 #ifndef TPETRA_BLOCKMULTIVECTOR_DECL_HPP
43 #define TPETRA_BLOCKMULTIVECTOR_DECL_HPP
47 #include "Tpetra_MultiVector.hpp"
141 template<
class Scalar,
150 using STS = Teuchos::ScalarTraits<Scalar>;
201 Kokkos::MemoryTraits<Kokkos::Unmanaged> >
212 Kokkos::MemoryTraits<Kokkos::Unmanaged> >
241 const Teuchos::DataAccess copyOrView);
309 const size_t offset = 0);
317 const size_t offset = 0);
365 void scale (
const Scalar& val);
374 update (
const Scalar& alpha,
465 template<
class TargetMemorySpace>
467 mv_.template sync<typename TargetMemorySpace::memory_space> ();
481 template<
class TargetMemorySpace>
483 return mv_.template need_sync<typename TargetMemorySpace::memory_space> ();
501 template<
class TargetMemorySpace>
503 mv_.template modify<typename TargetMemorySpace::memory_space> ();
536 bool replaceLocalValues (
const LO localRowIndex,
const LO colIndex,
const Scalar vals[])
const;
548 bool replaceGlobalValues (
const GO globalRowIndex,
const LO colIndex,
const Scalar vals[])
const;
560 bool sumIntoLocalValues (
const LO localRowIndex,
const LO colIndex,
const Scalar vals[])
const;
572 bool sumIntoGlobalValues (
const GO globalRowIndex,
const LO colIndex,
const Scalar vals[])
const;
584 bool getLocalRowView (
const LO localRowIndex,
const LO colIndex, Scalar*& vals)
const;
596 bool getGlobalRowView (
const GO globalRowIndex,
const LO colIndex, Scalar*& vals)
const;
609 typename little_vec_type::HostMirror
610 getLocalBlock (
const LO localRowIndex,
const LO colIndex)
const;
624 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
630 const size_t numSameIDs,
637 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
645 Kokkos::DualView<packet_type*,
647 Kokkos::DualView<
size_t*,
649 size_t& constantNumPackets,
653 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
660 Kokkos::DualView<packet_type*,
662 Kokkos::DualView<
size_t*,
664 const size_t constantNumPackets,
678 return static_cast<size_t> (1);
721 replaceLocalValuesImpl (
const LO localRowIndex,
723 const Scalar vals[])
const;
726 sumIntoLocalValuesImpl (
const LO localRowIndex,
728 const Scalar vals[])
const;
730 static Teuchos::RCP<const mv_type>
733 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.
BlockMultiVector< Scalar, LO, GO, Node > & operator=(const BlockMultiVector< Scalar, LO, GO, Node > &)=default
Copy assigment (shallow copy).
virtual bool checkSizes(const Tpetra::SrcDistObject &source)
Compare the source and target (this) objects for compatibility.
void putScalar(const Scalar &val)
Fill all entries with the given value val.
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.
bool need_sync_host() const
Whether this object needs synchronization to the host.
bool need_sync() const
Whether this object needs synchronization to the given memory space.
GO global_ordinal_type
The type of global indices.
void sync_device()
Update data to the device.
impl_scalar_type * getRawPtr() const
Raw pointer to the MultiVector's data.
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).
Kokkos::View< const impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > const_little_vec_type
"Const block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector.
little_vec_type::HostMirror getLocalBlock(const LO localRowIndex, const LO colIndex) const
Get a host view of the degrees of freedom at the given mesh point.
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, Distributor &distor, const CombineMode combineMode)
Perform any unpacking and combining after communication.
bool replaceGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a global index.
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 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.
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, Distributor &distor)
Pack data and metadata for communication (sends).
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.
void sync_host()
Update data to the host.
BlockMultiVector()
Default constructor.
mv_type getMultiVectorView() const
Get a Tpetra::MultiVector that views this BlockMultiVector's data.
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.
void modify_host()
Mark data as modified on the host.
LO getNumVectors() const
Get the number of columns (vectors) in the BlockMultiVector.
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.
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.
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::Map< LO, GO, Node > map_type
The specialization of Tpetra::Map that this class uses.
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.
void modify_device()
Mark data as modified on the device.
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.
map_type getPointMap() const
Get this BlockMultiVector's (previously computed) point Map.
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)
Perform copies and permutations that are local to the calling (MPI) process.
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.
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.
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.
void modify()
Mark data as modified on the given memory space.
Scalar scalar_type
The type of entries in the object.
void sync()
Update data to the given target memory space, only if data in the "other" space have been marked as m...
mv_type mv_
The Tpetra::MultiVector used to represent the data.
Base class for distributed Tpetra objects that support data redistribution.
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
Kokkos::Device specialization for communication buffers.
typename ::Kokkos::Details::ArithTraits< Scalar >::val_type packet_type
The type of each datum being sent or received in an Import or Export.
Sets up and executes a communication plan for a Tpetra DistObject.
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.
void modify_device()
Mark data as modified on the device side.
typename Kokkos::Details::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
typename map_type::device_type device_type
This class' preferred Kokkos device type.
size_t getNumVectors() const
Number of columns in the multivector.
void sync_host()
Synchronize to Host.
void sync_device()
Synchronize to Device.
void modify_host()
Mark data as modified on the host side.
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.