40 #ifndef TPETRA_BLOCKMULTIVECTOR_DECL_HPP
41 #define TPETRA_BLOCKMULTIVECTOR_DECL_HPP
45 #include "Tpetra_MultiVector.hpp"
139 template<
class Scalar,
148 using STS = Teuchos::ScalarTraits<Scalar>;
199 Kokkos::MemoryTraits<Kokkos::Unmanaged> >
210 Kokkos::MemoryTraits<Kokkos::Unmanaged> >
239 const Teuchos::DataAccess copyOrView);
307 const size_t offset = 0);
315 const size_t offset = 0);
363 void scale (
const Scalar& val);
372 update (
const Scalar& alpha,
463 template<
class TargetMemorySpace>
465 mv_.template sync<typename TargetMemorySpace::memory_space> ();
479 template<
class TargetMemorySpace>
481 return mv_.template need_sync<typename TargetMemorySpace::memory_space> ();
499 template<
class TargetMemorySpace>
501 mv_.template modify<typename TargetMemorySpace::memory_space> ();
534 bool replaceLocalValues (
const LO localRowIndex,
const LO colIndex,
const Scalar vals[])
const;
546 bool replaceGlobalValues (
const GO globalRowIndex,
const LO colIndex,
const Scalar vals[])
const;
558 bool sumIntoLocalValues (
const LO localRowIndex,
const LO colIndex,
const Scalar vals[])
const;
570 bool sumIntoGlobalValues (
const GO globalRowIndex,
const LO colIndex,
const Scalar vals[])
const;
582 bool getLocalRowView (
const LO localRowIndex,
const LO colIndex, Scalar*& vals)
const;
594 bool getGlobalRowView (
const GO globalRowIndex,
const LO colIndex, Scalar*& vals)
const;
607 typename little_vec_type::HostMirror
608 getLocalBlock (
const LO localRowIndex,
const LO colIndex)
const;
624 const size_t numSameIDs,
635 Kokkos::DualView<packet_type*,
637 Kokkos::DualView<
size_t*,
639 size_t& constantNumPackets,
646 Kokkos::DualView<packet_type*,
648 Kokkos::DualView<
size_t*,
650 const size_t constantNumPackets,
664 return static_cast<size_t> (1);
707 replaceLocalValuesImpl (
const LO localRowIndex,
709 const Scalar vals[])
const;
712 sumIntoLocalValuesImpl (
const LO localRowIndex,
714 const Scalar vals[])
const;
716 static Teuchos::RCP<const mv_type>
719 static Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
725 #endif // TPETRA_BLOCKMULTIVECTOR_DECL_HPP
Node node_type
The Kokkos Node type; a legacy thing that will go away at some point.
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.
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.
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 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()
Update data to the given target memory space, only if data in the "other" space have been marked as m...
void putScalar(const Scalar &val)
Fill all entries with the given value val.
bool need_sync_device() const
Whether this object needs synchronization to the device.
LO local_ordinal_type
The type of local indices.
size_t getNumVectors() const
Number of columns in the multivector.
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 .
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.
size_t getStrideX() const
Stride between consecutive local entries in the same column.
Forward declaration of Tpetra::BlockCrsMatrix.
GO global_ordinal_type
The type of global indices.
virtual bool checkSizes(const Tpetra::SrcDistObject &source)
Compare the source and target (this) objects for compatibility.
void modify_device()
Mark data as modified on the device.
size_t getStrideY() const
Stride between consecutive local entries in the same row.
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.
bool need_sync_device() const
Whether this MultiVector needs synchronization to the device.
MultiVector for multiple degrees of freedom per mesh point.
bool need_sync_host() const
Whether this object needs synchronization to the host.
BlockMultiVector()
Default constructor.
LO getBlockSize() const
Get the number of degrees of freedom per mesh point.
typename::Kokkos::Details::ArithTraits< Scalar >::val_type packet_type
The type of each datum being sent or received in an Import or Export.
void sync_host()
Synchronize to Host.
bool isValidLocalMeshIndex(const LO meshLocalIndex) const
True if and only if meshLocalIndex is a valid local index in the mesh Map.
void modify()
Mark data as modified on the given memory space.
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.
Sets up and executes a communication plan for a Tpetra DistObject.
typename mv_type::device_type device_type
The Kokkos Device type.
CombineMode
Rule for combining data in an Import or Export.
map_type getPointMap() const
Get this BlockMultiVector's (previously computed) point Map.
void modify_host()
Mark data as modified on the host.
Tpetra::MultiVector< Scalar, LO, GO, Node > mv_type
The specialization of Tpetra::MultiVector that this class uses.
Tpetra::Map< LO, GO, Node > map_type
The specialization of Tpetra::Map that this class uses.
bool need_sync_host() const
Whether this MultiVector needs synchronization to the host.
void modify_device()
Mark data as modified on the device side.
mv_type mv_
The Tpetra::MultiVector used to represent the data.
Abstract base class for objects that can be the source of an Import or Export operation.
Forward declaration of Tpetra::BlockMultiVector.
void modify_host()
Mark data as modified on the host side.
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...
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.
Scalar scalar_type
The type of entries in the object.
typename map_type::device_type device_type
This class' preferred Kokkos device type.
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 need_sync() const
Whether this object needs synchronization to the given memory space.
typename Kokkos::Details::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
LO getNumVectors() const
Get the number of columns (vectors) in the BlockMultiVector.
mv_type getMultiVectorView() const
Get a Tpetra::MultiVector that views this BlockMultiVector's data.
typename dist_object_type::buffer_device_type buffer_device_type
Kokkos::Device specialization used for communication buffers.
size_t getStride() const
Stride between columns in the multivector.
impl_scalar_type * getRawPtr() const
Raw pointer to the MultiVector's data.
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
void scale(const Scalar &val)
Multiply all entries in place by the given value val.
Base class for distributed Tpetra objects that support data redistribution.
typename mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the object.
void sync_device()
Update data to the device.
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.
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...
void sync_host()
Update data to the host.
void sync_device()
Synchronize to Device.