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

Sparse matrix whose entries are small dense square blocks, all of the same dimensions. More...

#include <Tpetra_Experimental_BlockCrsMatrix_decl.hpp>

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

Public Types

typedef Kokkos::Device
< typename
device_type::execution_space,
buffer_memory_space
buffer_device_type
 Kokkos::Device specialization for communication buffers. More...
 
Public typedefs
using scalar_type = Scalar
 The type of entries in the matrix (that is, of each entry in each block). More...
 
using impl_scalar_type = typename BMV::impl_scalar_type
 The implementation type of entries in the matrix. More...
 
typedef LO local_ordinal_type
 The type of local indices. More...
 
typedef GO global_ordinal_type
 The type of global indices. More...
 
typedef Node node_type
 The Node type. More...
 
typedef Node::device_type device_type
 The Kokkos::Device specialization that this class uses. More...
 
typedef
device_type::execution_space 
execution_space
 The Kokkos execution space that this class uses. More...
 
typedef device_type::memory_space memory_space
 The Kokkos memory space that this class uses. More...
 
typedef ::Tpetra::Map< LO, GO,
node_type
map_type
 The implementation of Map that this class uses. More...
 
typedef ::Tpetra::MultiVector
< Scalar, LO, GO, node_type
mv_type
 The implementation of MultiVector that this class uses. More...
 
typedef ::Tpetra::CrsGraph< LO,
GO, node_type
crs_graph_type
 The implementation of CrsGraph that this class uses. More...
 
typedef Kokkos::View
< impl_scalar_type
**, Kokkos::LayoutRight,
device_type,
Kokkos::MemoryTraits
< Kokkos::Unmanaged > > 
little_block_type
 The type used to access nonconst matrix blocks. More...
 
typedef Kokkos::View< const
impl_scalar_type
**, Kokkos::LayoutRight,
device_type,
Kokkos::MemoryTraits
< Kokkos::Unmanaged > > 
const_little_block_type
 The type used to access const matrix blocks. More...
 
typedef Kokkos::View
< impl_scalar_type
*, Kokkos::LayoutRight,
device_type,
Kokkos::MemoryTraits
< Kokkos::Unmanaged > > 
little_vec_type
 The type used to access nonconst vector blocks. More...
 
typedef Kokkos::View< const
impl_scalar_type
*, Kokkos::LayoutRight,
device_type,
Kokkos::MemoryTraits
< Kokkos::Unmanaged > > 
const_little_vec_type
 The type used to access const vector blocks. More...
 
Typedefs
typedef MultiVector< Scalar,
LocalOrdinal, GlobalOrdinal,
Node >::mag_type 
mag_type
 Type of a norm result. More...
 

Public Member Functions

virtual Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm () const
 The communicator over which this matrix is distributed. More...
 
virtual Teuchos::RCP< Node > getNode () const
 The Kokkos Node instance. More...
 
virtual global_size_t getGlobalNumCols () const
 The global number of columns of this matrix. More...
 
virtual size_t getNodeNumCols () const
 The number of columns needed to apply the forward operator on this node. More...
 
virtual GO getIndexBase () const
 The index base for global indices in this matrix. More...
 
virtual global_size_t getGlobalNumEntries () const
 The global number of stored (structurally nonzero) entries. More...
 
virtual size_t getNodeNumEntries () const
 The local number of stored (structurally nonzero) entries. More...
 
virtual size_t getNumEntriesInGlobalRow (GO globalRow) const
 The current number of entries on the calling process in the specified global row. More...
 
virtual global_size_t
TPETRA_DEPRECATED 
getGlobalNumDiags () const
 Number of diagonal entries in the matrix's graph, over all processes in the matrix's communicator. More...
 
virtual size_t TPETRA_DEPRECATED getNodeNumDiags () const
 Number of diagonal entries in the matrix's graph, on the calling process. More...
 
virtual size_t getGlobalMaxNumRowEntries () const
 The maximum number of entries across all rows/columns on all nodes. More...
 
virtual bool hasColMap () const
 Whether this matrix has a well-defined column map. More...
 
virtual bool TPETRA_DEPRECATED isLowerTriangular () const
 Whether the matrix's graph is locally lower triangular. More...
 
virtual bool TPETRA_DEPRECATED isUpperTriangular () const
 Whether the matrix's graph is locally upper triangular. More...
 
virtual bool isLocallyIndexed () const
 Whether matrix indices are locally indexed. More...
 
virtual bool isGloballyIndexed () const
 Whether matrix indices are globally indexed. More...
 
virtual bool isFillComplete () const
 Whether fillComplete() has been called. More...
 
virtual bool supportsRowViews () const
 Whether this object implements getLocalRowView() and getGlobalRowView(). More...
 
Constructors and destructor
 BlockCrsMatrix ()
 Default constructor: Makes an empty block matrix. More...
 
 BlockCrsMatrix (const crs_graph_type &graph, const LO blockSize)
 Constructor that takes a graph and a block size. More...
 
 BlockCrsMatrix (const crs_graph_type &graph, const map_type &domainPointMap, const map_type &rangePointMap, const LO blockSize)
 Constructor that takes a graph, domain and range point Maps, and a block size. More...
 
virtual ~BlockCrsMatrix ()
 Destructor (declared virtual for memory safety). More...
 
Implementation of Tpetra::Operator
Teuchos::RCP< const map_typegetDomainMap () const
 Get the (point) domain Map of this matrix. More...
 
Teuchos::RCP< const map_typegetRangeMap () const
 Get the (point) range Map of this matrix. More...
 
Teuchos::RCP< const map_typegetRowMap () const
 get the (mesh) map for the rows of this block matrix. More...
 
Teuchos::RCP< const map_typegetColMap () const
 get the (mesh) map for the columns of this block matrix. More...
 
global_size_t getGlobalNumRows () const
 get the global number of block rows More...
 
size_t getNodeNumRows () const
 get the local number of block rows More...
 
size_t getNodeMaxNumRowEntries () const
 Maximum number of entries in any row of the matrix, on this process. More...
 
void apply (const mv_type &X, mv_type &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
 For this matrix A, compute Y := beta * Y + alpha * Op(A) * X. More...
 
bool hasTransposeApply () const
 Whether it is valid to apply the transpose or conjugate transpose of this matrix. More...
 
void setAllToScalar (const Scalar &alpha)
 Set all matrix entries equal to alpha. More...
 
Implementation of Teuchos::Describable
std::string description () const
 One-line description of this object. More...
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
 Print a description of this object to the given output stream. More...
 
Implementation of "dual view semantics"
template<class MemorySpace >
void modify ()
 Mark the matrix's values as modified in the given memory space. More...
 
template<class MemorySpace >
bool need_sync () const
 Whether the matrix's values need sync'ing to the given memory space. More...
 
template<class MemorySpace >
void sync ()
 Sync the matrix's values to the given memory space. More...
 
template<class MemorySpace >
auto getValues () -> decltype(val_.template view< typename MemorySpace::memory_space >())
 Get the host or device View of the matrix's values (val_). More...
 
Extraction Methods
virtual void getGlobalRowCopy (GO GlobalRow, const Teuchos::ArrayView< GO > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const
 Get a copy of the given global row's entries. More...
 
virtual void getGlobalRowView (GO GlobalRow, Teuchos::ArrayView< const GO > &indices, Teuchos::ArrayView< const Scalar > &values) const
 Get a constant, nonpersisting, globally indexed view of the given row of the matrix. More...
 
virtual void getLocalDiagCopy (::Tpetra::Vector< Scalar, LO, GO, Node > &diag) const
 Get a copy of the diagonal entries, distributed by the row Map. More...
 
Mathematical methods
virtual void leftScale (const ::Tpetra::Vector< Scalar, LO, GO, Node > &x)
 Scale the RowMatrix on the left with the given Vector x. More...
 
virtual void rightScale (const ::Tpetra::Vector< Scalar, LO, GO, Node > &x)
 Scale the RowMatrix on the right with the given Vector x. More...
 
virtual
typename::Tpetra::RowMatrix
< Scalar, LO, GO, Node >
::mag_type 
getFrobeniusNorm () const
 The Frobenius norm of the matrix. More...
 
Matrix query methods
virtual size_t getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const =0
 The current number of entries on the calling process in the specified global row. More...
 
virtual size_t getNumEntriesInLocalRow (LocalOrdinal localRow) const =0
 The current number of entries on the calling process in the specified local row. More...
 
Extraction Methods
virtual void getGlobalRowCopy (GlobalOrdinal GlobalRow, const Teuchos::ArrayView< GlobalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const =0
 Get a copy of the given global row's entries. More...
 
virtual void getLocalRowCopy (LocalOrdinal LocalRow, const Teuchos::ArrayView< LocalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const =0
 Get a copy of the given local row's entries. More...
 
virtual void getGlobalRowView (GlobalOrdinal GlobalRow, Teuchos::ArrayView< const GlobalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const =0
 Get a constant, nonpersisting, globally indexed view of the given row of the matrix. More...
 
virtual void getLocalRowView (LocalOrdinal LocalRow, Teuchos::ArrayView< const LocalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const =0
 Get a constant, nonpersisting, locally indexed view of the given row of the matrix. More...
 
virtual LocalOrdinal getLocalRowViewRaw (const LocalOrdinal lclRow, LocalOrdinal &numEnt, const LocalOrdinal *&lclColInds, const Scalar *&vals) const
 Get a constant, nonpersisting, locally indexed view of the given row of the matrix, using "raw" pointers instead of Teuchos::ArrayView. More...
 
virtual void getLocalDiagCopy (Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const =0
 Get a copy of the diagonal entries, distributed by the row Map. More...
 
Mathematical methods
virtual void leftScale (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
 Scale the matrix on the left with the given Vector. More...
 
virtual void rightScale (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
 Scale the matrix on the right with the given Vector. More...
 
virtual Teuchos::RCP
< RowMatrix< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > > 
add (const Scalar &alpha, const RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) const
 Return a new RowMatrix which is the result of beta*this + alpha*A. More...
 
Implementation of Packable interface
virtual void pack (const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Teuchos::Array< char > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets, Distributor &distor) const
 Pack this object's data for an Import or Export. More...
 
Pure virtual functions to be overridden by subclasses.
virtual void apply (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const =0
 Computes the operator-multivector application. More...
 
Public methods for redistributing data
void doImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
 Import data into this object using an Import object ("forward mode"). More...
 
void doImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
 Import data into this object using an Export object ("reverse mode"). More...
 
void doExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
 Export data into this object using an Export object ("forward mode"). More...
 
void doExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
 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...
 
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 Types

typedef char packet_type
 Implementation detail; tells. More...
 
enum  ReverseOption
 Whether the data transfer should be performed in forward or reverse mode. More...
 
typedef device_type::memory_space buffer_memory_space
 Kokkos memory space for communication buffers. More...
 

Protected Member Functions

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)
 Redistribute data across memory images. More...
 
virtual bool reallocArraysForNumPacketsPerLid (const size_t numExportLIDs, const size_t numImportLIDs)
 Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary. More...
 
virtual void createViews () const
 Hook for creating a const view. More...
 
virtual void createViewsNonConst (KokkosClassic::ReadWriteOption rwo)
 Hook for creating a nonconst view. More...
 
virtual void releaseViews () const
 Hook for releasing views. More...
 
bool reallocImportsIfNeeded (const size_t newSize, const bool debug=false)
 Reallocate imports_ if needed. 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)
 
virtual void copyAndPermute (const ::Tpetra::SrcDistObject &source, size_t numSameIDs, const Teuchos::ArrayView< const LO > &permuteToLIDs, const Teuchos::ArrayView< const LO > &permuteFromLIDs)
 
virtual void packAndPrepare (const ::Tpetra::SrcDistObject &source, const Teuchos::ArrayView< const LO > &exportLIDs, Teuchos::Array< packet_type > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets,::Tpetra::Distributor &distor)
 
virtual void unpackAndCombine (const Teuchos::ArrayView< const LO > &importLIDs, const Teuchos::ArrayView< const packet_type > &imports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t constantNumPackets,::Tpetra::Distributor &distor,::Tpetra::CombineMode CM)
 
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.")

virtual bool checkSizes (const SrcDistObject &source)=0
 Compare the source and target (this) objects for compatibility. More...
 
virtual bool useNewInterface ()
 Whether the subclass implements the "old" or "new" (Kokkos-friendly) interface. More...
 
virtual void copyAndPermute (const SrcDistObject &source, size_t numSameIDs, const Teuchos::ArrayView< const local_ordinal_type > &permuteToLIDs, const Teuchos::ArrayView< const local_ordinal_type > &permuteFromLIDs)
 Perform copies and permutations that are local to this process. More...
 
virtual void copyAndPermuteNew (const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, device_type > &permuteFromLIDs)
 
virtual void packAndPrepare (const SrcDistObject &source, const Teuchos::ArrayView< const local_ordinal_type > &exportLIDs, Teuchos::Array< packet_type > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets, Distributor &distor)
 Perform any packing or preparation required for communication. More...
 
virtual void packAndPrepareNew (const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, const Kokkos::DualView< size_t *, buffer_device_type > &numPacketsPerLID, size_t &constantNumPackets, Distributor &distor)
 
virtual void unpackAndCombine (const Teuchos::ArrayView< const local_ordinal_type > &importLIDs, const Teuchos::ArrayView< const packet_type > &imports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t constantNumPackets, Distributor &distor, CombineMode CM)
 Perform any unpacking and combining after communication (old version that uses Teuchos memory management classes to hold data). More...
 
virtual void unpackAndCombineNew (const Kokkos::DualView< const local_ordinal_type *, device_type > &importLIDs, const Kokkos::DualView< const packet_type *, buffer_device_type > &imports, const Kokkos::DualView< const size_t *, buffer_device_type > &numPacketsPerLID, const size_t constantNumPackets, Distributor &distor, const CombineMode CM)
 Perform any unpacking and combining after communication (new version that uses Kokkos data structures to hold data). More...
 

Protected Attributes

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

Block operations

LO getBlockSize () const
 The number of degrees of freedom per mesh point. More...
 
virtual Teuchos::RCP< const
::Tpetra::RowGraph< LO, GO,
Node > > 
getGraph () const
 Get the (mesh) graph. More...
 
const crs_graph_typegetCrsGraph () const
 
void applyBlock (const BlockMultiVector< Scalar, LO, GO, Node > &X, BlockMultiVector< Scalar, LO, GO, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, const Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), const Scalar beta=Teuchos::ScalarTraits< Scalar >::zero())
 Version of apply() that takes BlockMultiVector input and output. More...
 
void gaussSeidelCopy (MultiVector< Scalar, LO, GO, Node > &X, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &B, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &D, const Scalar &dampingFactor, const ESweepDirection direction, const int numSweeps, const bool zeroInitialGuess) const
 Version of gaussSeidel(), with fewer requirements on X. More...
 
void reorderedGaussSeidelCopy (::Tpetra::MultiVector< Scalar, LO, GO, Node > &X, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &B, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &D, const Teuchos::ArrayView< LO > &rowIndices, const Scalar &dampingFactor, const ESweepDirection direction, const int numSweeps, const bool zeroInitialGuess) const
 Version of reorderedGaussSeidel(), with fewer requirements on X. More...
 
void localGaussSeidel (const BlockMultiVector< Scalar, LO, GO, Node > &Residual, BlockMultiVector< Scalar, LO, GO, Node > &Solution, const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &D_inv, const Scalar &omega, const ESweepDirection direction) const
 Local Gauss-Seidel solve, given a factorized diagonal. More...
 
LO replaceLocalValues (const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
 Replace values at the given (mesh, i.e., block) column indices, in the given (mesh, i.e., block) row. More...
 
LO sumIntoLocalValues (const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
 Sum into values at the given (mesh, i.e., block) column indices, in the given (mesh, i.e., block) row. More...
 
LO getLocalRowView (const LO localRowInd, const LO *&colInds, Scalar *&vals, LO &numInds) const
 Get a view of the (mesh, i.e., block) row, using local (mesh, i.e., block) indices. More...
 
void getLocalRowView (LO LocalRow, Teuchos::ArrayView< const LO > &indices, Teuchos::ArrayView< const Scalar > &values) const
 Not implemented. More...
 
void getLocalRowCopy (LO LocalRow, const Teuchos::ArrayView< LO > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const
 Not implemented. More...
 
little_block_type getLocalBlock (const LO localRowInd, const LO localColInd) const
 
LO getLocalRowOffsets (const LO localRowInd, ptrdiff_t offsets[], const LO colInds[], const LO numColInds) const
 Get relative offsets corresponding to the given rows, given by local row index. More...
 
LO replaceLocalValuesByOffsets (const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
 Like replaceLocalValues, but avoids computing row offsets. More...
 
LO sumIntoLocalValuesByOffsets (const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
 Like sumIntoLocalValues, but avoids computing row offsets. More...
 
size_t getNumEntriesInLocalRow (const LO localRowInd) const
 Return the number of entries in the given row on the calling process. More...
 
bool localError () const
 Whether this object had an error on the calling process. More...
 
std::string errorMessages () const
 The current stream of error messages. More...
 
void getLocalDiagOffsets (const Kokkos::View< size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
 Get offsets of the diagonal entries in the matrix. More...
 
void TPETRA_DEPRECATED getLocalDiagOffsets (Teuchos::ArrayRCP< size_t > &offsets) const
 DEPRECATED overload of this method that writes offsets to a Teuchos::ArrayRCP instead of a Kokkos::View. More...
 
void getLocalDiagCopy (const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &diag, const Kokkos::View< const size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
 Variant of getLocalDiagCopy() that uses precomputed offsets and puts diagonal blocks in a 3-D Kokkos::View. More...
 
void getLocalDiagCopy (const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &diag, const Teuchos::ArrayView< const size_t > &offsets) const
 Variant of getLocalDiagCopy() that uses precomputed offsets and puts diagonal blocks in a 3-D Kokkos::View. More...
 
void TPETRA_DEPRECATED getLocalDiagCopy (BlockCrsMatrix< Scalar, LO, GO, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const
 Variant of getLocalDiagCopy() that uses precomputed offsets. More...
 
LO absMaxLocalValues (const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
 Like sumIntoLocalValues, but for the ABSMAX combine mode. More...
 
LO absMaxLocalValuesByOffsets (const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
 Like sumIntoLocalValuesByOffsets, but for the ABSMAX combine mode. More...
 

Detailed Description

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
class Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >

Sparse matrix whose entries are small dense square blocks, all of the same dimensions.

Template Parameters
ScalarThe type of the numerical entries of the matrix. (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.

Please read the documentation of BlockMultiVector first.

This class implements a sparse matrix whose entries are small dense square blocks, all of the same dimensions. The intended application is to store the discretization of a partial differential equation with multiple degrees of freedom per mesh point, where all mesh points have the same number of degrees of freedom. This class stores values associated with the degrees of freedom of a single mesh point contiguously, in a getBlockSize() by getBlockSize() block, in row-major format. The matrix's graph represents the mesh points, with one entry per mesh point. This saves storage over using a CrsMatrix, which requires one graph entry per matrix entry.

This class requires a fill-complete Tpetra::CrsGraph for construction. Thus, it has a row Map and a column Map already. As a result, BlockCrsMatrix only needs to provide access using local indices. Access using local indices is faster anyway, since conversion from global to local indices requires a hash table lookup per index. Users are responsible for converting from global to local indices if necessary. Please be aware that the row Map and column Map may differ, so you may not use local row and column indices interchangeably.

We reserve the right to change the block layout in the future. Best practice is to use this class' little_block_type and const_little_block_type typedefs to access blocks, since their types offer access to entries in a layout-independent way. These two typedefs are both Kokkos::View specializations.

Here is an example of how to fill into this object using raw-pointer views.

int err = 0;
// At least one entry, so &offsets[0] always makes sense.
Teuchos::Array<ptrdiff_t> offsets (1);
for (LO localRowInd = 0; localRowInd < localNumRows; ++localRowInd) {
// Get a view of the current row.
// You may modify the values, but not the column indices.
const LO* localColInds;
Scalar* vals;
LO numEntries;
err = A.getLocalRowView (localRowInd, localColInds, vals, numEntries);
if (err != 0) {
break;
}
// Modify the entries in the current row.
for (LO k = 0; k < numEntries; ++k) {
Scalar* const curBlock = vals[blockSize * blockSize * k];
// Blocks are stored in row-major format.
for (LO j = 0; j < blockSize; ++j) {
for (LO i = 0; i < blockSize; ++i) {
const Scalar curVal = &curBlock[i + j * blockSize];
// Some function f of the current value and mesh point
curBlock[i + j * blockSize] = f (curVal, localColInds[k], ...);
}
}
}
}

Definition at line 138 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

Member Typedef Documentation

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef char Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::packet_type
protected

Implementation detail; tells.

Definition at line 149 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
using Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::scalar_type = Scalar

The type of entries in the matrix (that is, of each entry in each block).

Definition at line 156 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
using Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::impl_scalar_type = typename BMV::impl_scalar_type

The implementation type of entries in the matrix.

Letting scalar_type and impl_scalar_type differ helps this class work correctly for Scalar types like std::complex<T>, which lack the necessary CUDA device macros and volatile overloads to work correctly with Kokkos.

Definition at line 164 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef LO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::local_ordinal_type

The type of local indices.

Definition at line 167 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef GO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::global_ordinal_type

The type of global indices.

Definition at line 169 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef Node Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::node_type

The Node type.

Prefer device_type, execution_space, and memory_space (see below), which relate directly to Kokkos.

Definition at line 174 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef Node::device_type Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::device_type

The Kokkos::Device specialization that this class uses.

Definition at line 177 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef device_type::execution_space Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::execution_space

The Kokkos execution space that this class uses.

Definition at line 179 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef device_type::memory_space Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::memory_space

The Kokkos memory space that this class uses.

Definition at line 181 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef ::Tpetra::Map<LO, GO, node_type> Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::map_type

The implementation of Map that this class uses.

Definition at line 184 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef ::Tpetra::MultiVector<Scalar, LO, GO, node_type> Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::mv_type

The implementation of MultiVector that this class uses.

Definition at line 186 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef ::Tpetra::CrsGraph<LO, GO, node_type> Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::crs_graph_type

The implementation of CrsGraph that this class uses.

Definition at line 188 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef Kokkos::View<impl_scalar_type**, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::little_block_type

The type used to access nonconst matrix blocks.

Definition at line 195 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef Kokkos::View<const impl_scalar_type**, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::const_little_block_type

The type used to access const matrix blocks.

Definition at line 201 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef Kokkos::View<impl_scalar_type*, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::little_vec_type

The type used to access nonconst vector blocks.

Definition at line 207 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef Kokkos::View<const impl_scalar_type*, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::const_little_vec_type

The type used to access const vector blocks.

Definition at line 213 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
typedef MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::mag_type Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type
inherited

Type of a norm result.

This is usually the same as the type of the magnitude (absolute value) of Scalar, but may differ for certain Scalar types.

Definition at line 107 of file Tpetra_RowMatrix_decl.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::buffer_memory_space
protectedinherited

Kokkos memory space for communication buffers.

See #1088 for why this is not just device_type::memory_space.

Definition at line 693 of file Tpetra_DistObject_decl.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::buffer_device_type
inherited

Kokkos::Device specialization for communication buffers.

See #1088 for why this is not just device_type::device_type.

This needs to be public so that I can declare functions like packAndPrepareWithOwningPIDs.

Warning
This is an implementation detail. DO NOT DEPEND ON IT. It may disappear or change at any time.

Definition at line 710 of file Tpetra_DistObject_decl.hpp.

Member Enumeration Documentation

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
enum Tpetra::Classes::DistObject::ReverseOption
protectedinherited

Whether the data transfer should be performed in forward or reverse mode.

"Reverse mode" means calling doExport() with an Import object, or calling doImport() with an Export object. "Forward mode" means calling doExport() with an Export object, or calling doImport() with an Import object.

Definition at line 605 of file Tpetra_DistObject_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar , class LO , class GO , class Node >
Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::BlockCrsMatrix ( )

Default constructor: Makes an empty block matrix.

Definition at line 658 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO, class GO , class Node >
Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::BlockCrsMatrix ( const crs_graph_type graph,
const LO  blockSize 
)

Constructor that takes a graph and a block size.

The graph represents the mesh. This constructor computes the point Maps corresponding to the given graph's domain and range Maps. If you already have those point Maps, it is better to call the four-argument constructor.

Parameters
graph[in] A fill-complete graph.
blockSize[in] Number of degrees of freedom per mesh point.

Definition at line 673 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO, class GO , class Node >
Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::BlockCrsMatrix ( const crs_graph_type graph,
const map_type domainPointMap,
const map_type rangePointMap,
const LO  blockSize 
)

Constructor that takes a graph, domain and range point Maps, and a block size.

The graph represents the mesh. This constructor uses the given domain and range point Maps, rather than computing them. The given point Maps must be the same as the above two-argument constructor would have computed.

Definition at line 731 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::~BlockCrsMatrix ( )
inlinevirtual

Destructor (declared virtual for memory safety).

Definition at line 246 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

Member Function Documentation

template<class Scalar , class LO , class GO , class Node >
Teuchos::RCP< const typename BlockCrsMatrix< Scalar, LO, GO, Node >::map_type > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getDomainMap ( ) const
virtual

Get the (point) domain Map of this matrix.

Implements Tpetra::Classes::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 791 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
Teuchos::RCP< const typename BlockCrsMatrix< Scalar, LO, GO, Node >::map_type > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getRangeMap ( ) const
virtual

Get the (point) range Map of this matrix.

Implements Tpetra::Classes::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 800 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
Teuchos::RCP< const typename BlockCrsMatrix< Scalar, LO, GO, Node >::map_type > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getRowMap ( ) const
virtual

get the (mesh) map for the rows of this block matrix.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 809 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
Teuchos::RCP< const typename BlockCrsMatrix< Scalar, LO, GO, Node >::map_type > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getColMap ( ) const
virtual

get the (mesh) map for the columns of this block matrix.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 817 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
global_size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getGlobalNumRows ( ) const
virtual

get the global number of block rows

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 825 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getNodeNumRows ( ) const
virtual

get the local number of block rows

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 833 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getNodeMaxNumRowEntries ( ) const
virtual

Maximum number of entries in any row of the matrix, on this process.

Precondition
Subclasses reserve the right to impose preconditions on the matrix's state.

This method only uses the matrix's graph. Explicitly stored zeros count as "entries."

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 841 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO , class GO , class Node >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::apply ( const mv_type X,
mv_type Y,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
Scalar  alpha = Teuchos::ScalarTraits<Scalar>::one (),
Scalar  beta = Teuchos::ScalarTraits<Scalar>::zero () 
) const

For this matrix A, compute Y := beta * Y + alpha * Op(A) * X.

Op(A) is A if mode is Teuchos::NO_TRANS, the transpose of A if mode is Teuchos::TRANS, and the conjugate transpose of A if mode is Teuchos::CONJ_TRANS.

If alpha is zero, ignore X's entries on input; if beta is zero, ignore Y's entries on input. This follows the BLAS convention, and only matters if X resp. Y have Inf or NaN entries.

Definition at line 849 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
bool Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::hasTransposeApply ( ) const
inlinevirtual

Whether it is valid to apply the transpose or conjugate transpose of this matrix.

Reimplemented from Tpetra::Classes::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 290 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar, class LO , class GO , class Node >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::setAllToScalar ( const Scalar &  alpha)

Set all matrix entries equal to alpha.

Definition at line 942 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
std::string Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::description ( ) const
virtual

One-line description of this object.

Reimplemented from Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3261 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel 
) const
virtual

Print a description of this object to the given output stream.

Parameters
out[out] Output stream to which to print. Valid values include Teuchos::VERB_DEFAULT, Teuchos::VERB_NONE, Teuchos::VERB_LOW, Teuchos::VERB_MEDIUM, Teuchos::VERB_HIGH, and Teuchos::VERB_EXTREME.
verbLevel[in] Verbosity level at which to print.
Warning
If verbLevel is Teuchos::VERB_EXTREME, this method has collective semantics over the matrix's communicator.

The following pseudocode shows how to wrap your std::ostream object in a Teuchos::FancyOStream, and pass it into this method:

// ...
std::ostream& yourObject = ...;
Teuchos::RCP<Teuchos::FancyOStream> wrappedStream =
Teuchos::getFancyOStream (Teuchos::rcpFromRef (yourObject));
const Teuchos::EVerbosityLevel verbLevel = ...;
A.describe (*wrappedStream, verbLevel);

Reimplemented from Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3285 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
LO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getBlockSize ( ) const
inline

The number of degrees of freedom per mesh point.

Definition at line 338 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar , class LO , class GO , class Node >
Teuchos::RCP< const ::Tpetra::RowGraph< LO, GO, Node > > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getGraph ( ) const
virtual
template<class Scalar, class LO, class GO, class Node>
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::applyBlock ( const BlockMultiVector< Scalar, LO, GO, Node > &  X,
BlockMultiVector< Scalar, LO, GO, Node > &  Y,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
const Scalar  alpha = Teuchos::ScalarTraits<Scalar>::one (),
const Scalar  beta = Teuchos::ScalarTraits<Scalar>::zero () 
)

Version of apply() that takes BlockMultiVector input and output.

This method is deliberately not marked const, because it may do lazy initialization of temporary internal block multivectors.

Definition at line 914 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::gaussSeidelCopy ( MultiVector< Scalar, LO, GO, Node > &  X,
const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &  B,
const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &  D,
const Scalar &  dampingFactor,
const ESweepDirection  direction,
const int  numSweeps,
const bool  zeroInitialGuess 
) const

Version of gaussSeidel(), with fewer requirements on X.

Not Implemented

Definition at line 1216 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::reorderedGaussSeidelCopy ( ::Tpetra::MultiVector< Scalar, LO, GO, Node > &  X,
const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &  B,
const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &  D,
const Teuchos::ArrayView< LO > &  rowIndices,
const Scalar &  dampingFactor,
const ESweepDirection  direction,
const int  numSweeps,
const bool  zeroInitialGuess 
) const

Version of reorderedGaussSeidel(), with fewer requirements on X.

Not Implemented

Definition at line 1234 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::localGaussSeidel ( const BlockMultiVector< Scalar, LO, GO, Node > &  Residual,
BlockMultiVector< Scalar, LO, GO, Node > &  Solution,
const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &  D_inv,
const Scalar &  omega,
const ESweepDirection  direction 
) const

Local Gauss-Seidel solve, given a factorized diagonal.

Parameters
Residual[in] The "residual" (right-hand side) block (multi)vector
Solution[in/out] On input: the initial guess / current approximate solution. On output: the new approximate solution.
D_inv[in] Block diagonal, the explicit inverse of this matrix's block diagonal (possibly modified for algorithmic reasons).
omega[in] (S)SOR relaxation coefficient
direction[in] Forward, Backward, or Symmetric.

One may access block i in D_inv using the following code:

auto D_ii = Kokkos::subview(D_inv, i, Kokkos::ALL(), Kokkos::ALL());

The resulting block is b x b, where b = this->getBlockSize().

Definition at line 1105 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO , class Node >
LO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::replaceLocalValues ( const LO  localRowInd,
const LO  colInds[],
const Scalar  vals[],
const LO  numColInds 
) const

Replace values at the given (mesh, i.e., block) column indices, in the given (mesh, i.e., block) row.

Parameters
localRowInd[in] Local mesh (i.e., block) index of the row in which to replace.
colInds[in] Local mesh (i.e., block) column ind{ex,ices} at which to replace values. colInds[k] is the local column index whose new values start at vals[getBlockSize() * getBlockSize() * k], and colInds has length at least numColInds. This method will only access the first numColInds entries of colInds.
vals[in] The new values to use at the given column indices. Values for each block are stored contiguously, in row major layout, with no padding between rows or between blocks. Thus, if b = getBlockSize(), then vals[k*b*b] .. vals[(k+1)*b*b-1] are the values to use for block colInds[k].
numColInds[in] The number of entries of colInds.
Returns
The number of valid entries of colInds. colInds[k] is valid if and only if it is a valid local mesh (i.e., block) column index. This method succeeded if and only if the return value equals the input argument numColInds.

Definition at line 981 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO , class Node >
LO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::sumIntoLocalValues ( const LO  localRowInd,
const LO  colInds[],
const Scalar  vals[],
const LO  numColInds 
) const

Sum into values at the given (mesh, i.e., block) column indices, in the given (mesh, i.e., block) row.

Parameters
localRowInd[in] Local mesh (i.e., block) index of the row in which to sum.
colInds[in] Local mesh (i.e., block) column ind{ex,ices} at which to sum. colInds[k] is the local column index whose new values start at vals[getBlockSize() * getBlockSize() * k], and colInds has length at least numColInds. This method will only access the first numColInds entries of colInds.
vals[in] The new values to sum in at the given column indices. Values for each block are stored contiguously, in row major layout, with no padding between rows or between blocks. Thus, if b = getBlockSize(), then vals[k*b*b] .. vals[(k+1)*b*b-1] are the values to use for block colInds[k].
numColInds[in] The number of entries of colInds.
Returns
The number of valid entries of colInds. colInds[k] is valid if and only if it is a valid local mesh (i.e., block) column index. This method succeeded if and only if the return value equals the input argument numColInds.

Definition at line 1414 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO , class Node >
LO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getLocalRowView ( const LO  localRowInd,
const LO *&  colInds,
Scalar *&  vals,
LO &  numInds 
) const

Get a view of the (mesh, i.e., block) row, using local (mesh, i.e., block) indices.

This matrix has a graph, and we assume that the graph is fill complete on input to the matrix's constructor. Thus, the matrix has a column Map, and it stores column indices as local indices. This means you can view the column indices as local indices directly. However, you may not view them as global indices directly, since the column indices are not stored as global indices in the graph.

Parameters
localRowInd[in] Local (mesh, i.e., block) row index.
colInds[out] If localRowInd is valid on the calling process, then on output, this is a pointer to the local (mesh, i.e., block) column indices in the given (mesh, i.e., block) row. If localRowInd is not valid, then this is undefined. (Please check the return value of this method.)
vals[out] If localRowInd is valid on the calling process, then on output, this is a pointer to the row's values. If localRowInd is not valid, then this is undefined. (Please check the return value of this method.)
numInds[in] The number of (mesh, i.e., block) indices in colInds on output.
Returns
0 if localRowInd is valid, else Teuchos::OrdinalTraits<LO>::invalid().

Definition at line 1491 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO , class Node >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getLocalRowView ( LO  LocalRow,
Teuchos::ArrayView< const LO > &  indices,
Teuchos::ArrayView< const Scalar > &  values 
) const

Not implemented.

Definition at line 3715 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO , class Node >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getLocalRowCopy ( LO  LocalRow,
const Teuchos::ArrayView< LO > &  Indices,
const Teuchos::ArrayView< Scalar > &  Values,
size_t &  NumEntries 
) const

Not implemented.

Definition at line 1539 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO, class GO , class Node >
LO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getLocalRowOffsets ( const LO  localRowInd,
ptrdiff_t  offsets[],
const LO  colInds[],
const LO  numColInds 
) const

Get relative offsets corresponding to the given rows, given by local row index.

The point of this method is to precompute the results of searching for the offsets corresponding to the given column indices. You may then reuse these search results in replaceLocalValuesByOffsets or sumIntoLocalValuesByOffsets.

Offsets are block offsets; they are for column indices, not for values.

Parameters
localRowInd[in] Local index of the row.
offsets[out] On output: relative offsets corresponding to the given column indices. Must have at least numColInds entries.
colInds[in] The local column indices for which to compute offsets. Must have at least numColInds entries. This method will only read the first numColsInds entries.
numColInds[in] Number of entries in colInds to read.
Returns
The number of valid column indices in colInds. This method succeeded if and only if the return value equals the input argument numColInds.

Definition at line 1565 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO , class Node >
LO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::replaceLocalValuesByOffsets ( const LO  localRowInd,
const ptrdiff_t  offsets[],
const Scalar  vals[],
const LO  numOffsets 
) const

Like replaceLocalValues, but avoids computing row offsets.

Returns
The number of valid column indices in colInds. This method succeeded if and only if the return value equals the input argument numColInds.

Definition at line 1598 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO , class Node >
LO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::sumIntoLocalValuesByOffsets ( const LO  localRowInd,
const ptrdiff_t  offsets[],
const Scalar  vals[],
const LO  numOffsets 
) const

Like sumIntoLocalValues, but avoids computing row offsets.

Returns
The number of valid column indices in colInds. This method succeeded if and only if the return value equals the input argument numColInds.

Definition at line 1676 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO, class GO , class Node >
size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getNumEntriesInLocalRow ( const LO  localRowInd) const

Return the number of entries in the given row on the calling process.

If the given local row index is invalid, this method (sensibly) returns zero, since the calling process trivially does not own any entries in that row.

Definition at line 1717 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
bool Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::localError ( ) const
inline

Whether this object had an error on the calling process.

Import and Export operations using this object as the target of the Import or Export may incur local errors, if some process encounters an LID in its list which is not a valid mesh row local index on that process. In that case, we don't want to throw an exception, because not all processes may throw an exception; this can result in deadlock or put Tpetra in an incorrect state, due to lack of consistency across processes. Instead, we set a local error flag and ignore the incorrect data. When unpacking, we do the same with invalid column indices. If you want to check whether some process experienced an error, you must do a reduction or all-reduce over this flag. Every time you initiate a new Import or Export with this object as the target, we clear this flag. (Note to developers: we clear it at the beginning of checkSizes().)

Definition at line 598 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
std::string Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::errorMessages ( ) const
inline

The current stream of error messages.

This is only nonempty on the calling process if localError() returns true. In that case, it stores a stream of human-readable, endline-separated error messages encountered during an Import or Export cycle. Every time you initiate a new Import or Export with this object as the target, we clear this stream. (Note to developers: we clear it at the beginning of checkSizes().)

If you want to print this, you are responsible for ensuring that it is valid for the calling MPI process to print to whatever output stream you use. On some MPI implementations, you may need to send the string to Process 0 for printing.

Definition at line 616 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getLocalDiagOffsets ( const Kokkos::View< size_t *, device_type, Kokkos::MemoryUnmanaged > &  offsets) const

Get offsets of the diagonal entries in the matrix.

Warning
This method is only for expert users.
We make no promises about backwards compatibility for this method. It may disappear or change at any time.
This method must be called collectively. We reserve the right to do extra checking in a debug build that will require collectives.
Precondition
The matrix must be locally indexed (which means that it has a column Map).
All diagonal entries of the matrix's graph must be populated on this process. Results are undefined otherwise.
Postcondition
offsets.extent(0) == getNodeNumRows()

This method creates an array of offsets of the local diagonal entries in the matrix. This array is suitable for use in the two-argument version of getLocalDiagCopy(). However, its contents are not defined in any other context. For example, you should not rely on offsets(i) being the index of the diagonal entry in the views returned by getLocalRowView(). This may be the case, but it need not be. (For example, we may choose to optimize the lookups down to the optimized storage level, in which case the offsets will be computed with respect to the underlying storage format, rather than with respect to the views.)

If the matrix has a const ("static") graph, and if that graph is fill complete, then the offsets array remains valid through calls to fillComplete() and resumeFill(). "Invalidates" means that you must call this method again to recompute the offsets.

Definition at line 1058 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void TPETRA_DEPRECATED Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getLocalDiagOffsets ( Teuchos::ArrayRCP< size_t > &  offsets) const

DEPRECATED overload of this method that writes offsets to a Teuchos::ArrayRCP instead of a Kokkos::View.

Please use the version of this method directly above, that writes offsets a Kokkos::View instead of to a Teuchos::ArrayRCP.

Definition at line 1067 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getLocalDiagCopy ( const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &  diag,
const Kokkos::View< const size_t *, device_type, Kokkos::MemoryUnmanaged > &  offsets 
) const

Variant of getLocalDiagCopy() that uses precomputed offsets and puts diagonal blocks in a 3-D Kokkos::View.

Parameters
diag[out] On input: Must be preallocated, with dimensions at least (number of diagonal blocks on the calling process) x getBlockSize() x getBlockSize(). On output: the diagonal blocks. Leftmost index is "which block," then the row index within a block, then the column index within a block.

This method uses the offsets of the diagonal entries, as precomputed by getLocalDiagOffsets(), to speed up copying the diagonal of the matrix.

Definition at line 1284 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getLocalDiagCopy ( const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &  diag,
const Teuchos::ArrayView< const size_t > &  offsets 
) const

Variant of getLocalDiagCopy() that uses precomputed offsets and puts diagonal blocks in a 3-D Kokkos::View.

Parameters
diag[out] On input: Must be preallocated, with dimensions at least (number of diagonal blocks on the calling process) x getBlockSize() x getBlockSize(). On output: the diagonal blocks. Leftmost index is "which block," then the row index within a block, then the column index within a block.

This method uses the offsets of the diagonal entries, as precomputed by getLocalDiagOffsets(), to speed up copying the diagonal of the matrix.

Definition at line 1334 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO, class Node>
void TPETRA_DEPRECATED Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getLocalDiagCopy ( BlockCrsMatrix< Scalar, LO, GO, Node > &  diag,
const Teuchos::ArrayView< const size_t > &  offsets 
) const

Variant of getLocalDiagCopy() that uses precomputed offsets.

Warning
This overload of the method is DEPRECATED. Call the overload above that returns the diagonal blocks as a 3-D Kokkos::View.

This method uses the offsets of the diagonal entries, as precomputed by getLocalDiagOffsets(), to speed up copying the diagonal of the matrix.

If the matrix has a const ("static") graph, and if that graph is fill complete, then the offsets array remains valid through calls to fillComplete() and resumeFill().

Definition at line 1253 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO , class Node >
LO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::absMaxLocalValues ( const LO  localRowInd,
const LO  colInds[],
const Scalar  vals[],
const LO  numColInds 
) const
protected

Like sumIntoLocalValues, but for the ABSMAX combine mode.

Definition at line 1371 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO , class Node >
LO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::absMaxLocalValuesByOffsets ( const LO  localRowInd,
const ptrdiff_t  offsets[],
const Scalar  vals[],
const LO  numOffsets 
) const
protected

Like sumIntoLocalValuesByOffsets, but for the ABSMAX combine mode.

Definition at line 1637 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
template<class MemorySpace >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::modify ( )
inline

Mark the matrix's values as modified in the given memory space.

Definition at line 896 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
template<class MemorySpace >
bool Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::need_sync ( ) const
inline

Whether the matrix's values need sync'ing to the given memory space.

Definition at line 909 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
template<class MemorySpace >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::sync ( )
inline

Sync the matrix's values to the given memory space.

Definition at line 922 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LO = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
template<class MemorySpace >
auto Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getValues ( ) -> decltype (val_.template view<typename MemorySpace::memory_space> ())
inline

Get the host or device View of the matrix's values (val_).

Warning
This is for EXPERT USE ONLY. We make NO PROMISES OF BACKWARDS COMPATIBILITY. YOU HAVE BEEN WARNED. CAVEAT LECTOR!!!
Template Parameters
MemorySpaceMemory space for which to get the Kokkos::View.

This is not const, because we reserve the right to do lazy allocation on device / "fast" memory. Host / "slow" memory allocations generally are not lazy; that way, the host fill interface always works in a thread-parallel context without needing to synchronize on the allocation.

Definition at line 948 of file Tpetra_Experimental_BlockCrsMatrix_decl.hpp.

template<class Scalar , class LO , class GO , class Node >
Teuchos::RCP< const Teuchos::Comm< int > > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getComm ( ) const
virtual

The communicator over which this matrix is distributed.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3538 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
Teuchos::RCP< Node > Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getNode ( ) const
virtual
template<class Scalar , class LO , class GO , class Node >
global_size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getGlobalNumCols ( ) const
virtual

The global number of columns of this matrix.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3555 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getNodeNumCols ( ) const
virtual

The number of columns needed to apply the forward operator on this node.

This is the same as the number of elements listed in the column Map. It is not necessarily the same as the number of domain Map elements owned by the calling process.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3563 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
GO Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getIndexBase ( ) const
virtual

The index base for global indices in this matrix.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3571 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
global_size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getGlobalNumEntries ( ) const
virtual

The global number of stored (structurally nonzero) entries.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3579 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getNodeNumEntries ( ) const
virtual

The local number of stored (structurally nonzero) entries.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3587 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO, class Node >
size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getNumEntriesInGlobalRow ( GO  globalRow) const
virtual

The current number of entries on the calling process in the specified global row.

Note that if the row Map is overlapping, then the calling process might not necessarily store all the entries in the row. Some other process might have the rest of the entries.

Returns
Teuchos::OrdinalTraits<size_t>::invalid() if the specified global row does not belong to this graph, else the number of entries.

Definition at line 3595 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
global_size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getGlobalNumDiags ( ) const
virtual

Number of diagonal entries in the matrix's graph, over all processes in the matrix's communicator.

Precondition
! this->isFillActive()
Warning
This method is DEPRECATED. DO NOT CALL IT. It may go away at any time.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3603 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getNodeNumDiags ( ) const
virtual

Number of diagonal entries in the matrix's graph, on the calling process.

Precondition
! this->isFillActive()
Warning
This method is DEPRECATED. DO NOT CALL IT. It may go away at any time.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3612 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
size_t Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getGlobalMaxNumRowEntries ( ) const
virtual

The maximum number of entries across all rows/columns on all nodes.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3621 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
bool Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::hasColMap ( ) const
virtual

Whether this matrix has a well-defined column map.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3629 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
bool Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::isLowerTriangular ( ) const
virtual

Whether the matrix's graph is locally lower triangular.

Warning
DO NOT CALL THIS METHOD! This method is DEPRECATED and will DISAPPEAR VERY SOON per #2630.
Precondition
! isFillActive(). If fill is active, this method's behavior is undefined.
Note
This is entirely a local property. That means this method may return different results on different processes.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3637 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
bool Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::isUpperTriangular ( ) const
virtual

Whether the matrix's graph is locally upper triangular.

Warning
DO NOT CALL THIS METHOD! This method is DEPRECATED and will DISAPPEAR VERY SOON per #2630.
Precondition
! isFillActive(). If fill is active, this method's behavior is undefined.
Note
This is entirely a local property. That means this method may return different results on different processes.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3646 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
bool Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::isLocallyIndexed ( ) const
virtual

Whether matrix indices are locally indexed.

A RowMatrix may store column indices either as global indices (of type GO), or as local indices (of type LO). In some cases (for example, if the column Map has not been computed), it is not possible to switch from global to local indices without extra work. Furthermore, some operations only work for one or the other case.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3655 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
bool Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::isGloballyIndexed ( ) const
virtual

Whether matrix indices are globally indexed.

A RowMatrix may store column indices either as global indices (of type GO), or as local indices (of type LO). In some cases (for example, if the column Map has not been computed), it is not possible to switch from global to local indices without extra work. Furthermore, some operations only work for one or the other case.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3663 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
bool Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::isFillComplete ( ) const
virtual

Whether fillComplete() has been called.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3671 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
bool Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::supportsRowViews ( ) const
virtual
template<class Scalar, class LO , class GO, class Node >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getGlobalRowCopy ( GO  GlobalRow,
const Teuchos::ArrayView< GO > &  Indices,
const Teuchos::ArrayView< Scalar > &  Values,
size_t &  NumEntries 
) const
virtual

Get a copy of the given global row's entries.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Parameters
GlobalRow[in] Global index of the row.
Indices[out] Global indices of the columns corresponding to values.
Values[out] Matrix values.
NumEntries[out] Number of stored entries on the calling process; length of Indices and Values.

This method throws std::runtime_error if either Indices or Values is not large enough to hold the data associated with row GlobalRow. If GlobalRow does not belong to the calling process, then the method sets NumIndices to Teuchos::OrdinalTraits<size_t>::invalid(), and does not modify Indices or Values.

Definition at line 3688 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO , class GO, class Node >
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getGlobalRowView ( GO  GlobalRow,
Teuchos::ArrayView< const GO > &  indices,
Teuchos::ArrayView< const Scalar > &  values 
) const
virtual

Get a constant, nonpersisting, globally indexed view of the given row of the matrix.

The returned views of the column indices and values are not guaranteed to persist beyond the lifetime of this. Furthermore, some RowMatrix implementations allow changing the values, or the indices and values. Any such changes invalidate the returned views.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Precondition
isGloballyIndexed () && supportsRowViews ()
Postcondition
indices.size () == getNumEntriesInGlobalRow (GlobalRow)
Parameters
GlobalRow[in] Global index of the row.
Indices[out] Global indices of the columns corresponding to values.
Values[out] Matrix values.

If GlobalRow does not belong to this node, then indices is set to null.

Definition at line 3702 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getLocalDiagCopy ( ::Tpetra::Vector< Scalar, LO, GO, Node > &  diag) const
virtual

Get a copy of the diagonal entries, distributed by the row Map.

On input, the Vector's Map must be the same as the row Map of the matrix. (That is, this->getRowMap ()->isSameAs (* (diag.getMap ())) == true.)

On return, the entries of diag are filled with the diagonal entries of the matrix stored on this process. Note that if the row Map is overlapping, multiple processes may own the same diagonal element. You may combine these overlapping diagonal elements by doing an Export from the row Map Vector to a range Map Vector.

Definition at line 3727 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::leftScale ( const ::Tpetra::Vector< Scalar, LO, GO, Node > &  x)
virtual

Scale the RowMatrix on the left with the given Vector x.

On return, for all entries i,j in the matrix, $A(i,j) = x(i)*A(i,j)$.

Definition at line 3785 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar, class LO, class GO, class Node>
void Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::rightScale ( const ::Tpetra::Vector< Scalar, LO, GO, Node > &  x)
virtual

Scale the RowMatrix on the right with the given Vector x.

On return, for all entries i,j in the matrix, $A(i,j) = x(j)*A(i,j)$.

Definition at line 3796 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
typename::Tpetra::RowMatrix< Scalar, LO, GO, Node >::mag_type Tpetra::Experimental::Classes::BlockCrsMatrix< Scalar, LO, GO, Node >::getFrobeniusNorm ( ) const
virtual

The Frobenius norm of the matrix.

This method computes and returns the Frobenius norm of the matrix. The Frobenius norm $\|A\|_F$ for the matrix $A$ is defined as $\|A\|_F = \sqrt{ \sum_{i,j} |A(i,j)|^2 }$. It has the same value as the Euclidean norm of a vector made by stacking the columns of $A$.

Implements Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3815 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual size_t Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInGlobalRow ( GlobalOrdinal  globalRow) const
pure virtualinherited

The current number of entries on the calling process in the specified global row.

Note that if the row Map is overlapping, then the calling process might not necessarily store all the entries in the row. Some other process might have the rest of the entries.

Returns
Teuchos::OrdinalTraits<size_t>::invalid() if the specified global row does not belong to this graph, else the number of entries.

Implemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual size_t Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInLocalRow ( LocalOrdinal  localRow) const
pure virtualinherited

The current number of entries on the calling process in the specified local row.

Note that if the row Map is overlapping, then the calling process might not necessarily store all the entries in the row. Some other process might have the rest of the entries.

Returns
Teuchos::OrdinalTraits<size_t>::invalid() if the specified local row is not valid for this graph, else the number of entries.

Implemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowCopy ( GlobalOrdinal  GlobalRow,
const Teuchos::ArrayView< GlobalOrdinal > &  Indices,
const Teuchos::ArrayView< Scalar > &  Values,
size_t &  NumEntries 
) const
pure virtualinherited

Get a copy of the given global row's entries.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Parameters
GlobalRow[in] Global index of the row.
Indices[out] Global indices of the columns corresponding to values.
Values[out] Matrix values.
NumEntries[out] Number of stored entries on the calling process; length of Indices and Values.

This method throws std::runtime_error if either Indices or Values is not large enough to hold the data associated with row GlobalRow. If GlobalRow does not belong to the calling process, then the method sets NumIndices to Teuchos::OrdinalTraits<size_t>::invalid(), and does not modify Indices or Values.

Implemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowCopy ( LocalOrdinal  LocalRow,
const Teuchos::ArrayView< LocalOrdinal > &  Indices,
const Teuchos::ArrayView< Scalar > &  Values,
size_t &  NumEntries 
) const
pure virtualinherited

Get a copy of the given local row's entries.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Parameters
LocalRow[in] Local index of the row.
Indices[out] Local indices of the columns corresponding to values.
Values[out] Matrix values.
NumEntries[out] Number of stored entries on the calling process; length of Indices and Values.

This method throws std::runtime_error if either Indices or Values is not large enough to hold the data associated with row LocalRow. If LocalRow does not belong to the calling process, then the method sets NumIndices to Teuchos::OrdinalTraits<size_t>::invalid(), and does not modify Indices or Values.

Implemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowView ( GlobalOrdinal  GlobalRow,
Teuchos::ArrayView< const GlobalOrdinal > &  indices,
Teuchos::ArrayView< const Scalar > &  values 
) const
pure virtualinherited

Get a constant, nonpersisting, globally indexed view of the given row of the matrix.

The returned views of the column indices and values are not guaranteed to persist beyond the lifetime of this. Furthermore, some RowMatrix implementations allow changing the values, or the indices and values. Any such changes invalidate the returned views.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Precondition
isGloballyIndexed () && supportsRowViews ()
Postcondition
indices.size () == getNumEntriesInGlobalRow (GlobalRow)
Parameters
GlobalRow[in] Global index of the row.
Indices[out] Global indices of the columns corresponding to values.
Values[out] Matrix values.

If GlobalRow does not belong to this node, then indices is set to null.

Implemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowView ( LocalOrdinal  LocalRow,
Teuchos::ArrayView< const LocalOrdinal > &  indices,
Teuchos::ArrayView< const Scalar > &  values 
) const
pure virtualinherited

Get a constant, nonpersisting, locally indexed view of the given row of the matrix.

The returned views of the column indices and values are not guaranteed to persist beyond the lifetime of this. Furthermore, some RowMatrix implementations allow changing the values, or the indices and values. Any such changes invalidate the returned views.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Precondition
isLocallyIndexed () && supportsRowViews ()
Postcondition
indices.size () == getNumEntriesInGlobalRow (LocalRow)
Parameters
LocalRow[in] Local index of the row.
Indices[out] Local indices of the columns corresponding to values.
Values[out] Matrix values.

If LocalRow does not belong to this node, then indices is set to null.

Implemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
LocalOrdinal Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowViewRaw ( const LocalOrdinal  lclRow,
LocalOrdinal &  numEnt,
const LocalOrdinal *&  lclColInds,
const Scalar *&  vals 
) const
virtualinherited

Get a constant, nonpersisting, locally indexed view of the given row of the matrix, using "raw" pointers instead of Teuchos::ArrayView.

The returned views of the column indices and values are not guaranteed to persist beyond the lifetime of this. Furthermore, some RowMatrix implementations allow changing the values, or the indices and values. Any such changes invalidate the returned views.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Precondition
isLocallyIndexed () && supportsRowViews ()
Postcondition
numEnt == getNumEntriesInGlobalRow (LocalRow)
Parameters
lclRow[in] Local index of the row.
numEnt[out] Number of entries in the row that are stored on the calling process.
lclColInds[out] Local indices of the columns corresponding to values.
vals[out] Matrix values.
Returns
Error code; zero on no error.

Reimplemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 598 of file Tpetra_RowMatrix_def.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagCopy ( Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  diag) const
pure virtualinherited

Get a copy of the diagonal entries, distributed by the row Map.

On input, the Vector's Map must be the same as the row Map of the matrix. (That is, this->getRowMap ()->isSameAs (* (diag.getMap ())) == true.)

On return, the entries of diag are filled with the diagonal entries of the matrix stored on this process. Note that if the row Map is overlapping, multiple processes may own the same diagonal element. You may combine these overlapping diagonal elements by doing an Export from the row Map Vector to a range Map Vector.

Implemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::leftScale ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  x)
pure virtualinherited

Scale the matrix on the left with the given Vector.

On return, for all entries i,j in the matrix, $A(i,j) = x(i)*A(i,j)$.

Implemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::rightScale ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  x)
pure virtualinherited

Scale the matrix on the right with the given Vector.

On return, for all entries i,j in the matrix, $A(i,j) = x(j)*A(i,j)$.

Implemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::add ( const Scalar &  alpha,
const RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Scalar &  beta,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  domainMap = Teuchos::null,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
) const
virtualinherited

Return a new RowMatrix which is the result of beta*this + alpha*A.

The new RowMatrix is actually a CrsMatrix (which see). Note that RowMatrix is a read-only interface (not counting the left and right scale methods), so it is impossible to implement an in-place add using just that interface.

For brevity, call this matrix B, and the result matrix C. C's row Map will be identical to B's row Map. It is correct, though less efficient, for A and B not to have the same row Maps. We could make C's row Map the union of the two row Maps in that case. However, we don't want row Maps to grow for a repeated sequence of additions with matrices with different row Maps. Furthermore, the fact that the user called this method on B, rather than on A, suggests a preference for using B's distribution. The most reasonable thing to do, then, is to use B's row Map for C.

A and B must have identical or congruent communicators. This method must be called as a collective over B's communicator.

The parameters are optional and may be null. Here are the parameters that this function accepts:

  • "Call fillComplete" (bool): If true, call fillComplete on the result matrix C. This is true by default.
  • "Constructor parameters" (sublist): If provided, give these parameters to C's constructor.
  • "fillComplete parameters" (sublist): If provided, and if "Call fillComplete" is true, then give these parameters to C's fillComplete call.

It is not strictly necessary that a RowMatrix always have a domain and range Map. For example, a CrsMatrix does not have a domain and range Map until after its first fillComplete call. Neither A nor B need to have a domain and range Map in order to call add(). If at least one of them has a domain and range Map, you need not supply a domain and range Map to this method. If you ask this method to call fillComplete on C (it does by default), it will supply the any missing domain or range Maps from either B's or A's (in that order) domain and range Maps. If neither A nor B have a domain and range Map, and if you ask this method to call fillComplete, then you must supply both a domain Map and a range Map to this method.

This method comes with a default implementation, since the RowMatrix interface suffices for implementing it. Subclasses (like CrsMatrix) may override this implementation, for example to improve its performance, given additional knowledge about the subclass. Subclass implementations may need to do a dynamic cast on A in order to know its type.

Reimplemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 58 of file Tpetra_RowMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::Classes::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::pack ( const Teuchos::ArrayView< const LocalOrdinal > &  exportLIDs,
Teuchos::Array< char > &  exports,
const Teuchos::ArrayView< size_t > &  numPacketsPerLID,
size_t &  constantNumPackets,
Distributor distor 
) const
virtualinherited

Pack this object's data for an Import or Export.

Warning
To be called only by the packAndPrepare method of appropriate classes of DistObject.

Subclasses may override this method to speed up or otherwise improve the implementation by exploiting more specific details of the subclass.

Implements Tpetra::Classes::Packable< char, LocalOrdinal >.

Definition at line 298 of file Tpetra_RowMatrix_def.hpp.

template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >::apply ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  X,
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  Y,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
Scalar  alpha = Teuchos::ScalarTraits< Scalar >::one(),
Scalar  beta = Teuchos::ScalarTraits< Scalar >::zero() 
) const
pure virtualinherited

Computes the operator-multivector application.

Loosely, performs $Y = \alpha \cdot A^{\textrm{mode}} \cdot X + \beta \cdot Y$. However, the details of operation vary according to the values of alpha and beta. Specifically

  • if beta == 0, apply() must overwrite Y, so that any values in Y (including NaNs) are ignored.
  • if alpha == 0, apply() may short-circuit the operator, so that any values in X (including NaNs) are ignored.

Implemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Packet , class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::doImport ( const SrcDistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
CombineMode  CM 
)
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.

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.

Definition at line 252 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::doImport ( const SrcDistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter,
CombineMode  CM 
)
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.

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.

Definition at line 336 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::doExport ( const SrcDistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter,
CombineMode  CM 
)
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.

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.

Definition at line 294 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::doExport ( const SrcDistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
CombineMode  CM 
)
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.

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.

Definition at line 378 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, 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.

Definition at line 420 of file Tpetra_DistObject_def.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual Teuchos::RCP<const map_type> Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, 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 510 of file Tpetra_DistObject_decl.hpp.

template<class Packet , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, 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.

Definition at line 1609 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, 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::Classes::MultiVector< Scalar, LO, GO, Node >.

Definition at line 214 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, 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::Classes::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::Classes::MultiVector< Scalar, LO, GO, Node >, and Tpetra::Details::CooMatrix< SC, LO, GO, NT >.

Definition at line 427 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::doTransfer ( const SrcDistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &  src,
const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &  transfer,
const char  modeString[],
const ReverseOption  revOp,
const CombineMode  CM 
)
protectedvirtualinherited

Redistribute data across memory images.

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.

Definition at line 434 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, 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 doTransferOld() and doTransferNew(). This needs to be protected, but that doesn't mean users should call this method.

Definition at line 617 of file Tpetra_DistObject_def.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual bool Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::checkSizes ( const SrcDistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &  source)
protectedpure virtualinherited
template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual bool Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::useNewInterface ( )
inlineprotectedvirtualinherited

Whether the subclass implements the "old" or "new" (Kokkos-friendly) interface.

The "old" interface consists of copyAndPermute, packAndPrepare, and unpackAndCombine. The "new" interface consists of copyAndPermuteNew, packAndPrepareNew, and unpackAndCombineNew. We prefer the new interface, because it facilitates thread parallelization using Kokkos data structures.

At some point, we will remove the old interface, and rename the "new" interface (by removing "New" from the methods' names), so that it becomes the only interface.

Reimplemented in Tpetra::Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::Classes::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::Classes::MultiVector< Scalar, LO, GO, Node >, and Tpetra::Details::CooMatrix< SC, LO, GO, NT >.

Definition at line 761 of file Tpetra_DistObject_decl.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::copyAndPermute ( const SrcDistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &  source,
size_t  numSameIDs,
const Teuchos::ArrayView< const local_ordinal_type > &  permuteToLIDs,
const Teuchos::ArrayView< const local_ordinal_type > &  permuteFromLIDs 
)
inlineprotectedvirtualinherited

Perform copies and permutations that are local to this process.

Parameters
source[in] On entry, the source object, from which we are distributing. We distribute to the destination object, which is *this object.
numSameIDs[in] The umber of elements that are the same on the source and destination (this) objects. These elements are owned by the same process in both the source and destination objects. No permutation occurs.
numPermuteIDs[in] The number of elements that are locally permuted between the source and destination objects.
permuteToLIDs[in] List of the elements that are permuted. They are listed by their LID in the destination object.
permuteFromLIDs[in] List of the elements that are permuted. They are listed by their LID in the source object.

Reimplemented in Tpetra::Classes::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Classes::CrsGraph< LO, GO, node_type >.

Definition at line 781 of file Tpetra_DistObject_decl.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::packAndPrepare ( const SrcDistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Teuchos::ArrayView< const local_ordinal_type > &  exportLIDs,
Teuchos::Array< packet_type > &  exports,
const Teuchos::ArrayView< size_t > &  numPacketsPerLID,
size_t &  constantNumPackets,
Distributor distor 
)
inlineprotectedvirtualinherited

Perform any packing or preparation required for communication.

Parameters
source[in] Source object for the redistribution.
exportLIDs[in] List of the entries (as local IDs in the source object) we will be sending to other images.
exports[out] On exit, the buffer for data to send.
numPacketsPerLID[out] On exit, the implementation of this method must do one of two things: set numPacketsPerLID[i] to contain the number of packets to be exported for exportLIDs[i] and set constantNumPackets to zero, or set constantNumPackets to a nonzero value. If the latter, the implementation need not fill numPacketsPerLID.
constantNumPackets[out] On exit, 0 if numPacketsPerLID has variable contents (different size for each LID). If nonzero, then it is expected that the number of packets per LID is constant, and that constantNumPackets is that value.
distor[in] The Distributor object we are using.

Definition at line 816 of file Tpetra_DistObject_decl.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::unpackAndCombine ( const Teuchos::ArrayView< const local_ordinal_type > &  importLIDs,
const Teuchos::ArrayView< const packet_type > &  imports,
const Teuchos::ArrayView< size_t > &  numPacketsPerLID,
size_t  constantNumPackets,
Distributor distor,
CombineMode  CM 
)
inlineprotectedvirtualinherited

Perform any unpacking and combining after communication (old version that uses Teuchos memory management classes to hold data).

Parameters
importLIDs[in] List of the entries (as LIDs in the destination object) we received from other images.
imports[in] Buffer containing data we received.
numPacketsPerLID[in] If constantNumPackets is zero, then numPacketsPerLID[i] contains the number of packets imported for importLIDs[i].
constantNumPackets[in] If nonzero, then numPacketsPerLID is constant (same value in all entries) and constantNumPackets is that value. If zero, then numPacketsPerLID[i] is the number of packets imported for importLIDs[i].
distor[in] The Distributor object we are using.
CM[in] The combine mode to use when combining the imported entries with existing entries.

Definition at line 857 of file Tpetra_DistObject_decl.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
virtual void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::unpackAndCombineNew ( const Kokkos::DualView< const local_ordinal_type *, device_type > &  importLIDs,
const Kokkos::DualView< const packet_type *, buffer_device_type > &  imports,
const Kokkos::DualView< const size_t *, buffer_device_type > &  numPacketsPerLID,
const size_t  constantNumPackets,
Distributor distor,
const CombineMode  CM 
)
inlineprotectedvirtualinherited

Perform any unpacking and combining after communication (new version that uses Kokkos data structures to hold data).

The imports input argument controls whether this method should unpack on host or unpack on device.

Parameters
importLIDs[in] List of the entries (as LIDs in the destination object) we received from other images.
imports[in] Buffer containing data we received.
numPacketsPerLID[in] If constantNumPackets is zero, then numPacketsPerLID[i] contains the number of packets imported for importLIDs[i].
constantNumPackets[in] If nonzero, then numPacketsPerLID is constant (same value in all entries) and constantNumPackets is that value. If zero, then numPacketsPerLID[i] is the number of packets imported for importLIDs[i].
distor[in] The Distributor object we are using.
CM[in] The combine mode to use when combining the imported entries with existing entries.

Definition at line 891 of file Tpetra_DistObject_decl.hpp.

template<class Packet , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::createViews ( ) const
protectedvirtualinherited

Hook for creating a const view.

doTransfer() calls this on the source object. By default, it does nothing, but the source object can use this as a hint to fetch data from a compute buffer on an off-CPU device (such as a GPU) into host memory.

Definition at line 1624 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::createViewsNonConst ( KokkosClassic::ReadWriteOption  rwo)
protectedvirtualinherited

Hook for creating a nonconst view.

doTransfer() calls this on the destination (*this) object. By default, it does nothing, but the destination object can use this as a hint to fetch data from a compute buffer on an off-CPU device (such as a GPU) into host memory.

Parameters
rwo[in] Whether to create a write-only or a read-and-write view. For Kokkos Node types where compute buffers live in a separate memory space (e.g., in the device memory of a discrete accelerator like a GPU), a write-only view only requires copying from host memory to the compute buffer, whereas a read-and-write view requires copying both ways (once to read, from the compute buffer to host memory, and once to write, back to the compute buffer).

Definition at line 1630 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::releaseViews ( ) const
protectedvirtualinherited

Hook for releasing views.

Note
This is no longer called (and is therefore no longer needed) for subclasses for which useNewInterface() returns true.

doTransfer() calls this on both the source and destination objects, once it no longer needs to access that object's data. By default, this method does nothing. Implementations may use this as a hint to free host memory which is a view of a compute buffer, once the host memory view is no longer needed. Some implementations may prefer to mirror compute buffers in host memory; for these implementations, releaseViews() may do nothing.

Definition at line 1636 of file Tpetra_DistObject_def.hpp.

template<class Packet , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::reallocImportsIfNeeded ( const size_t  newSize,
const bool  debug = false 
)
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_.
debug[in] Whether to print (copious) debug output to stderr.
Returns
Whether we actually reallocated.

We don't need a "reallocExportsIfNeeded" method, because exports_ always gets passed into packAndPrepareNew() by nonconst reference. Thus, that method can resize the DualView without needing to call other DistObject methods.

Definition at line 592 of file Tpetra_DistObject_def.hpp.

Member Data Documentation

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
Teuchos::RCP<const map_type> Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node >::map_
protectedinherited

The Map over which this object is distributed.

Definition at line 942 of file Tpetra_DistObject_decl.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
Kokkos::DualView<packet_type*, buffer_device_type> Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, 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 951 of file Tpetra_DistObject_decl.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
Kokkos::DualView<size_t*, buffer_device_type> Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, 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 983 of file Tpetra_DistObject_decl.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
Kokkos::DualView<packet_type*, buffer_device_type> Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, 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 990 of file Tpetra_DistObject_decl.hpp.

template<class Packet, class LocalOrdinal = ::Tpetra::Details::DefaultTypes::local_ordinal_type, class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type, class Node = ::Tpetra::Details::DefaultTypes::node_type>
Kokkos::DualView<size_t*, buffer_device_type> Tpetra::Classes::DistObject< Packet, LocalOrdinal, GlobalOrdinal, 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 1005 of file Tpetra_DistObject_decl.hpp.


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