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

Public Member Functions | |
| 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. | |
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. | |
Operator Methods | |
| const Teuchos::RCP< const Map < LocalOrdinal, GlobalOrdinal, Node > > & | getDomainMap () const |
| Returns the Map associated with the domain of this operator. | |
| const Teuchos::RCP< const Map < LocalOrdinal, GlobalOrdinal, Node > > & | getRangeMap () const |
| Returns the Map associated with the range of this operator, which must be compatible with Y.getMap(). | |
| 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< GlobalOrdinal, Scalar > &blockEntry) |
| Copy the contents of the input block-entry into the matrix. | |
| void | setLocalBlockEntry (LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, const Teuchos::SerialDenseMatrix< LocalOrdinal, Scalar > &blockEntry) |
| Copy the contents of the input block-entry into the matrix. | |
| void | sumIntoGlobalBlockEntry (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, const Teuchos::SerialDenseMatrix< GlobalOrdinal, Scalar > &blockEntry) |
| Add the contents of the input block-entry into the matrix. | |
| void | sumIntoLocalBlockEntry (LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, const Teuchos::SerialDenseMatrix< LocalOrdinal, 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. | |
Extraction Methods | |
| 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. | |
Overridden from Teuchos::Describable | |
| std::string | description () const |
| 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. | |
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 81 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 58 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() != true, 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 82 of file Tpetra_VbrMatrix_def.hpp.
| Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::~VbrMatrix | ( | ) | [virtual] |
Destructor.
Definition at line 119 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 151 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 166 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 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 126 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 Map associated with the range of this operator, which must be compatible with Y.getMap().
Note that this is a point-entry map, not a block-map.
Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 134 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 61 of file Tpetra_CrsMatrixMultiplyOp_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 242 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.
Reimplemented from Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 265 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 273 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 289 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 297 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 305 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 281 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 313 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 142 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 559 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::setGlobalBlockEntry | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| const Teuchos::SerialDenseMatrix< GlobalOrdinal, 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.
This method may be called any time (before or after fillComplete()).
Definition at line 606 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::setLocalBlockEntry | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| const Teuchos::SerialDenseMatrix< LocalOrdinal, 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 622 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::sumIntoGlobalBlockEntry | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| const Teuchos::SerialDenseMatrix< GlobalOrdinal, 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.
This method may be called any time (before or after fillComplete()).
Definition at line 638 of file Tpetra_VbrMatrix_def.hpp.
| void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::sumIntoLocalBlockEntry | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| const Teuchos::SerialDenseMatrix< LocalOrdinal, 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 654 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.
This method may be called any time (before or after fillComplete()).
Definition at line 670 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 687 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.
This method may be called any time (before or after fillComplete()).
Definition at line 704 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 721 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.
Definition at line 894 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 912 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 386 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 321 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 507 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 434 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 481 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 |
Print the object with some verbosity level to a FancyOStream object.
Definition at line 939 of file Tpetra_VbrMatrix_def.hpp.
1.7.4