|
Tpetra Matrix/Vector Services
Version of the Day
|
VbrMatrix: Variable block row matrix. More...
#include <Tpetra_VbrMatrix_decl.hpp>

Public Types | |
| typedef Scalar | scalar_type |
| The type of the entries of the input and output multivectors. | |
| typedef LocalOrdinal | local_ordinal_type |
| The local index type. | |
| typedef GlobalOrdinal | global_ordinal_type |
| The global index type. | |
| typedef Node | node_type |
| The Kokkos Node type. | |
Public Member Functions | |
Constructor/Destructor Methods | |
| VbrMatrix (const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > &blkRowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile) | |
| Constructor specifying the row-map and the max number of (block) non-zeros for all rows. | |
| VbrMatrix (const Teuchos::RCP< const BlockCrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &blkGraph) | |
| Constructor specifying a pre-filled block-graph. | |
| virtual | ~VbrMatrix () |
| Destructor. | |
Advanced Mathematical operations | |
| template<class DomainScalar , class RangeScalar > | |
| void | multiply (const MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp trans, RangeScalar alpha, RangeScalar beta) const |
| Multiply this matrix by a MultiVector. | |
| template<class DomainScalar , class RangeScalar > | |
| void | solve (const MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &Y, MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > &X, Teuchos::ETransp trans) const |
| Triangular Solve -- Matrix must be triangular. | |
Operator Methods | |
| const Teuchos::RCP< const Map < LocalOrdinal, GlobalOrdinal, Node > > & | getDomainMap () const |
| Returns the (point-entry) Map associated with the domain of this operator. | |
| const Teuchos::RCP< const Map < LocalOrdinal, GlobalOrdinal, Node > > & | getRangeMap () const |
| Returns the (point-entry) Map associated with the range of this operator. | |
| void | apply (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp trans=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const |
| Computes the operator-multivector application. | |
| void | applyInverse (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Teuchos::ETransp trans) const |
| Triangular Solve -- Matrix must be triangular. | |
| bool | hasTransposeApply () const |
| Indicates whether this operator supports applying the adjoint operator. | |
Attribute Query Methods | |
| const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & | getBlockRowMap () const |
| Returns the block-row map. | |
| const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & | getBlockColMap () const |
| Returns the block-column map. | |
| const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & | getBlockDomainMap () const |
| Returns the block-domain map. | |
| const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & | getBlockRangeMap () const |
| Returns the block-range map. | |
| const Teuchos::RCP< const Map < LocalOrdinal, GlobalOrdinal, Node > > & | getPointRowMap () const |
| Returns the point-row map. | |
| const Teuchos::RCP< const Map < LocalOrdinal, GlobalOrdinal, Node > > & | getPointColMap () const |
| Returns the point-column map. | |
| bool | isFillComplete () const |
| Return true if fillComplete has been called, false otherwise. | |
Insertion Methods | |
| void | putScalar (Scalar s) |
| Set the specified scalar throughout the matrix. | |
| void | setGlobalBlockEntry (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, const Teuchos::SerialDenseMatrix< int, Scalar > &blockEntry) |
| Copy the contents of the input block-entry into the matrix. | |
| void | setLocalBlockEntry (LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, const Teuchos::SerialDenseMatrix< int, Scalar > &blockEntry) |
| Copy the contents of the input block-entry into the matrix. | |
| void | sumIntoGlobalBlockEntry (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, const Teuchos::SerialDenseMatrix< int, Scalar > &blockEntry) |
| Add the contents of the input block-entry into the matrix. | |
| void | sumIntoLocalBlockEntry (LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, const Teuchos::SerialDenseMatrix< int, Scalar > &blockEntry) |
| Add the contents of the input block-entry into the matrix. | |
| void | setGlobalBlockEntry (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, LocalOrdinal blkRowSize, LocalOrdinal blkColSize, LocalOrdinal LDA, const Teuchos::ArrayView< const Scalar > &blockEntry) |
| Copy the contents of the input block-entry into the matrix. | |
| void | setLocalBlockEntry (LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, LocalOrdinal blkRowSize, LocalOrdinal blkColSize, LocalOrdinal LDA, const Teuchos::ArrayView< const Scalar > &blockEntry) |
| Copy the contents of the input block-entry into the matrix. | |
| void | sumIntoGlobalBlockEntry (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, LocalOrdinal blkRowSize, LocalOrdinal blkColSize, LocalOrdinal LDA, const Teuchos::ArrayView< const Scalar > &blockEntry) |
| Add the contents of the input block-entry into the matrix. | |
| void | sumIntoLocalBlockEntry (LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, LocalOrdinal blkRowSize, LocalOrdinal blkColSize, LocalOrdinal LDA, const Teuchos::ArrayView< const Scalar > &blockEntry) |
| Add the contents of the input block-entry into the matrix. | |
Transformational Methods | |
| void | fillComplete (const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > &blockDomainMap, const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > &blockRangeMap) |
| Transition the matrix to the packed, optimized-storage state. | |
| void | fillComplete () |
| Transition the matrix to the packed, optimized-storage state. | |
| void | globalAssemble () |
| Communicate non-local contributions to the processors that own those contributions. | |
Extraction Methods | |
| void | getGlobalBlockRowView (GlobalOrdinal globalBlockRow, LocalOrdinal &numPtRows, Teuchos::ArrayView< const GlobalOrdinal > &blockCols, Teuchos::Array< LocalOrdinal > &ptColsPerBlockCol, Teuchos::Array< Teuchos::ArrayRCP< const Scalar > > &blockEntries) const |
| Returns a const read-only view of the block-entries in the specified row. | |
| void | getLocalBlockRowView (LocalOrdinal localBlockRow, LocalOrdinal &numPtRows, Teuchos::ArrayView< const LocalOrdinal > &blockCols, Teuchos::Array< LocalOrdinal > &ptColsPerBlockCol, Teuchos::ArrayRCP< const Scalar > &blockEntries) const |
| Returns a const read-only view of the block-entries in the specified row. | |
| void | getGlobalBlockEntryView (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, LocalOrdinal &numPtRows, LocalOrdinal &numPtCols, Teuchos::ArrayRCP< const Scalar > &blockEntry) const |
| Returns a const read-only view of a block-entry. | |
| void | getGlobalBlockEntryViewNonConst (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, LocalOrdinal &numPtRows, LocalOrdinal &numPtCols, Teuchos::ArrayRCP< Scalar > &blockEntry) |
| Returns a non-const read-write view of a block-entry. | |
| void | getLocalBlockEntryView (LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, LocalOrdinal &numPtRows, LocalOrdinal &numPtCols, Teuchos::ArrayRCP< const Scalar > &blockEntry) const |
| Returns a const read-only view of a block-entry. | |
| void | getLocalBlockEntryViewNonConst (LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, LocalOrdinal &numPtRows, LocalOrdinal &numPtCols, Teuchos::ArrayRCP< Scalar > &blockEntry) |
| Returns a non-const read-write view of a block-entry. | |
| void | getLocalDiagCopy (Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const |
| Return a copy of the (point-entry) diagonal values. | |
|
const Teuchos::RCP< const BlockCrsGraph< LocalOrdinal, GlobalOrdinal, Node > > & | getBlockCrsGraph () |
Overridden from Teuchos::DistObject | |
| bool | checkSizes (const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source) |
| Compare the source and target (this) objects for compatibility. | |
| void | copyAndPermute (const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, size_t numSameIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteToLIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteFromLIDs) |
| Perform copies and permutations that are local to this process. | |
| void | packAndPrepare (const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Teuchos::Array< char > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets, Distributor &distor) |
| Perform any packing or preparation required for communication. | |
| void | unpackAndCombine (const Teuchos::ArrayView< const LocalOrdinal > &importLIDs, const Teuchos::ArrayView< const char > &imports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t constantNumPackets, Distributor &distor, CombineMode CM) |
| Perform any unpacking and combining after communication. | |
Overridden from Teuchos::Describable | |
| std::string | description () const |
| One-line descriptiion of this object. | |
| void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const |
| Print the object with some verbosity level to a FancyOStream object. | |
Public methods for redistributing data | |
| void | doImport (const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM) |
| Import using an Import object ("forward mode"). | |
| void | doImport (const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM) |
| Import using an Export object ("reverse mode"). | |
| void | doExport (const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM) |
| Export using an Export object ("forward mode"). | |
| void | doExport (const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM) |
| Export using an Import object ("reverse mode"). | |
Attribute accessor methods | |
| bool | isDistributed () const |
| Whether this is a globally distributed object. | |
| const Teuchos::RCP< const Map < LocalOrdinal, GlobalOrdinal, Node > > & | getMap () const |
| The Map with which this DistObject was constructed. | |
I/O methods | |
| void | print (std::ostream &os) const |
| Print this object to the given output stream. | |
Protected Member Functions | |
| virtual void | doTransfer (const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, CombineMode CM, size_t numSameIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteToLIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteFromLIDs, const Teuchos::ArrayView< const LocalOrdinal > &remoteLIDs, const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Distributor &distor, ReverseOption revOp) |
| Redistribute data across memory images. | |
| virtual void | createViews () const |
| Hook for creating a const view. | |
| virtual void | createViewsNonConst (Kokkos::ReadWriteOption rwo) |
| Hook for creating a nonconst view. | |
| virtual void | releaseViews () const |
| Hook for releasing views. | |
Protected Attributes | |
| Teuchos::RCP< const Map < LocalOrdinal, GlobalOrdinal, Node > > | map_ |
| The Map over which this object is distributed. | |
VbrMatrix: Variable block row matrix.
The VbrMatrix class has two significant 'states', distinguished by whether or not storage has been optimized (packed) or not.
When the matrix is in the non-optimized-storage state, internal data storage is in a non-contiguous data-structure that allows for convenient insertion of data.
When the matrix is in the optimized-storage state, internal data is stored in contiguous (packed) arrays. When in this state, existing entries may be updated and replaced, but no new entries (indices and/or coefficients) may be inserted. In other words, the sparsity pattern or structure of the matrix may not be changed.
Use of the matrix as an Operator (performing matrix-vector multiplication) is only allowed when it is in the optimized-storage state.
VbrMatrix has two constructors, one which leaves the matrix in the optimized- storage state, and another which leaves the matrix in the non-optimized-storage state.
When the VbrMatrix is constructed in the non-optimized-storage state, (and then filled using methods such as setGlobalBlockEntry etc.), it can then be transformed to the optimized-storage state by calling the method fillComplete().
Once in the optimized-storage state, the VbrMatrix can not be returned to the non-optimized-storage state.
Definition at line 96 of file Tpetra_VbrMatrix_decl.hpp.
| typedef Scalar Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::scalar_type |
The type of the entries of the input and output multivectors.
Reimplemented from Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 98 of file Tpetra_VbrMatrix_decl.hpp.
| typedef LocalOrdinal Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::local_ordinal_type |
The local index type.
Reimplemented from Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 99 of file Tpetra_VbrMatrix_decl.hpp.
| typedef GlobalOrdinal Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::global_ordinal_type |
The global index type.
Reimplemented from Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 100 of file Tpetra_VbrMatrix_decl.hpp.
| typedef Node Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::node_type |
The Kokkos Node type.
Reimplemented from Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 101 of file Tpetra_VbrMatrix_decl.hpp.
| Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::VbrMatrix | ( | const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & | blkRowMap, |
| size_t | maxNumEntriesPerRow, | ||
| ProfileType | pftype = DynamicProfile |
||
| ) |
Constructor specifying the row-map and the max number of (block) non-zeros for all rows.
After this constructor completes, the VbrMatrix is in the non-packed, non-optimized-storage, isFillComplete()==false state. Block-entries (rectangular, dense submatrices) may be inserted using class methods such as setGlobalBlockEntry(...), declared below.
Definition at line 71 of file Tpetra_VbrMatrix_def.hpp.
| Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::VbrMatrix | ( | const Teuchos::RCP< const BlockCrsGraph< LocalOrdinal, GlobalOrdinal, Node > > & | blkGraph | ) |
Constructor specifying a pre-filled block-graph.
Constructing a VbrMatrix with a pre-filled graph means that the matrix will start out in the optimized-storage state, i.e., isFillComplete()==true. The graph provided to this constructor must be already filled. (If blkGraph->isFillComplete() == false, an exception is thrown.)
Entries in the input BlockCrsGraph correspond to block-entries in the VbrMatrix. In other words, the VbrMatrix will have a block-row corresponding to each row in the graph, and a block-entry corresponding to each column- index in the graph.
Definition at line 96 of file Tpetra_VbrMatrix_def.hpp.
| Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::~VbrMatrix | ( | ) | [virtual] |
Destructor.
Definition at line 134 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::multiply | ( | const MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > & | X, |
| MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > & | Y, | ||
| Teuchos::ETransp | trans, | ||
| RangeScalar | alpha, | ||
| RangeScalar | beta | ||
| ) | const |
Multiply this matrix by a MultiVector.
X is required to be post-imported, i.e., described by the column map of the matrix. Y is required to be pre-exported, i.e., described by the row map of the matrix. See also the Operator::apply method which is implemented below.
Definition at line 166 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::solve | ( | const MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > & | Y, |
| MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > & | X, | ||
| Teuchos::ETransp | trans | ||
| ) | const |
Triangular Solve -- Matrix must be triangular.
Find X such that A*X = Y. X is required to be post-imported, i.e., described by the column map of the matrix. Y is required to be pre-exported, i.e., described by the row map of the matrix.
Both X and Y are required to have constant stride.
Note that if the diagonal block-entries are stored, they must be triangular. I.e., the matrix structure must be block-triangular, and any diagonal blocks must be "point"-triangular, meaning that coefficients on the "wrong" side of the point-diagonal must be zero.
Definition at line 181 of file Tpetra_VbrMatrix_def.hpp.
| const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getDomainMap | ( | ) | const [virtual] |
Returns the (point-entry) Map associated with the domain of this operator.
Note that this is a point-entry map, not a block-map.
Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 141 of file Tpetra_VbrMatrix_def.hpp.
| const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getRangeMap | ( | ) | const [virtual] |
Returns the (point-entry) Map associated with the range of this operator.
Note that this is a point-entry map, not a block-map.
Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 149 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::apply | ( | const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | X, |
| MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | Y, | ||
| Teuchos::ETransp | trans = Teuchos::NO_TRANS, |
||
| Scalar | alpha = Teuchos::ScalarTraits<Scalar>::one(), |
||
| Scalar | beta = Teuchos::ScalarTraits<Scalar>::zero() |
||
| ) | const [virtual] |
Computes the operator-multivector application.
Loosely, performs
. However, the details of operation vary according to the values of alpha and beta. Specifically
beta == 0, apply() must overwrite Y, so that any values in Y (including NaNs) are ignored.alpha == 0, apply() may short-circuit the operator, so that any values in X (including NaNs) are ignored. Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 232 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::applyInverse | ( | const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | Y, |
| MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | X, | ||
| Teuchos::ETransp | trans | ||
| ) | const |
Triangular Solve -- Matrix must be triangular.
Find X such that A*X = Y. Both X and Y are required to have constant stride.
Definition at line 257 of file Tpetra_VbrMatrix_def.hpp.
| bool Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::hasTransposeApply | ( | ) | const [virtual] |
Indicates whether this operator supports applying the adjoint operator.
VbrMatrix does support transpose-apply. (This method returns true.)
Reimplemented from Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 280 of file Tpetra_VbrMatrix_def.hpp.
| const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getBlockRowMap | ( | ) | const |
Returns the block-row map.
Definition at line 288 of file Tpetra_VbrMatrix_def.hpp.
| const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getBlockColMap | ( | ) | const |
Returns the block-column map.
Definition at line 304 of file Tpetra_VbrMatrix_def.hpp.
| const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getBlockDomainMap | ( | ) | const |
Returns the block-domain map.
Definition at line 312 of file Tpetra_VbrMatrix_def.hpp.
| const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getBlockRangeMap | ( | ) | const |
Returns the block-range map.
Definition at line 320 of file Tpetra_VbrMatrix_def.hpp.
| const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getPointRowMap | ( | ) | const |
Returns the point-row map.
Definition at line 296 of file Tpetra_VbrMatrix_def.hpp.
| const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getPointColMap | ( | ) | const |
Returns the point-column map.
Definition at line 328 of file Tpetra_VbrMatrix_def.hpp.
| bool Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::isFillComplete | ( | ) | const |
Return true if fillComplete has been called, false otherwise.
Definition at line 157 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::putScalar | ( | Scalar | s | ) |
Set the specified scalar throughout the matrix.
This method may be called any time (before or after fillComplete()).
Definition at line 922 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::setGlobalBlockEntry | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| const Teuchos::SerialDenseMatrix< int, Scalar > & | blockEntry | ||
| ) |
Copy the contents of the input block-entry into the matrix.
This method will create the specified block-entry if it doesn't already exist, but only if fillComplete() has not yet been called.
If the specified block-entry already exists in the matrix, it will be over-written (replaced) by the input block-entry.
Note that if globalBlockRow is not owned by the local processor (as indicated by getBlockRowMap()) then the block-entry is held in temporary storage until globalAssemble() is called (which is called internally by fillComplete()) and then globalAssemble() performans the communication needed to move the data to the owning processor.
This method may be called any time (before or after fillComplete()).
Definition at line 973 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::setLocalBlockEntry | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| const Teuchos::SerialDenseMatrix< int, Scalar > & | blockEntry | ||
| ) |
Copy the contents of the input block-entry into the matrix.
This method will throw an exception if fillComplete() has not yet been called, or if the specified block-entry doesn't already exist in the matrix.
The coefficients of the specified block-entry will be over-written (replaced) by the input block-entry.
Definition at line 989 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::sumIntoGlobalBlockEntry | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| const Teuchos::SerialDenseMatrix< int, Scalar > & | blockEntry | ||
| ) |
Add the contents of the input block-entry into the matrix.
This method will create the specified block-entry if it doesn't already exist, but only if fillComplete() has not yet been called.
If the specified block-entry already exists in the matrix, the contents of the input block-entry will be added to the values that are already present.
Note that if globalBlockRow is not owned by the local processor (as indicated by getBlockRowMap()) then the block-entry is held in temporary storage until globalAssemble() is called (which is called internally by fillComplete()) and then globalAssemble() performans the communication needed to move the data to the owning processor.
This method may be called any time (before or after fillComplete()).
Definition at line 1005 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::sumIntoLocalBlockEntry | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| const Teuchos::SerialDenseMatrix< int, Scalar > & | blockEntry | ||
| ) |
Add the contents of the input block-entry into the matrix.
This method will throw an exception if fillComplete() has not yet been called, or if the specified block-entry doesn't already exist in the matrix.
The contents of the input block-entry will be added to the values that are already present in the matrix.
Definition at line 1021 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::setGlobalBlockEntry | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| LocalOrdinal | blkRowSize, | ||
| LocalOrdinal | blkColSize, | ||
| LocalOrdinal | LDA, | ||
| const Teuchos::ArrayView< const Scalar > & | blockEntry | ||
| ) |
Copy the contents of the input block-entry into the matrix.
This method will create the specified block-entry if it doesn't already exist, but only if fillComplete() has not yet been called.
If the specified block-entry already exists in the matrix, it will be over-written (replaced) by the input block-entry.
Note that if globalBlockRow is not owned by the local processor (as indicated by getBlockRowMap()) then the block-entry is held in temporary storage until globalAssemble() is called (which is called internally by fillComplete()) and then globalAssemble() performans the communication needed to move the data to the owning processor.
This method may be called any time (before or after fillComplete()).
Definition at line 1037 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::setLocalBlockEntry | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| LocalOrdinal | blkRowSize, | ||
| LocalOrdinal | blkColSize, | ||
| LocalOrdinal | LDA, | ||
| const Teuchos::ArrayView< const Scalar > & | blockEntry | ||
| ) |
Copy the contents of the input block-entry into the matrix.
This method will throw an exception if fillComplete() has not yet been called, or if the specified block-entry doesn't already exist in the matrix.
The coefficients for the specified block-entry will be over-written (replaced) by the input block-entry.
Definition at line 1054 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::sumIntoGlobalBlockEntry | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| LocalOrdinal | blkRowSize, | ||
| LocalOrdinal | blkColSize, | ||
| LocalOrdinal | LDA, | ||
| const Teuchos::ArrayView< const Scalar > & | blockEntry | ||
| ) |
Add the contents of the input block-entry into the matrix.
This method will create the specified block-entry if it doesn't already exist, but only if fillComplete() has not yet been called.
If the specified block-entry already exists in the matrix, the contents of the input block-entry will be added to the values that are already present.
Note that if globalBlockRow is not owned by the local processor (as indicated by getBlockRowMap()) then the block-entry is held in temporary storage until globalAssemble() is called (which is called internally by fillComplete()) and then globalAssemble() performans the communication needed to move the data to the owning processor.
This method may be called any time (before or after fillComplete()).
Definition at line 1071 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::sumIntoLocalBlockEntry | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| LocalOrdinal | blkRowSize, | ||
| LocalOrdinal | blkColSize, | ||
| LocalOrdinal | LDA, | ||
| const Teuchos::ArrayView< const Scalar > & | blockEntry | ||
| ) |
Add the contents of the input block-entry into the matrix.
This method will throw an exception if fillComplete() has not yet been called, or if the specified block-entry doesn't already exist in the matrix.
The contents of the input block-entry will be added to the values that are already present in the matrix.
Definition at line 1088 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::fillComplete | ( | const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & | blockDomainMap, |
| const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & | blockRangeMap | ||
| ) |
Transition the matrix to the packed, optimized-storage state.
This method also sets the domain and range maps. This method internally calls globalAssemble().
Definition at line 1260 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::fillComplete | ( | ) |
Transition the matrix to the packed, optimized-storage state.
This method internally calls fillComplete(getBlockRowMap(),getBlockRowMap()).
Definition at line 1280 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::globalAssemble | ( | ) |
Communicate non-local contributions to the processors that own those contributions.
Definition at line 1290 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getGlobalBlockRowView | ( | GlobalOrdinal | globalBlockRow, |
| LocalOrdinal & | numPtRows, | ||
| Teuchos::ArrayView< const GlobalOrdinal > & | blockCols, | ||
| Teuchos::Array< LocalOrdinal > & | ptColsPerBlockCol, | ||
| Teuchos::Array< Teuchos::ArrayRCP< const Scalar > > & | blockEntries | ||
| ) | const |
Returns a const read-only view of the block-entries in the specified row.
Can only be called if isFillComplete()==false
Definition at line 336 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getLocalBlockRowView | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal & | numPtRows, | ||
| Teuchos::ArrayView< const LocalOrdinal > & | blockCols, | ||
| Teuchos::Array< LocalOrdinal > & | ptColsPerBlockCol, | ||
| Teuchos::ArrayRCP< const Scalar > & | blockEntries | ||
| ) | const |
Returns a const read-only view of the block-entries in the specified row.
Can only be called if isFillComplete()==true
Definition at line 362 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getGlobalBlockEntryView | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| LocalOrdinal & | numPtRows, | ||
| LocalOrdinal & | numPtCols, | ||
| Teuchos::ArrayRCP< const Scalar > & | blockEntry | ||
| ) | const |
Returns a const read-only view of a block-entry.
The arguments numPtRows and numPtCols are set to the dimensions of the block- entry on output. The stride (LDA in Blas terminology) is equal to numPtRows.
This method may be called any time (before or after fillComplete()), but will throw an exception if the specified block-entry doesn't already exist.
Definition at line 462 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getGlobalBlockEntryViewNonConst | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| LocalOrdinal & | numPtRows, | ||
| LocalOrdinal & | numPtCols, | ||
| Teuchos::ArrayRCP< Scalar > & | blockEntry | ||
| ) |
Returns a non-const read-write view of a block-entry.
Creates the block-entry if it doesn't already exist, and if:
Important Note: Be very careful managing the lifetime of this view. If fillComplete() has been called, and if you are running on a GPU, this view may be a copy of memory from the GPU, and your changes to the view won't be copied back to the GPU until your ArrayRCP is destroyed or set to Teuchos::null.
Definition at line 396 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getLocalBlockEntryView | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| LocalOrdinal & | numPtRows, | ||
| LocalOrdinal & | numPtCols, | ||
| Teuchos::ArrayRCP< const Scalar > & | blockEntry | ||
| ) | const |
Returns a const read-only view of a block-entry.
The arguments numPtRows and numPtCols are set to the dimensions of the block- entry on output. The stride (LDA in Blas terminology) is equal to numPtRows. Throws an exception if fillComplete() has not yet been called, or if the specified block-entry doesn't exist.
This method may only be called after fillComplete() has been called, and will throw an exception if the specified block-entry doesn't already exist.
Definition at line 870 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getLocalBlockEntryViewNonConst | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| LocalOrdinal & | numPtRows, | ||
| LocalOrdinal & | numPtCols, | ||
| Teuchos::ArrayRCP< Scalar > & | blockEntry | ||
| ) |
Returns a non-const read-write view of a block-entry.
The arguments numPtRows and numPtCols are set to the dimensions of the block- entry on output. The stride (LDA in Blas terminology) is equal to numPtRows. Throws an exception if fillComplete() has not yet been called, or if the specified block-entry doesn't exist.
Important Note: Be very careful managing the lifetime of this view. If fillComplete() has been called, and if you are running on a GPU, this view may be a copy of memory from the GPU, and your changes to the view won't be copied back to the GPU until your ArrayRCP is destroyed or set to Teuchos::null.
This method may only be called after fillComplete() has been called, and will throw an exception if the specified block-entry doesn't already exist.
Definition at line 510 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getLocalDiagCopy | ( | Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | diag | ) | const |
Return a copy of the (point-entry) diagonal values.
Throws an exception if the input-vector's map is not the same as getBlockRowMap()->getPointMap().
Definition at line 557 of file Tpetra_VbrMatrix_def.hpp.
| bool Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::checkSizes | ( | const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > & | source | ) | [virtual] |
Compare the source and target (this) objects for compatibility.
Implements Tpetra::DistObject< char, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 583 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::copyAndPermute | ( | const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > & | source, |
| size_t | numSameIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | permuteToLIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | permuteFromLIDs | ||
| ) | [virtual] |
Perform copies and permutations that are local to this process.
| 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. |
Implements Tpetra::DistObject< char, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 593 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::packAndPrepare | ( | const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > & | source, |
| const Teuchos::ArrayView< const LocalOrdinal > & | exportLIDs, | ||
| Teuchos::Array< char > & | exports, | ||
| const Teuchos::ArrayView< size_t > & | numPacketsPerLID, | ||
| size_t & | constantNumPackets, | ||
| Distributor & | distor | ||
| ) | [virtual] |
Perform any packing or preparation required for communication.
| 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, numPacketsPerLID[i] contains the number of packets to be exported for exportLIDs[i]. If constantNumPackets is nonzero, you should use that instead, and not rely on numPacketsPerLID[i] being filled. |
| constantNumPackets | [out] On exit, 0 if numPacketsPerLID has variable contents (different size for each LID). If nonzero, then it is expected that num-packets-per-LID is constant, and constantNumPackets holds that value. |
| distor | [in] The Distributor object we are using. |
Implements Tpetra::DistObject< char, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 670 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::unpackAndCombine | ( | const Teuchos::ArrayView< const LocalOrdinal > & | importLIDs, |
| const Teuchos::ArrayView< const char > & | imports, | ||
| const Teuchos::ArrayView< size_t > & | numPacketsPerLID, | ||
| size_t | constantNumPackets, | ||
| Distributor & | distor, | ||
| CombineMode | CM | ||
| ) | [virtual] |
Perform any unpacking and combining after communication.
| 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] 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, use numPacketsPerLID[i] instead. |
| distor | [in] The Distributor object we are using. |
| CM | [in] The combine mode to use when combining the imported entries with existing entries. |
Implements Tpetra::DistObject< char, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 818 of file Tpetra_VbrMatrix_def.hpp.
| std::string Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::description | ( | ) | const [virtual] |
One-line descriptiion of this object.
We declare this method virtual so that subclasses of DistObject may override it.
Reimplemented from Tpetra::DistObject< char, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1329 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::describe | ( | Teuchos::FancyOStream & | out, |
| const Teuchos::EVerbosityLevel | verbLevel = Teuchos::Describable::verbLevel_default |
||
| ) | const [virtual] |
Print the object with some verbosity level to a FancyOStream object.
Reimplemented from Tpetra::DistObject< char, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1345 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::doImport | ( | const DistObject< char , LocalOrdinal, GlobalOrdinal, Node > & | source, |
| const Import< LocalOrdinal, GlobalOrdinal, Node > & | importer, | ||
| CombineMode | CM | ||
| ) | [inherited] |
Import using an Import object ("forward mode").
| void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::doImport | ( | const DistObject< char , LocalOrdinal, GlobalOrdinal, Node > & | source, |
| const Export< LocalOrdinal, GlobalOrdinal, Node > & | exporter, | ||
| CombineMode | CM | ||
| ) | [inherited] |
Import using an Export object ("reverse mode").
| void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::doExport | ( | const DistObject< char , LocalOrdinal, GlobalOrdinal, Node > & | dest, |
| const Export< LocalOrdinal, GlobalOrdinal, Node > & | exporter, | ||
| CombineMode | CM | ||
| ) | [inherited] |
Export using an Export object ("forward mode").
| void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::doExport | ( | const DistObject< char , LocalOrdinal, GlobalOrdinal, Node > & | dest, |
| const Import< LocalOrdinal, GlobalOrdinal, Node > & | importer, | ||
| CombineMode | CM | ||
| ) | [inherited] |
Export using an Import object ("reverse mode").
| bool Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::isDistributed | ( | ) | const [inline, 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.
| const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::getMap | ( | ) | const [inline, inherited] |
The Map with which this DistObject was constructed.
Definition at line 172 of file Tpetra_DistObject.hpp.
| void Tpetra::DistObject< char , 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.
| virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::doTransfer | ( | const DistObject< char , LocalOrdinal, GlobalOrdinal, Node > & | source, |
| CombineMode | CM, | ||
| size_t | numSameIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | permuteToLIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | permuteFromLIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | remoteLIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | exportLIDs, | ||
| Distributor & | distor, | ||
| ReverseOption | revOp | ||
| ) | [protected, virtual, inherited] |
Redistribute data across memory images.
| source | [in] The source object, to redistribute into the destination object, which is *this object. |
| CM | [in] The combine mode that describes how to combine values that map to the same global ID on the same process. |
| permuteToLIDs | [in] See copyAndPermute(). |
| permuteFromLIDs | [in] See copyAndPermute(). |
| remoteLIDs | [in] List of entries (as local IDs) in the destination object to receive from other processes. |
| exportLIDs | [in] See packAndPrepare(). |
| distor | [in/out] The Distributor object that knows how to redistribute data. |
| revOp | [in] Whether to do a forward or reverse mode redistribution. |
| virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::createViews | ( | ) | const [inline, protected, virtual, inherited] |
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 358 of file Tpetra_DistObject.hpp.
| virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::createViewsNonConst | ( | Kokkos::ReadWriteOption | rwo | ) | [inline, protected, virtual, inherited] |
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.
| 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 375 of file Tpetra_DistObject.hpp.
| virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::releaseViews | ( | ) | const [inline, protected, virtual, inherited] |
Hook for releasing views.
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 387 of file Tpetra_DistObject.hpp.
Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::map_ [protected, inherited] |
The Map over which this object is distributed.
Definition at line 390 of file Tpetra_DistObject.hpp.
1.7.6.1