Tpetra parallel linear algebra
Version of the Day
|
MultiVector for multiple degrees of freedom per mesh point. More...
#include <Tpetra_Experimental_BlockMultiVector_decl.hpp>
Public Types | |
Typedefs to facilitate template metaprogramming. | |
using | map_type = Tpetra::Map< LO, GO, Node > |
The specialization of Tpetra::Map that this class uses. More... | |
using | mv_type = Tpetra::MultiVector< Scalar, LO, GO, Node > |
The specialization of Tpetra::MultiVector that this class uses. More... | |
using | scalar_type = Scalar |
The type of entries in the object. More... | |
using | impl_scalar_type = typename mv_type::impl_scalar_type |
The implementation type of entries in the object. More... | |
using | local_ordinal_type = LO |
The type of local indices. More... | |
using | global_ordinal_type = GO |
The type of global indices. More... | |
using | device_type = typename mv_type::device_type |
The Kokkos Device type. More... | |
using | node_type = Node |
The Kokkos Node type; a legacy thing that will go away at some point. More... | |
using | buffer_device_type = typename dist_object_type::buffer_device_type |
Kokkos::Device specialization used for communication buffers. More... | |
typedef 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. More... | |
typedef 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. More... | |
Typedefs | |
using | execution_space = typename device_type::execution_space |
The Kokkos execution space. More... | |
Public Member Functions | |
Constructors | |
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). More... | |
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). More... | |
BlockMultiVector (const mv_type &X_mv, const map_type &meshMap, const LO blockSize) | |
View an existing MultiVector. More... | |
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. More... | |
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. More... | |
BlockMultiVector () | |
Default constructor. More... | |
Coarse-grained operations | |
void | putScalar (const Scalar &val) |
Fill all entries with the given value val . More... | |
void | scale (const Scalar &val) |
Multiply all entries in place by the given value val . More... | |
void | update (const Scalar &alpha, const BlockMultiVector< Scalar, LO, GO, Node > &X, const Scalar &beta) |
Update: this = beta*this + alpha*X . More... | |
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. More... | |
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 . More... | |
Implementation of "dual view semantics" | |
template<class TargetMemorySpace > | |
void | sync () |
Update data to the given target memory space, only if data in the "other" space have been marked as modified. More... | |
void | sync_host () |
Update data to the host. More... | |
void | sync_device () |
Update data to the device. More... | |
template<class TargetMemorySpace > | |
bool | need_sync () const |
Whether this object needs synchronization to the given memory space. More... | |
bool | need_sync_host () const |
Whether this object needs synchronization to the host. More... | |
bool | need_sync_device () const |
Whether this object needs synchronization to the device. More... | |
template<class TargetMemorySpace > | |
void | modify () |
Mark data as modified on the given memory space. More... | |
void | modify_host () |
Mark data as modified on the host. More... | |
void | modify_device () |
Mark data as modified on the device. More... | |
Fine-grained data access | |
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. More... | |
bool | replaceGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const |
Replace all values at the given mesh point, using a global index. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
Public methods for redistributing data | |
void | doImport (const SrcDistObject &source, const Import< LO, GO, Node > &importer, const CombineMode CM, const bool restrictedMode=false) |
Import data into this object using an Import object ("forward mode"). More... | |
void | doImport (const SrcDistObject &source, const Export< LO, GO, Node > &exporter, const CombineMode CM, const bool restrictedMode=false) |
Import data into this object using an Export object ("reverse mode"). More... | |
void | doExport (const SrcDistObject &source, const Export< LO, GO, Node > &exporter, const CombineMode CM, const bool restrictedMode=false) |
Export data into this object using an Export object ("forward mode"). More... | |
void | doExport (const SrcDistObject &source, const Import< LO, GO, Node > &importer, const CombineMode CM, const bool restrictedMode=false) |
Export data into this object using an Import object ("reverse mode"). More... | |
Attribute accessor methods | |
bool | isDistributed () const |
Whether this is a globally distributed object. More... | |
virtual Teuchos::RCP< const map_type > | getMap () const |
The Map describing the parallel distribution of this object. More... | |
I/O methods | |
void | print (std::ostream &os) const |
Print this object to the given output stream. More... | |
Implementation of Teuchos::Describable | |
virtual std::string | description () const |
One-line descriptiion of this object. More... | |
virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const |
Print a descriptiion of this object to the given output stream. More... | |
Methods for use only by experts | |
virtual void | removeEmptyProcessesInPlace (const Teuchos::RCP< const map_type > &newMap) |
Remove processes which contain no entries in this object's Map. More... | |
Protected Member Functions | |
impl_scalar_type * | getRawPtr () const |
Raw pointer to the MultiVector's data. More... | |
size_t | getStrideX () const |
Stride between consecutive local entries in the same column. More... | |
size_t | getStrideY () const |
Stride between consecutive local entries in the same row. More... | |
bool | isValidLocalMeshIndex (const LO meshLocalIndex) const |
True if and only if meshLocalIndex is a valid local index in the mesh Map. More... | |
virtual size_t | constantNumberOfPackets () const |
Whether the implementation's instance promises always to have a constant number of packets per LID (local index), and if so, how many packets per LID there are. More... | |
virtual void | doTransfer (const SrcDistObject &src, const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &transfer, const char modeString[], const ReverseOption revOp, const CombineMode CM, const bool restrictedMode) |
Redistribute data across (MPI) processes. More... | |
virtual bool | reallocArraysForNumPacketsPerLid (const size_t numExportLIDs, const size_t numImportLIDs) |
Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary. More... | |
virtual void | doTransferNew (const SrcDistObject &src, const CombineMode CM, 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 Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &remoteLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Distributor &distor, const ReverseOption revOp, const bool commOnHost, const bool restrictedMode) |
Implementation detail of doTransfer. More... | |
Implementation of Tpetra::DistObject. | |
The methods here implement Tpetra::DistObject. They let BlockMultiVector participate in Import and Export operations. Users don't have to worry about these methods. | |
virtual bool | checkSizes (const Tpetra::SrcDistObject &source) |
Compare the source and target (this) objects for compatibility. More... | |
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) |
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) |
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) |
Protected Attributes | |
mv_type | mv_ |
The Tpetra::MultiVector used to represent the data. More... | |
Methods implemented by subclasses and used by doTransfer(). | |
The doTransfer() method uses the subclass' implementations of these methods to implement data transfer. Subclasses of DistObject must implement these methods. This is an instance of the Template Method Pattern. ("Template" here doesn't mean "C++ template"; it means "pattern with holes that are filled in by the subclass' method implementations.") | |
Teuchos::RCP< const map_type > | map_ |
The Map over which this object is distributed. More... | |
Kokkos::DualView< packet_type *, buffer_device_type > | imports_ |
Buffer into which packed data are imported (received from other processes). More... | |
Kokkos::DualView< size_t *, buffer_device_type > | numImportPacketsPerLID_ |
Number of packets to receive for each receive operation. More... | |
Kokkos::DualView< packet_type *, buffer_device_type > | exports_ |
Buffer from which packed data are exported (sent to other processes). More... | |
Kokkos::DualView< size_t *, buffer_device_type > | numExportPacketsPerLID_ |
Number of packets to send for each send operation. More... | |
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. More... | |
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). More... | |
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. More... | |
bool | reallocImportsIfNeeded (const size_t newSize, const bool verbose, const std::string *prefix) |
Reallocate imports_ if needed. More... | |
Access to Maps, the block size, and a MultiVector view. | |
map_type | getPointMap () const |
Get this BlockMultiVector's (previously computed) point Map. More... | |
LO | getBlockSize () const |
Get the number of degrees of freedom per mesh point. More... | |
LO | getNumVectors () const |
Get the number of columns (vectors) in the BlockMultiVector. More... | |
mv_type | getMultiVectorView () const |
Get a Tpetra::MultiVector that views this BlockMultiVector's data. More... | |
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. More... | |
MultiVector for multiple degrees of freedom per mesh point.
Scalar | The type of each entry of the block multivector. (You can use real-valued or complex-valued types here, unlike in Epetra, where the scalar type is always double .) |
LO | The type of local indices. See the documentation of the first template parameter of Map for requirements. |
GO | The type of global indices. See the documentation of the second template parameter of Map for requirements. |
Node | The Kokkos Node type. See the documentation of the third template parameter of Map for requirements. |
BlockMultiVector is like Tpetra::MultiVector, but its interface supports multiple degrees of freedom per mesh point. You can specify a mesh point by its local or global index, and read or write the values at that point. Every mesh point must have the same number of degrees of freedom. We call the number of degrees of freedom per mesh point the block size.
BlockMultiVector has view semantics. This means that the copy constructor and assignment operator (operator=) all do shallow copies. If you want to create a deep copy, call createCopy(); if you want to copy deeply between two BlockMultiVector objects, call deep_copy().
The replace*, sumInto*, and get* methods all are meant to be thread-safe. They don't throw exceptions on error; instead, they return an error code. They are marked "const" because Kokkos requires this in order to use them in parallel kernels. "Const" is legitimate here, even though some of these methods may modify entries of the vector. This is because they do not change any constants or pointers (e.g., they do not reallocate memory).
BlockMultiVector stores entries in a column corresponding to degrees of freedom for the same mesh point contiguously. This is not strictly necessary; we could generalize so that they are stored in a strided fashion, or even take away the layout assumption and only allow access to entries by copy (e.g., getLocalRowCopy() instead of getLocalRowView()).
Here is how you are supposed to use this class:
Note that BlockMultiVector is not a subclass of Tpetra::MultiVector. I did this on purpose, for simplicity. You should really think of BlockMultiVector as an in-place reinterpretation of a Tpetra::MultiVector.
This is the first Tpetra class to assume view semantics of all Tpetra objects. In particular, it assumes that Map has have view semantics, so that its copy constructor and assignment operator do shallow copies, and its default constructor (that takes no arguments) compiles and succeeds, creating a meaningful empty object.
Design philosophy of the new block objects:
Point 1 means that users fill by mesh points, not degrees of freedom. Thus, they mainly care about the distribution of mesh points, not so much about what GID we assign to which degree of freedom. The latter is just another Map (the "point Map") from their perspective, if they care at all. Preconditioners that respect the block structure also just want the mesh Map. Iterative linear solvers treat the matrix and preconditioner as black-box operators. This means that they only care about the domain and range point Maps.
The latter motivates Point 2. BlockMultiVector views a MultiVector, and iterative solvers (and users) can access the MultiVector and work with it directly, along with its (point) Map. It doesn't make sense for BlockMultiVector to implement MultiVector, because the desired fill interfaces of the two classes are different.
Definition at line 146 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
using Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::map_type = Tpetra::Map<LO, GO, Node> |
The specialization of Tpetra::Map that this class uses.
Definition at line 159 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
using Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::mv_type = Tpetra::MultiVector<Scalar, LO, GO, Node> |
The specialization of Tpetra::MultiVector that this class uses.
Definition at line 161 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
using Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::scalar_type = Scalar |
The type of entries in the object.
Definition at line 164 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
using Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::impl_scalar_type = typename mv_type::impl_scalar_type |
The implementation type of entries in the object.
Letting scalar_type and impl_scalar_type differ addresses Tpetra's work-around to deal with missing device macros and volatile overloads in types like std::complex<T>.
Definition at line 170 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
using Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::local_ordinal_type = LO |
The type of local indices.
Definition at line 172 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
using Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::global_ordinal_type = GO |
The type of global indices.
Definition at line 174 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
using Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::device_type = typename mv_type::device_type |
The Kokkos Device type.
Definition at line 176 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
using Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::node_type = Node |
The Kokkos Node type; a legacy thing that will go away at some point.
Definition at line 179 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
using Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::buffer_device_type = typename dist_object_type::buffer_device_type |
Kokkos::Device specialization used for communication buffers.
Definition at line 182 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
typedef Kokkos::View<impl_scalar_type*, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::little_vec_type |
"Block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector.
A "block view" lets you address all degrees of freedom at a mesh point. You don't have to use this class to access the degrees of freedom. If you do choose to use this class, it implements operator()(LO i), so you can access and modify its entries.
The preferred way to refer to the little_vec_type and const_little_vec_type types, is to get them from the typedefs below. This is because different specializations of BlockVector reserve the right to use different types to implement little_vec_type or const_little_vec_type. This was our porting strategy circa 2014 to move from "classic" Tpetra to the Kokkos refactor version.
Definition at line 203 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
typedef Kokkos::View<const impl_scalar_type*, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::const_little_vec_type |
"Const block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector.
This is just like little_vec_type, except that you can't modify its entries.
Definition at line 214 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inherited |
The Kokkos execution space.
Definition at line 349 of file Tpetra_DistObject_decl.hpp.
Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::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).
meshMap | [in] Map that describes the distribution of mesh points (rather than the distribution of unknowns for those mesh points). |
blockSize | [in] The number of degrees of freedom per mesh point. We assume that this is the same for all mesh points in the above Map. |
numVecs | [in] Number of vectors (columns) in the MultiVector. |
The mesh Map describes the distribution of mesh points. Its corresponding point Map describes the distribution of degrees of freedom corresponding to those mesh points. If you have already computed the point Map corresponding to the above mesh Map, then it is more efficient to call the four-argument constructor below, that takes both the mesh Map and the point Map.
There are two ways to get the point Map corresponding to a given mesh Map and block size. You may either call the class method makePointMap(), or you may call this three-argument constructor, and then call getPointMap().
The point Map enables reinterpretation of a BlockMultiVector as a standard Tpetra::MultiVector. This lets users solve linear systems with Trilinos' solvers and preconditioners, that expect vectors as Tpetra::MultiVector instances.
Definition at line 126 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::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).
See the documentation of the three-argument constructor above.
Definition at line 142 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::BlockMultiVector | ( | const mv_type & | X_mv, |
const map_type & | meshMap, | ||
const LO | blockSize | ||
) |
View an existing MultiVector.
X_mv | [in/out] The MultiVector to view. It MUST have view semantics; otherwise this constructor throws. Its Map must be the same (in the sense of isSameAs) as the point Map corresponding to the given mesh Map and block size. |
meshMap | [in] The mesh Map to use for interpreting the given MultiVector (in place) as a BlockMultiVector. |
blockSize | [in] The number of degrees of freedom per mesh point. We assume that this is the same for all mesh points. |
Definition at line 159 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::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.
This method corresponds to MultiVector's "offset view" constructor.
Definition at line 218 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::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.
This method corresponds to MultiVector's "offset view" constructor.
Definition at line 235 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::BlockMultiVector | ( | ) |
Default constructor.
Creates an empty BlockMultiVector. An empty BlockMultiVector has zero rows, and block size zero.
Definition at line 251 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
|
static |
Create and return the point Map corresponding to the given mesh Map and block size.
This is a class ("static") method so that you can make and reuse a point Map for creating different BlockMultiVector instances, using the more efficient four-argument constructor.
Definition at line 263 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
|
inline |
Get this BlockMultiVector's (previously computed) point Map.
It is always valid to call this method. A BlockMultiVector always has a point Map. We do not compute the point Map lazily.
Definition at line 319 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inline |
Get the number of degrees of freedom per mesh point.
Definition at line 324 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inline |
Get the number of columns (vectors) in the BlockMultiVector.
Definition at line 329 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
BlockMultiVector< Scalar, LO, GO, Node >::mv_type Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::getMultiVectorView | ( | ) | const |
Get a Tpetra::MultiVector that views this BlockMultiVector's data.
This is how you can give a BlockMultiVector to Trilinos' solvers and preconditioners.
Definition at line 105 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
void Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::putScalar | ( | const Scalar & | val | ) |
Fill all entries with the given value val
.
Definition at line 897 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
void Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::scale | ( | const Scalar & | val | ) |
Multiply all entries in place by the given value val
.
Definition at line 904 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
void Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::update | ( | const Scalar & | alpha, |
const BlockMultiVector< Scalar, LO, GO, Node > & | X, | ||
const Scalar & | beta | ||
) |
Update: this = beta*this + alpha*X
.
Update this BlockMultiVector with scaled values of X. If beta is zero, overwrite *this
unconditionally, even if it contains NaN entries. It is legal for the input X to alias this MultiVector.
Definition at line 911 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
void Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::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.
Compute *this := alpha * D * X
, where D is a block diagonal matrix, stored as a 3-D Kokkos::View. This method is the block analog of Tpetra::MultiVector::elementWiseMultiply, and is likewise useful for implementing (block) Jacobi.
alpha | [in] Coefficient by which to scale the result. We treat alpha = 0 as a special case, following the BLAS rules. That is, if alpha = 0, this method does this->putScalar(0) . |
D | [in] Block diagonal, as a 3-D Kokkos::View. The leftmost index indicates which block, the middle index the row within a block, and the rightmost index the column within a block. |
pivots | [in] Pivots (from LU factorization of the blocks) |
X | [in] Input Block(Multi)Vector; may alias *this . |
D is really the inverse of some BlockCrsMatrix's block diagonal. You may compute the inverse of each block however you like. One way is to use GETRF, then GETRI.
Definition at line 1104 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
void Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::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 .
This method computes the block Jacobi update , where Y is *this<>, D the (explicitly stored) inverse block diagonal of a BlockCrsMatrix A, and . The method may use Z as scratch space.
Folks who optimize sparse matrix-vector multiply kernels tend not to write special-purpose kernels like this one. Thus, this kernel consolidates all the other code that block Jacobi needs, while exploiting the existing sparse matrix-vector multiply kernel in BlockCrsMatrix. That consolidation minimizes thread-parallel kernel launch overhead.
alpha | [in] Coefficient of the "block scaled" term. We treat alpha = 0 as a special case, following the BLAS rules. That is, if alpha = 0, this method does Y = beta * Y. |
D | [in] Block diagonal, as a 3-D Kokkos::View. The leftmost index indicates which block, the middle index the row within a block, and the rightmost index the column within a block. |
X | [in] The first of two block (multi)vectors whose difference is "block scaled" |
Z | [in/out] On input: The second of two block (multi)vectors whose difference is "block scaled." This method may use Z as scratch space. |
beta | [in] Coefficient of Y. We treat beta = 0 as a special case, following the BLAS rules. That is, if beta = 0, the initial contents of Y are ignored. |
Definition at line 1136 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
|
inline |
Update data to the given target memory space, only if data in the "other" space have been marked as modified.
If TargetMemorySpace
is the same as this object's memory space, then copy data from host to device. Otherwise, copy data from device to host. In either case, only copy if the source of the copy has been modified.
This is a one-way synchronization only. If the target of the copy has been modified, this operation will discard those modifications. It will also reset both device and host modified flags.
Definition at line 449 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inline |
Update data to the host.
Definition at line 454 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inline |
Update data to the device.
Definition at line 459 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inline |
Whether this object needs synchronization to the given memory space.
Definition at line 465 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inline |
Whether this object needs synchronization to the host.
Definition at line 470 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inline |
Whether this object needs synchronization to the device.
Definition at line 475 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inline |
Mark data as modified on the given memory space.
If TargetDeviceType::memory_space
is the same as this object's memory space, then mark the device's data as modified. Otherwise, mark the host's data as modified.
Definition at line 485 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inline |
Mark data as modified on the host.
Definition at line 490 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inline |
Mark data as modified on the device.
Definition at line 495 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
bool Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::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.
localRowIndex | [in] Local index of the mesh point. |
colIndex | [in] Column (vector) to modify. |
vals | [in] Input values with which to replace whatever existing values are at the mesh point. |
Definition at line 322 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
bool Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::replaceGlobalValues | ( | const GO | globalRowIndex, |
const LO | colIndex, | ||
const Scalar | vals[] | ||
) | const |
Replace all values at the given mesh point, using a global index.
globalRowIndex | [in] Global index of the mesh point. |
colIndex | [in] Column (vector) to modify. |
vals | [in] Input values with which to sum into whatever existing values are at the mesh point. |
Definition at line 337 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
bool Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::sumIntoLocalValues | ( | const LO | localRowIndex, |
const LO | colIndex, | ||
const Scalar | vals[] | ||
) | const |
Sum into all values at the given mesh point, using a local index.
localRowIndex | [in] Local index of the mesh point. |
colIndex | [in] Column (vector) to modify. |
vals | [in] Input values with which to replace whatever existing values are at the mesh point. |
Definition at line 366 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
bool Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::sumIntoGlobalValues | ( | const GO | globalRowIndex, |
const LO | colIndex, | ||
const Scalar | vals[] | ||
) | const |
Sum into all values at the given mesh point, using a global index.
globalRowIndex | [in] Global index of the mesh point. |
colIndex | [in] Column (vector) to modify. |
vals | [in] Input values with which to replace whatever existing values are at the mesh point. |
Definition at line 381 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
bool Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::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.
localRowIndex | [in] Local index of the mesh point. |
colIndex | [in] Column (vector) to view. |
vals | [out] View of the entries at the given mesh point. |
Definition at line 397 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
bool Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::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.
globalRowIndex | [in] Global index of the mesh point. |
colIndex | [in] Column (vector) to view. |
vals | [out] View of the entries at the given mesh point. |
Definition at line 411 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
BlockMultiVector< Scalar, LO, GO, Node >::little_vec_type::HostMirror Tpetra::Experimental::BlockMultiVector< Scalar, LO, GO, Node >::getLocalBlock | ( | const LO | localRowIndex, |
const LO | colIndex | ||
) | const |
Get a host view of the degrees of freedom at the given mesh point.
Prefer using auto
to let the compiler compute the return type. This gives us the freedom to change this type in the future. If you insist not to use auto
, then please use the little_vec_type
typedef to deduce the correct return type; don't try to hard-code the return type yourself.
Definition at line 426 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
|
protectedvirtual |
Compare the source and target (this) objects for compatibility.
Implements Tpetra::DistObject< Scalar, LO, GO, Node >.
Definition at line 483 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
|
inlineprotected |
Raw pointer to the MultiVector's data.
Definition at line 655 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inlineprotected |
Stride between consecutive local entries in the same column.
Definition at line 660 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
inlineprotected |
Stride between consecutive local entries in the same row.
Definition at line 665 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
protected |
True if and only if meshLocalIndex
is a valid local index in the mesh Map.
Definition at line 889 of file Tpetra_Experimental_BlockMultiVector_def.hpp.
|
inherited |
Import data into this object using an Import object ("forward mode").
The input DistObject is always the source of the data redistribution operation, and the *this
object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use this method with your precomputed Import object if you want to do an Import, else use doExport() with a precomputed Export object.
"Restricted Mode" does two things:
*this
, in a "locallyFitted" sense. This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.
source | [in] The "source" object for redistribution. |
importer | [in] Precomputed data redistribution plan. Its source Map must be the same as the input DistObject's Map, and its target Map must be the same as this->getMap() . |
CM | [in] How to combine incoming data with the same global index. |
|
inherited |
Import data into this object using an Export object ("reverse mode").
The input DistObject is always the source of the data redistribution operation, and the *this
object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use the version of doImport() that takes a precomputed Import object in that case.
"Restricted Mode" does two things:
*this
, in a "locallyFitted" sense. This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.
source | [in] The "source" object for redistribution. |
exporter | [in] Precomputed data redistribution plan. Its target Map must be the same as the input DistObject's Map, and its source Map must be the same as this->getMap() . (Note the difference from forward mode.) |
CM | [in] How to combine incoming data with the same global index. |
|
inherited |
Export data into this object using an Export object ("forward mode").
The input DistObject is always the source of the data redistribution operation, and the *this
object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use this method with your precomputed Export object if you want to do an Export, else use doImport() with a precomputed Import object.
"Restricted Mode" does two things:
*this
, in a "locallyFitted" sense. This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.
source | [in] The "source" object for redistribution. |
exporter | [in] Precomputed data redistribution plan. Its source Map must be the same as the input DistObject's Map, and its target Map must be the same as this->getMap() . |
CM | [in] How to combine incoming data with the same global index. |
|
inherited |
Export data into this object using an Import object ("reverse mode").
The input DistObject is always the source of the data redistribution operation, and the *this
object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use the version of doExport() that takes a precomputed Export object in that case.
"Restricted Mode" does two things:
*this
, in a "locallyFitted" sense. This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.
source | [in] The "source" object for redistribution. |
importer | [in] Precomputed data redistribution plan. Its target Map must be the same as the input DistObject's Map, and its source Map must be the same as this->getMap() . (Note the difference from forward mode.) |
CM | [in] How to combine incoming data with the same global index. |
|
inherited |
Whether this is a globally distributed object.
For a definition of "globally distributed" (and its opposite, "locally replicated"), see the documentation of Map's isDistributed() method.
|
inlinevirtualinherited |
The Map describing the parallel distribution of this object.
Note that some Tpetra objects might be distributed using multiple Map objects. For example, CrsMatrix has both a row Map and a column Map. It is up to the subclass to decide which Map to use when invoking the DistObject constructor.
Definition at line 541 of file Tpetra_DistObject_decl.hpp.
|
inherited |
Print this object to the given output stream.
We generally assume that all MPI processes can print to the given stream.
|
virtualinherited |
One-line descriptiion of this object.
We declare this method virtual so that subclasses of DistObject may override it.
Reimplemented in Tpetra::MultiVector< Scalar, LO, GO, Node >.
|
virtualinherited |
Print a descriptiion of this object to the given output stream.
We declare this method virtual so that subclasses of Distobject may override it.
Reimplemented in Tpetra::MultiVector< Scalar, LO, GO, Node >.
|
virtualinherited |
Remove processes which contain no entries in this object's Map.
On input, this object is distributed over the Map returned by getMap() (the "original Map," with its communicator, the "original communicator"). The input newMap
of this method must be the same as the result of calling getMap()->removeEmptyProcesses()
. On processes in the original communicator which contain zero entries ("excluded processes," as opposed to "included processes"), the input newMap
must be Teuchos::null
(which is what getMap()->removeEmptyProcesses()
returns anyway).
On included processes, reassign this object's Map (that would be returned by getMap()) to the input newMap
, and do any work that needs to be done to restore correct semantics. On excluded processes, free any data that needs freeing, and do any other work that needs to be done to restore correct semantics.
This method has collective semantics over the original communicator. On exit, the only method of this object which is safe to call on excluded processes is the destructor. This implies that subclasses' destructors must not contain communication operations.
Reimplemented in Tpetra::MultiVector< Scalar, LO, GO, Node >.
|
protectedvirtualinherited |
Whether the implementation's instance promises always to have a constant number of packets per LID (local index), and if so, how many packets per LID there are.
If this method returns zero, the instance says that it might possibly have a different number of packets for each LID (local index) to send or receive. If it returns nonzero, the instance promises that the number of packets is the same for all LIDs, and that the return value is this number of packets per LID.
The default implementation of this method returns zero. This does not affect the behavior of doTransfer() in any way. If a nondefault implementation returns nonzero, doTransfer() will use this information to avoid unnecessary allocation and / or resizing of arrays.
Reimplemented in Tpetra::MultiVector< Scalar, LO, GO, Node >.
|
protectedvirtualinherited |
Redistribute data across (MPI) processes.
src | [in] The source object, to redistribute into the target object, which is *this object. |
transfer | [in] The Export or Import object representing the communication pattern. (Details::Transfer is the common base class of these two objects.) |
modeString | [in] Human-readable string, for verbose debugging output and error output, explaining what function called this method. Example: "doImport (forward)", "doExport (reverse)". |
revOp | [in] Whether to do a forward or reverse mode redistribution. |
CM | [in] The combine mode that describes how to combine values that map to the same global ID on the same process. |
|
protectedvirtualinherited |
Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary.
numExportLIDs | [in] Number of entries in the exportLIDs input array argument of doTransfer(). |
numImportLIDs | [in] Number of entries in the remoteLIDs input array argument of doTransfer(). |
|
protectedvirtualinherited |
Implementation detail of doTransfer.
LID DualViews come from the Transfer object given to doTransfer. They are always sync'd on both host and device. Users must never attempt to modify or sync them.
|
protectedvirtualinherited |
Perform copies and permutations that are local to the calling (MPI) process.
Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this
, packs the source object's data.
permuteToLIDs.need_sync_host()
, permuteToLIDs.need_sync_device()
, permuteFromLIDs.need_sync_host()
, and permuteFromLIDs.need_sync_device()
are all false.source | [in] On entry, the source object of the Export or Import operation. |
numSameIDs | [in] The number of elements that are the same on the source and target objects. These elements live on the same process in both the source and target objects. |
permuteToLIDs | [in] List of the elements that are permuted. They are listed by their local index (LID) in the destination object. |
permuteFromLIDs | [in] List of the elements that are permuted. They are listed by their local index (LID) in the source object. |
|
protectedvirtualinherited |
Pack data and metadata for communication (sends).
Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this
, packs the source object's data.
exportLIDs.need_sync_host ()
and exportLIDs.need_sync_device()
are both false.source | [in] Source object for the redistribution. |
exportLIDs | [in] List of the entries (as local IDs in the source object) that Tpetra will send to other processes. |
exports | [out] On exit, the packed data to send. Implementations must reallocate this as needed (prefer reusing the existing allocation if possible), and may modify and/or sync this wherever they like. |
numPacketsPerLID | [out] On exit, the implementation of this method must do one of two things: either set numPacketsPerLID[i] to the number of packets to be packed for exportLIDs[i] and set constantNumPackets to zero, or set constantNumPackets to a nonzero value. If the latter, the implementation must not modify the entries of numPacketsPerLID . If the former, the implementation may sync numPacketsPerLID this wherever it likes, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference. |
constantNumPackets | [out] On exit, 0 if the number of packets per LID could differ, else (if nonzero) the number of packets per LID (which must be constant). |
distor | [in] The Distributor object we are using. Most implementations will not use this. |
|
protectedvirtualinherited |
Perform any unpacking and combining after communication.
Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this
, unpacks the received data into itself, possibly modifying its entries.
importLIDs.need_sync_host ()
and importLIDs.need_sync_device()
are both false.importLIDs | [in] List of the entries (as LIDs in the destination object) we received from other processes. |
imports | [in/out] On input: Buffer of received data to unpack. DistObject promises nothing about where this is sync'd. Implementations may sync this wherever they like, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference. |
numPacketsPerLID | [in/out] On input: If constantNumPackets is zero, then numPacketsPerLID[i] contains the number of packets imported for importLIDs[i]. DistObject promises nothing about where this is sync'd. Implementations may sync this wherever they like, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference. |
constantNumPackets | [in] If nonzero, then the number of packets per LID is the same for all entries ("constant") and constantNumPackets is that number. If zero, then numPacketsPerLID[i] is the number of packets to unpack for LID importLIDs[i] . |
distor | [in] The Distributor object we are using. Most implementations will not use this. |
combineMode | [in] The CombineMode to use when combining the imported entries with existing entries. |
|
protectedinherited |
Reallocate imports_ if needed.
This unfortunately must be declared protected, for the same reason that imports_ is declared protected.
newSize | [in] New size of imports_. |
verbose | [in] Whether to print verbose debugging output to stderr on every (MPI) process in the communicator. |
prefix | [in] If verbose is true , then this is a nonnull prefix to print at the beginning of each line of verbose debugging output. Otherwise, not used. |
We don't need a "reallocExportsIfNeeded" method, because exports_
always gets passed into packAndPrepare() by nonconst reference. Thus, that method can resize the DualView without needing to call other DistObject methods.
|
protected |
The Tpetra::MultiVector used to represent the data.
Definition at line 686 of file Tpetra_Experimental_BlockMultiVector_decl.hpp.
|
protectedinherited |
The Map over which this object is distributed.
Definition at line 1140 of file Tpetra_DistObject_decl.hpp.
|
protectedinherited |
Buffer into which packed data are imported (received from other processes).
Unfortunately, I had to declare these protected, because CrsMatrix uses them at one point. Please, nobody else use them.
Definition at line 1149 of file Tpetra_DistObject_decl.hpp.
|
protectedinherited |
Number of packets to receive for each receive operation.
This array is used in Distributor::doPosts() (and doReversePosts()) when starting the ireceive operation.
This may be ignored in doTransfer() if constantNumPackets is nonzero, indicating a constant number of packets per LID. (For example, MultiVector sets the constantNumPackets output argument of packAndPrepare() to the number of columns in the multivector.)
Unfortunately, I had to declare this protected, because CrsMatrix uses it at one point. Please, nobody else use it.
Definition at line 1187 of file Tpetra_DistObject_decl.hpp.
|
protectedinherited |
Buffer from which packed data are exported (sent to other processes).
Unfortunately, I had to declare this protected, because CrsMatrix uses it at one point. Please, nobody else use it.
Definition at line 1194 of file Tpetra_DistObject_decl.hpp.
|
protectedinherited |
Number of packets to send for each send operation.
This array is used in Distributor::doPosts() (and doReversePosts()) for preparing for the send operation.
This may be ignored in doTransfer() if constantNumPackets is nonzero, indicating a constant number of packets per LID. (For example, MultiVector sets the constantNumPackets output argument of packAndPrepare() to the number of columns in the multivector.)
Unfortunately, I had to declare this protected, because CrsMatrix uses them at one point. Please, nobody else use it.
Definition at line 1209 of file Tpetra_DistObject_decl.hpp.