Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Protected Member Functions | Protected Attributes | List of all members
Tpetra::BlockVector< Scalar, LO, GO, Node > Class Template Reference

Vector for multiple degrees of freedom per mesh point. More...

#include <Tpetra_BlockVector_decl.hpp>

Inheritance diagram for Tpetra::BlockVector< Scalar, LO, GO, Node >:
Inheritance graph
[legend]

Public Types

Typedefs to facilitate template metaprogramming.
typedef base_type::scalar_type scalar_type
 The type of entries in the vector. More...
 
typedef base_type::impl_scalar_type impl_scalar_type
 The implementation type of entries in the vector. More...
 
typedef
base_type::local_ordinal_type 
local_ordinal_type
 The type of local indices. More...
 
typedef
base_type::global_ordinal_type 
global_ordinal_type
 The type of global indices. More...
 
typedef base_type::node_type node_type
 The Kokkos Node type. More...
 
typedef Node::device_type device_type
 The Kokkos Device type. More...
 
typedef Tpetra::Map< LO, GO, Node > map_type
 The specialization of Tpetra::Map that this class uses. More...
 
typedef Tpetra::MultiVector
< Scalar, LO, GO, Node > 
mv_type
 The specialization of Tpetra::MultiVector that this class uses. More...
 
typedef Tpetra::Vector< Scalar,
LO, GO, Node > 
vec_type
 The specialization of Tpetra::Vector that this class uses. 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. 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. More...
 
Typedefs to facilitate template metaprogramming.
using buffer_device_type = typename dist_object_type::buffer_device_type
 Kokkos::Device specialization used for communication buffers. More...
 
Typedefs
using execution_space = typename device_type::execution_space
 The Kokkos execution space. More...
 

Public Member Functions

Constructors
 BlockVector ()
 Default constructor. More...
 
 BlockVector (const BlockVector< Scalar, LO, GO, Node > &)=default
 Copy constructor (shallow copy). More...
 
 BlockVector (BlockVector< Scalar, LO, GO, Node > &&)=default
 Move constructor (shallow move). More...
 
BlockVector< Scalar, LO, GO,
Node > & 
operator= (const BlockVector< Scalar, LO, GO, Node > &)=default
 Copy assigment (shallow copy). More...
 
BlockVector< Scalar, LO, GO,
Node > & 
operator= (BlockVector< Scalar, LO, GO, Node > &&)=default
 Move assigment (shallow move). More...
 
 BlockVector (const BlockVector< Scalar, LO, GO, Node > &in, const Teuchos::DataAccess copyOrView)
 "Copy constructor" with option to deep copy. More...
 
 BlockVector (const map_type &meshMap, const LO blockSize)
 Constructor that takes a mesh Map and a block size. More...
 
 BlockVector (const map_type &meshMap, const map_type &pointMap, const LO blockSize)
 Constructor that takes a mesh Map, a point Map, and a block size. More...
 
 BlockVector (const mv_type &X_mv, const map_type &meshMap, const LO blockSize)
 View an existing Vector or MultiVector. More...
 
 BlockVector (const vec_type &X_vec, const map_type &meshMap, const LO blockSize)
 View an existing Vector. More...
 
 BlockVector (const BlockVector< Scalar, LO, GO, Node > &X, const map_type &newMeshMap, const map_type &newPointMap, const size_t offset=0)
 View an existing BlockVector using a different mesh Map, supplying the corresponding point Map. More...
 
 BlockVector (const BlockVector< Scalar, LO, GO, Node > &X, const map_type &newMeshMap, const size_t offset=0)
 View an existing BlockVector using a different mesh Map; compute the new point Map. More...
 
Access to Maps, the block size, and a Vector view.
vec_type getVectorView ()
 Get a Tpetra::Vector that views this BlockVector's data. More...
 
Fine-grained data access
bool replaceLocalValues (const LO localRowIndex, const Scalar vals[]) const
 Replace all values at the given mesh point, using a local index. More...
 
bool replaceGlobalValues (const GO globalRowIndex, const Scalar vals[]) const
 Replace all values at the given mesh point, using a global index. More...
 
bool sumIntoLocalValues (const LO localRowIndex, const Scalar vals[]) const
 Sum into all values at the given mesh point, using a local index. More...
 
bool sumIntoGlobalValues (const GO globalRowIndex, const Scalar vals[]) const
 Sum into all values at the given mesh point, using a global index. More...
 
bool getLocalRowView (const LO localRowIndex, 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, Scalar *&vals) const
 Get a writeable view of the entries at the given mesh point, using a global index. More...
 
little_vec_type getLocalBlock (const LO localRowIndex) const
 Get a view of the degrees of freedom at the given mesh point, using a local index. 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 $Y = \beta * Y + \alpha D (X - Z)$. 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_typegetRawPtr () 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...
 

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...
 

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_typemap_
 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...
 
std::unique_ptr< std::string > createPrefix (const char className[], const char methodName[]) const
 
bool reallocImportsIfNeeded (const size_t newSize, const bool verbose, const std::string *prefix)
 Reallocate imports_ if needed. More...
 

Detailed Description

template<class Scalar, class LO, class GO, class Node>
class Tpetra::BlockVector< Scalar, LO, GO, Node >

Vector for multiple degrees of freedom per mesh point.

Author
Mark Hoemmen
Template Parameters
ScalarThe type of each entry of the block vector. (You can use real-valued or complex-valued types here, unlike in Epetra, where the scalar type is always double.)
LOThe type of local indices. See the documentation of the first template parameter of Map for requirements.
GOThe type of global indices. See the documentation of the second template parameter of Map for requirements.
NodeThe Kokkos Node type. See the documentation of the third template parameter of Map for requirements.

BlockVector 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.

BlockVector is a special case of BlockMultiVector, for "multivectors" that are not "multi." That is, a BlockVector has a single vector (column). Please refer to the documentation of BlockMultiVector for details.

Definition at line 80 of file Tpetra_BlockVector_decl.hpp.

Member Typedef Documentation

template<class Scalar, class LO, class GO, class Node>
typedef base_type::scalar_type Tpetra::BlockVector< Scalar, LO, GO, Node >::scalar_type

The type of entries in the vector.

Definition at line 90 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
typedef base_type::impl_scalar_type Tpetra::BlockVector< Scalar, LO, GO, Node >::impl_scalar_type

The implementation type of entries in the vector.

Definition at line 92 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
typedef base_type::local_ordinal_type Tpetra::BlockVector< Scalar, LO, GO, Node >::local_ordinal_type

The type of local indices.

Definition at line 94 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
typedef base_type::global_ordinal_type Tpetra::BlockVector< Scalar, LO, GO, Node >::global_ordinal_type

The type of global indices.

Definition at line 96 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
typedef base_type::node_type Tpetra::BlockVector< Scalar, LO, GO, Node >::node_type

The Kokkos Node type.

Definition at line 98 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
typedef Node::device_type Tpetra::BlockVector< Scalar, LO, GO, Node >::device_type

The Kokkos Device type.

Definition at line 100 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
typedef Tpetra::Map<LO, GO, Node> Tpetra::BlockVector< Scalar, LO, GO, Node >::map_type

The specialization of Tpetra::Map that this class uses.

Definition at line 103 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
typedef Tpetra::MultiVector<Scalar, LO, GO, Node> Tpetra::BlockVector< Scalar, LO, GO, Node >::mv_type

The specialization of Tpetra::MultiVector that this class uses.

Definition at line 105 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
typedef Tpetra::Vector<Scalar, LO, GO, Node> Tpetra::BlockVector< Scalar, LO, GO, Node >::vec_type

The specialization of Tpetra::Vector that this class uses.

Definition at line 107 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
typedef Kokkos::View<impl_scalar_type*, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Tpetra::BlockVector< Scalar, LO, GO, Node >::little_vec_type

"Block view" of all degrees of freedom at a mesh point.

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 gives us a porting strategy to move from "classic" Tpetra to the Kokkos refactor version.

Definition at line 127 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
typedef Kokkos::View<const impl_scalar_type*, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Tpetra::BlockVector< Scalar, LO, GO, Node >::const_little_vec_type

"Const block view" of all degrees of freedom at a mesh point.

This is just like little_vec_type, except that you can't modify its entries.

Definition at line 137 of file Tpetra_BlockVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
using Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::buffer_device_type = typename dist_object_type::buffer_device_type
inherited

Kokkos::Device specialization used for communication buffers.

Definition at line 179 of file Tpetra_BlockMultiVector_decl.hpp.

using Tpetra::DistObject< Scalar , LO , GO , Node >::execution_space = typename device_type::execution_space
inherited

The Kokkos execution space.

Definition at line 348 of file Tpetra_DistObject_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar , class LO , class GO , class Node >
Tpetra::BlockVector< Scalar, LO, GO, Node >::BlockVector ( )

Default constructor.

Creates an empty BlockVector. An empty BlockVector has zero rows, and block size zero.

Definition at line 49 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
Tpetra::BlockVector< Scalar, LO, GO, Node >::BlockVector ( const BlockVector< Scalar, LO, GO, Node > &  )
default

Copy constructor (shallow copy).

template<class Scalar, class LO, class GO, class Node>
Tpetra::BlockVector< Scalar, LO, GO, Node >::BlockVector ( BlockVector< Scalar, LO, GO, Node > &&  )
default

Move constructor (shallow move).

template<class Scalar, class LO, class GO, class Node>
Tpetra::BlockVector< Scalar, LO, GO, Node >::BlockVector ( const BlockVector< Scalar, LO, GO, Node > &  in,
const Teuchos::DataAccess  copyOrView 
)

"Copy constructor" with option to deep copy.

Definition at line 55 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
Tpetra::BlockVector< Scalar, LO, GO, Node >::BlockVector ( const map_type meshMap,
const LO  blockSize 
)

Constructor that takes a mesh Map and a block size.

Parameters
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.

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 three-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() (inherited from the parent class BlockMultiVector), or you may call this two-argument constructor, and then call getPointMap().

The point Map enables reinterpretation of a BlockVector as a standard Tpetra::Vector, or as a Tpetra::MultiVector with one column. This lets users solve linear systems with Trilinos' solvers and preconditioners, that expect vectors as Tpetra::MultiVector or Tpetra::Vector instances.

Definition at line 62 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
Tpetra::BlockVector< Scalar, LO, GO, Node >::BlockVector ( const map_type meshMap,
const map_type pointMap,
const LO  blockSize 
)

Constructor that takes a mesh Map, a point Map, and a block size.

See the documentation of the two-argument constructor above.

Definition at line 68 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
Tpetra::BlockVector< Scalar, LO, GO, Node >::BlockVector ( const mv_type X_mv,
const map_type meshMap,
const LO  blockSize 
)

View an existing Vector or MultiVector.

Parameters
X_mv[in/out] The Vector or 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. If this is a MultiVector, it must have only one column.
meshMap[in] The mesh Map to use for interpreting the given MultiVector or Vector (in place) as a BlockVector.
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 76 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
Tpetra::BlockVector< Scalar, LO, GO, Node >::BlockVector ( const vec_type X_vec,
const map_type meshMap,
const LO  blockSize 
)

View an existing Vector.

Parameters
X_vec[in/out] The Vector 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 Vector (in place) as a BlockVector.
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 89 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
Tpetra::BlockVector< Scalar, LO, GO, Node >::BlockVector ( const BlockVector< Scalar, LO, GO, Node > &  X,
const map_type newMeshMap,
const map_type newPointMap,
const size_t  offset = 0 
)

View an existing BlockVector using a different mesh Map, supplying the corresponding point Map.

This method corresponds to MultiVector's "offset view" constructor.

Definition at line 97 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
Tpetra::BlockVector< Scalar, LO, GO, Node >::BlockVector ( const BlockVector< Scalar, LO, GO, Node > &  X,
const map_type newMeshMap,
const size_t  offset = 0 
)

View an existing BlockVector using a different mesh Map; compute the new point Map.

This method corresponds to MultiVector's "offset view" constructor.

Definition at line 106 of file Tpetra_BlockVector_def.hpp.

Member Function Documentation

template<class Scalar, class LO, class GO, class Node>
BlockVector<Scalar, LO, GO, Node>& Tpetra::BlockVector< Scalar, LO, GO, Node >::operator= ( const BlockVector< Scalar, LO, GO, Node > &  )
default

Copy assigment (shallow copy).

template<class Scalar, class LO, class GO, class Node>
BlockVector<Scalar, LO, GO, Node>& Tpetra::BlockVector< Scalar, LO, GO, Node >::operator= ( BlockVector< Scalar, LO, GO, Node > &&  )
default

Move assigment (shallow move).

template<class Scalar , class LO , class GO , class Node >
BlockVector< Scalar, LO, GO, Node >::vec_type Tpetra::BlockVector< Scalar, LO, GO, Node >::getVectorView ( )

Get a Tpetra::Vector that views this BlockVector's data.

This is how you can give a BlockVector to Trilinos' solvers and preconditioners.

Definition at line 114 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO, class GO , class Node >
bool Tpetra::BlockVector< Scalar, LO, GO, Node >::replaceLocalValues ( const LO  localRowIndex,
const Scalar  vals[] 
) const

Replace all values at the given mesh point, using a local index.

Parameters
localRowIndex[in] Local index of the mesh point.
vals[in] Input values with which to replace whatever existing values are at the mesh point.
Returns
true if successful, else false. This method will not succeed if the given local index of the mesh point is invalid on the calling process.
Note
This method, the other "replace" and "sumInto" methods, and the view methods, are marked const. This is because they do not change pointers. They do, of course, change the values in the BlockVector, but that does not require marking the methods as nonconst.

Definition at line 122 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO , class GO, class Node >
bool Tpetra::BlockVector< Scalar, LO, GO, Node >::replaceGlobalValues ( const GO  globalRowIndex,
const Scalar  vals[] 
) const

Replace all values at the given mesh point, using a global index.

Parameters
globalRowIndex[in] Global index of the mesh point.
vals[in] Input values with which to sum into whatever existing values are at the mesh point.
Returns
true if successful, else false. This method will not succeed if the given global index of the mesh point is invalid on the calling process.

Definition at line 129 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO, class GO , class Node >
bool Tpetra::BlockVector< Scalar, LO, GO, Node >::sumIntoLocalValues ( const LO  localRowIndex,
const Scalar  vals[] 
) const

Sum into all values at the given mesh point, using a local index.

Parameters
localRowIndex[in] Local index of the mesh point.
vals[in] Input values with which to replace whatever existing values are at the mesh point.
Returns
true if successful, else false. This method will not succeed if the given local index of the mesh point is invalid on the calling process.

Definition at line 136 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO , class GO, class Node >
bool Tpetra::BlockVector< Scalar, LO, GO, Node >::sumIntoGlobalValues ( const GO  globalRowIndex,
const Scalar  vals[] 
) const

Sum into all values at the given mesh point, using a global index.

Parameters
globalRowIndex[in] Global index of the mesh point.
vals[in] Input values with which to replace whatever existing values are at the mesh point.
Returns
true if successful, else false. This method will not succeed if the given global index of the mesh point is invalid on the calling process.

Definition at line 143 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO, class GO , class Node >
bool Tpetra::BlockVector< Scalar, LO, GO, Node >::getLocalRowView ( const LO  localRowIndex,
Scalar *&  vals 
) const

Get a writeable view of the entries at the given mesh point, using a local index.

Parameters
localRowIndex[in] Local index of the mesh point.
vals[in] Input values with which to replace whatever existing values are at the mesh point.
Returns
true if successful, else false. This method will not succeed if the given local index of the mesh point is invalid on the calling process.

Definition at line 150 of file Tpetra_BlockVector_def.hpp.

template<class Scalar, class LO , class GO, class Node >
bool Tpetra::BlockVector< Scalar, LO, GO, Node >::getGlobalRowView ( const GO  globalRowIndex,
Scalar *&  vals 
) const

Get a writeable view of the entries at the given mesh point, using a global index.

Parameters
globalRowIndex[in] Global index of the mesh point.
vals[in] Input values with which to replace whatever existing values are at the mesh point.
Returns
true if successful, else false. This method will not succeed if the given global index of the mesh point is invalid on the calling process.

Definition at line 157 of file Tpetra_BlockVector_def.hpp.

template<class Scalar , class LO, class GO , class Node >
BlockVector< Scalar, LO, GO, Node >::little_vec_type Tpetra::BlockVector< Scalar, LO, GO, Node >::getLocalBlock ( const LO  localRowIndex) const

Get a view of the degrees of freedom at the given mesh point, using a local index.

Warning
This method's interface may change or disappear at any time. Please do not rely on it in your code yet.

The preferred way to refer to little_vec_type is to get it from BlockVector's typedef. This is because different specializations of BlockVector reserve the right to use different types to implement little_vec_type. This gives us a porting strategy to move from "classic" Tpetra to the Kokkos refactor version.

Definition at line 164 of file Tpetra_BlockVector_def.hpp.

template<class Scalar , class LO, class GO , class Node >
BlockMultiVector< Scalar, LO, GO, Node >::map_type Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::makePointMap ( const map_type meshMap,
const LO  blockSize 
)
staticinherited

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 256 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
map_type Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::getPointMap ( ) const
inlineinherited

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 334 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
LO Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::getBlockSize ( ) const
inlineinherited

Get the number of degrees of freedom per mesh point.

Definition at line 339 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
LO Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::getNumVectors ( ) const
inlineinherited

Get the number of columns (vectors) in the BlockMultiVector.

Definition at line 344 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar , class LO , class GO , class Node >
BlockMultiVector< Scalar, LO, GO, Node >::mv_type Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::getMultiVectorView ( ) const
inherited

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 102 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO , class GO , class Node >
void Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::putScalar ( const Scalar &  val)
inherited

Fill all entries with the given value val.

Definition at line 545 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO , class GO , class Node >
void Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::scale ( const Scalar &  val)
inherited

Multiply all entries in place by the given value val.

Definition at line 552 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::update ( const Scalar &  alpha,
const BlockMultiVector< Scalar, LO, GO, Node > &  X,
const Scalar &  beta 
)
inherited

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 559 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::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 
)
inherited

*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.

Parameters
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 752 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::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 
)
inherited

Block Jacobi update $Y = \beta * Y + \alpha D (X - Z)$.

This method computes the block Jacobi update $Y = \beta * Y + \alpha D (X - Z)$, where Y is *this<>, D the (explicitly stored) inverse block diagonal of a BlockCrsMatrix A, and $Z = A*Y$. 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.

Parameters
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 784 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO, class GO, class Node>
template<class TargetMemorySpace >
void Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::sync ( )
inlineinherited

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.

Note
This method doesn't know on its own whether you modified the data in either memory space. You must manually mark the MultiVector as modified in the space in which you modified it, by calling the modify() method with the appropriate template parameter.

Definition at line 464 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::sync_host ( )
inlineinherited

Update data to the host.

Definition at line 469 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::sync_device ( )
inlineinherited

Update data to the device.

Definition at line 474 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
template<class TargetMemorySpace >
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::need_sync ( ) const
inlineinherited

Whether this object needs synchronization to the given memory space.

Definition at line 480 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::need_sync_host ( ) const
inlineinherited

Whether this object needs synchronization to the host.

Definition at line 485 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::need_sync_device ( ) const
inlineinherited

Whether this object needs synchronization to the device.

Definition at line 490 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
template<class TargetMemorySpace >
void Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::modify ( )
inlineinherited

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 500 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::modify_host ( )
inlineinherited

Mark data as modified on the host.

Definition at line 505 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::modify_device ( )
inlineinherited

Mark data as modified on the device.

Definition at line 510 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO , class Node >
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::replaceLocalValues ( const LO  localRowIndex,
const LO  colIndex,
const Scalar  vals[] 
) const
inherited

Replace all values at the given mesh point, using local row and column indices.

Parameters
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.
Returns
true if successful, else false. This method will not succeed if the given local index of the mesh point is invalid on the calling process.
Note
This method, the other "replace" and "sumInto" methods, and the view methods, are marked const. This is because they do not change pointers. They do, of course, change the values in the BlockMultiVector, but that does not require marking the methods as nonconst.

Definition at line 315 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO, class GO, class Node >
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::replaceGlobalValues ( const GO  globalRowIndex,
const LO  colIndex,
const Scalar  vals[] 
) const
inherited

Replace all values at the given mesh point, using a global index.

Parameters
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.
Returns
true if successful, else false. This method will not succeed if the given global index of the mesh point is invalid on the calling process.

Definition at line 330 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO, class GO , class Node >
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::sumIntoLocalValues ( const LO  localRowIndex,
const LO  colIndex,
const Scalar  vals[] 
) const
inherited

Sum into all values at the given mesh point, using a local index.

Parameters
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.
Returns
true if successful, else false. This method will not succeed if the given local index of the mesh point is invalid on the calling process.

Definition at line 359 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO, class GO, class Node >
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::sumIntoGlobalValues ( const GO  globalRowIndex,
const LO  colIndex,
const Scalar  vals[] 
) const
inherited

Sum into all values at the given mesh point, using a global index.

Parameters
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.
Returns
true if successful, else false. This method will not succeed if the given global index of the mesh point is invalid on the calling process.

Definition at line 374 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO, class GO , class Node >
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::getLocalRowView ( const LO  localRowIndex,
const LO  colIndex,
Scalar *&  vals 
) const
inherited

Get a writeable view of the entries at the given mesh point, using a local index.

Parameters
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.
Returns
true if successful, else false. This method will not succeed if the given local index of the mesh point is invalid on the calling process.

Definition at line 390 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar, class LO, class GO, class Node >
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::getGlobalRowView ( const GO  globalRowIndex,
const LO  colIndex,
Scalar *&  vals 
) const
inherited

Get a writeable view of the entries at the given mesh point, using a global index.

Parameters
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.
Returns
true if successful, else false. This method will not succeed if the given global index of the mesh point is invalid on the calling process.

Definition at line 404 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar , class LO, class GO , class Node >
BlockMultiVector< Scalar, LO, GO, Node >::little_vec_type::HostMirror Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::getLocalBlock ( const LO  localRowIndex,
const LO  colIndex 
) const
inherited

Get a host view of the degrees of freedom at the given mesh point.

Warning
This method's interface may change or disappear at any time. Please do not rely on it in your code yet.

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 419 of file Tpetra_BlockMultiVector_def.hpp.

template<class Scalar , class LO , class GO , class Node >
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::checkSizes ( const Tpetra::SrcDistObject source)
protectedvirtualinherited

Compare the source and target (this) objects for compatibility.

Returns
True if they are compatible, else false.

Implements Tpetra::DistObject< Scalar, LO, GO, Node >.

Definition at line 476 of file Tpetra_BlockMultiVector_def.hpp.

virtual void Tpetra::DistObject< Scalar , LO , GO , Node >::copyAndPermute ( const SrcDistObject< Scalar, LO, GO, Node > &  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 
)
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.

Precondition
permuteToLIDs and permuteFromLIDs are sync'd to both host and device. That is, permuteToLIDs.need_sync_host(), permuteToLIDs.need_sync_device(), permuteFromLIDs.need_sync_host(), and permuteFromLIDs.need_sync_device() are all false.
Parameters
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.
virtual void Tpetra::DistObject< Scalar , LO , GO , Node >::packAndPrepare ( const SrcDistObject< Scalar, LO, GO, Node > &  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 
)
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.

Precondition
exportLIDs is sync'd to both host and device. That is, exportLIDs.need_sync_host () and exportLIDs.need_sync_device() are both false.
Parameters
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.
virtual void Tpetra::DistObject< Scalar , LO , GO , Node >::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 
)
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.

Precondition
importLIDs is sync'd to both host and device. That is, importLIDs.need_sync_host () and importLIDs.need_sync_device() are both false.
Parameters
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.
template<class Scalar, class LO, class GO, class Node>
impl_scalar_type* Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::getRawPtr ( ) const
inlineprotectedinherited

Raw pointer to the MultiVector's data.

Definition at line 658 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
size_t Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::getStrideX ( ) const
inlineprotectedinherited

Stride between consecutive local entries in the same column.

Definition at line 663 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar, class LO, class GO, class Node>
size_t Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::getStrideY ( ) const
inlineprotectedinherited

Stride between consecutive local entries in the same row.

Definition at line 668 of file Tpetra_BlockMultiVector_decl.hpp.

template<class Scalar , class LO, class GO , class Node >
bool Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::isValidLocalMeshIndex ( const LO  meshLocalIndex) const
protectedinherited

True if and only if meshLocalIndex is a valid local index in the mesh Map.

Definition at line 537 of file Tpetra_BlockMultiVector_def.hpp.

void Tpetra::DistObject< Scalar , LO , GO , Node >::doImport ( const SrcDistObject< Scalar, LO, GO, Node > &  source,
const Import< LO , GO , Node > &  importer,
const CombineMode  CM,
const bool  restrictedMode = false 
)
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:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *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.

Parameters
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.
void Tpetra::DistObject< Scalar , LO , GO , Node >::doImport ( const SrcDistObject< Scalar, LO, GO, Node > &  source,
const Export< LO , GO , Node > &  exporter,
const CombineMode  CM,
const bool  restrictedMode = false 
)
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:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *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.

Parameters
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.
void Tpetra::DistObject< Scalar , LO , GO , Node >::doExport ( const SrcDistObject< Scalar, LO, GO, Node > &  source,
const Export< LO , GO , Node > &  exporter,
const CombineMode  CM,
const bool  restrictedMode = false 
)
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:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *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.

Parameters
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.
void Tpetra::DistObject< Scalar , LO , GO , Node >::doExport ( const SrcDistObject< Scalar, LO, GO, Node > &  source,
const Import< LO , GO , Node > &  importer,
const CombineMode  CM,
const bool  restrictedMode = false 
)
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:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *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.

Parameters
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.
bool Tpetra::DistObject< Scalar , LO , GO , Node >::isDistributed ( ) const
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.

virtual Teuchos::RCP<const map_type> Tpetra::DistObject< Scalar , LO , GO , Node >::getMap ( ) const
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 540 of file Tpetra_DistObject_decl.hpp.

void Tpetra::DistObject< Scalar , LO , GO , Node >::print ( std::ostream &  os) const
inherited

Print this object to the given output stream.

We generally assume that all MPI processes can print to the given stream.

virtual std::string Tpetra::DistObject< Scalar , LO , GO , Node >::description ( ) const
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 >.

virtual void Tpetra::DistObject< Scalar , LO , GO , Node >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
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 >.

virtual void Tpetra::DistObject< Scalar , LO , GO , Node >::removeEmptyProcessesInPlace ( const Teuchos::RCP< const map_type > &  newMap)
virtualinherited

Remove processes which contain no entries in this object's Map.

Warning
This method is ONLY for use by experts. We highly recommend using the nonmember function of the same name defined in this file.
We make NO promises of backwards compatibility. This method may change or disappear at any time.

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.

Returns
The object's new Map. Its communicator is a new communicator, distinct from the old Map's communicator, which contains a subset of the processes in the old communicator.
Note
The name differs from Map's method removeEmptyProcesses(), in order to emphasize that the operation on DistObject happens in place, modifying the input, whereas the operation removeEmptyProcess() on Map does not modify the input.
To implementers of DistObject subclasses: The default implementation of this class throws std::logic_error.

Reimplemented in Tpetra::MultiVector< Scalar, LO, GO, Node >.

virtual size_t Tpetra::DistObject< Scalar , LO , GO , Node >::constantNumberOfPackets ( ) const
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 >.

virtual void Tpetra::DistObject< Scalar , LO , GO , Node >::doTransfer ( const SrcDistObject< Scalar, LO, GO, Node > &  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 
)
protectedvirtualinherited

Redistribute data across (MPI) processes.

Parameters
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.
virtual bool Tpetra::DistObject< Scalar , LO , GO , Node >::reallocArraysForNumPacketsPerLid ( const size_t  numExportLIDs,
const size_t  numImportLIDs 
)
protectedvirtualinherited

Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary.

Parameters
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().
Returns
Whether we actually reallocated either of the arrays.
Warning
This is an implementation detail of doTransferNew(). This needs to be protected, but that doesn't mean users should call this method.
virtual void Tpetra::DistObject< Scalar , LO , GO , Node >::doTransferNew ( const SrcDistObject< Scalar, LO, GO, Node > &  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 
)
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.

bool Tpetra::DistObject< Scalar , LO , GO , Node >::reallocImportsIfNeeded ( const size_t  newSize,
const bool  verbose,
const std::string *  prefix 
)
protectedinherited

Reallocate imports_ if needed.

This unfortunately must be declared protected, for the same reason that imports_ is declared protected.

Parameters
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.
Returns
Whether we actually reallocated.

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.

Member Data Documentation

template<class Scalar, class LO, class GO, class Node>
mv_type Tpetra::BlockMultiVector< Scalar, LO, GO, Node >::mv_
protectedinherited

The Tpetra::MultiVector used to represent the data.

Definition at line 689 of file Tpetra_BlockMultiVector_decl.hpp.

Teuchos::RCP<const map_type> Tpetra::DistObject< Scalar , LO , GO , Node >::map_
protectedinherited

The Map over which this object is distributed.

Definition at line 905 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<packet_type*, buffer_device_type> Tpetra::DistObject< Scalar , LO , GO , Node >::imports_
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 918 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<size_t*, buffer_device_type> Tpetra::DistObject< Scalar , LO , GO , Node >::numImportPacketsPerLID_
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 956 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<packet_type*, buffer_device_type> Tpetra::DistObject< Scalar , LO , GO , Node >::exports_
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 963 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<size_t*, buffer_device_type> Tpetra::DistObject< Scalar , LO , GO , Node >::numExportPacketsPerLID_
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 978 of file Tpetra_DistObject_decl.hpp.


The documentation for this class was generated from the following files: