11#ifndef TPETRA_BLOCKMULTIVECTOR_DECL_HPP
12#define TPETRA_BLOCKMULTIVECTOR_DECL_HPP
16#include "Tpetra_MultiVector.hpp"
106template<
class Scalar,
115 using STS = Teuchos::ScalarTraits<Scalar>;
164 typedef typename little_vec_type::HostMirror little_host_vec_type;
171 typedef Kokkos::View<const impl_scalar_type *, device_type>
176 typedef Kokkos::View<const impl_scalar_type*, host_device_type>
177 const_little_host_vec_type;
205 const Teuchos::DataAccess copyOrView);
273 const size_t offset = 0);
281 const size_t offset = 0);
303 static Teuchos::RCP<const map_type>
321 return static_cast<LO
> (
mv_.getNumVectors ());
420 template<
class TargetMemorySpace>
427 return mv_.need_sync_host();
432 return mv_.need_sync_device();
495 const_little_host_vec_type getLocalBlockHost(
496 const LO localRowIndex,
498 const Access::ReadOnlyStruct)
const;
500 little_host_vec_type getLocalBlockHost(
501 const LO localRowIndex,
503 const Access::ReadWriteStruct);
509 const LO localRowIndex,
511 const Access::OverwriteAllStruct);
525 using dist_object_type::
533 const size_t numSameIDs,
551 Kokkos::DualView<packet_type*,
553 Kokkos::DualView<
size_t*,
555 size_t& constantNumPackets)
override;
567 Kokkos::DualView<packet_type*,
569 Kokkos::DualView<
size_t*,
571 const size_t constantNumPackets,
581 return static_cast<size_t> (1);
586 return mv_.getStride ();
602 Teuchos::RCP<const map_type> pointMap_;
615 replaceLocalValuesImpl (
const LO localRowIndex,
617 const Scalar vals[]);
620 sumIntoLocalValuesImpl (
const LO localRowIndex,
622 const Scalar vals[]);
624 static Teuchos::RCP<const mv_type>
627 static Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
Forward declaration of Tpetra::BlockCrsMatrix.
Forward declaration of Tpetra::BlockMultiVector.
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.
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.
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.
BlockMultiVector(const map_type &meshMap, const LO blockSize, const LO numVecs)
Constructor that takes a mesh Map, a block size, and a number of vectors (columns).
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).
little_host_vec_type getLocalBlockHost(const LO localRowIndex, const LO colIndex, const Access::OverwriteAllStruct)
Get a local block on host, with the intent to overwrite all blocks in the BlockMultiVector before acc...
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
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(const BlockMultiVector< Scalar, LO, GO, Node > &X, const map_type &newMeshMap, const size_t offset=0)
View an existing BlockMultiVector using a different mesh Map; compute the new point Map.
BlockMultiVector(const BlockMultiVector< Scalar, LO, GO, Node > &in, const Teuchos::DataAccess copyOrView)
"Copy constructor" with option to deep copy.
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
LO getNumVectors() const
Get the number of columns (vectors) in the BlockMultiVector.
Tpetra::MultiVector< Scalar, LO, GO, Node > mv_type
BlockMultiVector(const BlockMultiVector< Scalar, LO, GO, Node > &X, const map_type &newMeshMap, const map_type &newPointMap, const size_t offset=0)
View an existing BlockMultiVector using a different mesh Map, supplying the corresponding 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, const CombineMode CM) override
Kokkos::View< impl_scalar_type *, device_type > little_vec_type
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
Tpetra::Map< LO, GO, Node > map_type
typename mv_type::device_type 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.
BlockMultiVector(const mv_type &X_mv, const map_type &meshMap, const LO blockSize)
View an existing MultiVector.
Kokkos::View< const impl_scalar_type *, device_type > const_little_vec_type
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.
BlockMultiVector(const map_type &meshMap, const map_type &pointMap, const LO blockSize, const LO numVecs)
Constructor that takes a mesh Map, a point Map, a block size, and a number of vectors (columns).
Base class for distributed Tpetra objects that support data redistribution.
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
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)
typename ::Kokkos::ArithTraits< Scalar >::val_type packet_type
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)
A parallel distribution of indices over processes.
One or more distributed dense vectors.
typename map_type::device_type device_type
typename Kokkos::ArithTraits< Scalar >::val_type impl_scalar_type
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.