#include <Xpetra_TpetraVbrMatrix.hpp>
Public Member Functions | |
| RCP< const Tpetra::VbrMatrix < Scalar, LocalOrdinal, GlobalOrdinal, Node > > | getTpetra_VbrMatrix () const |
Private Attributes | |
| const RCP< const Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > | mtx_ |
Constructor/Destructor Methods | |
| TpetraVbrMatrix (const Teuchos::RCP< const Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mtx) | |
| virtual | ~TpetraVbrMatrix () |
| Destructor. | |
Matrix 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. | |
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. | |
Definition at line 67 of file Xpetra_TpetraVbrMatrix.hpp.
| Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::TpetraVbrMatrix | ( | const Teuchos::RCP< const Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > & | mtx | ) | [inline] |
Definition at line 74 of file Xpetra_TpetraVbrMatrix.hpp.
| virtual Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::~TpetraVbrMatrix | ( | ) | [virtual] |
Destructor.
| const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getDomainMap | ( | ) | const [inline] |
Returns the Map associated with the domain of this operator.
Note that this is a point-entry map, not a block-map.
Definition at line 87 of file Xpetra_TpetraVbrMatrix.hpp.
| const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRangeMap | ( | ) | const [inline] |
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.
Definition at line 92 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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 [inline] |
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. Definition at line 100 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::applyInverse | ( | const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | Y, |
| MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | X, | ||
| Teuchos::ETransp | trans | ||
| ) | const [inline] |
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 110 of file Xpetra_TpetraVbrMatrix.hpp.
| bool Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::hasTransposeApply | ( | ) | const [inline] |
Indicates whether this operator supports applying the adjoint operator.
Definition at line 115 of file Xpetra_TpetraVbrMatrix.hpp.
| const Teuchos::RCP<const BlockMap<LocalOrdinal,GlobalOrdinal,Node> >& Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getBlockRowMap | ( | ) | const [inline] |
Returns the block-row map.
Definition at line 123 of file Xpetra_TpetraVbrMatrix.hpp.
| const Teuchos::RCP<const BlockMap<LocalOrdinal,GlobalOrdinal,Node> >& Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getBlockColMap | ( | ) | const [inline] |
Returns the block-column map.
Definition at line 126 of file Xpetra_TpetraVbrMatrix.hpp.
| const Teuchos::RCP<const BlockMap<LocalOrdinal,GlobalOrdinal,Node> >& Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getBlockDomainMap | ( | ) | const [inline] |
Returns the block-domain map.
Definition at line 129 of file Xpetra_TpetraVbrMatrix.hpp.
| const Teuchos::RCP<const BlockMap<LocalOrdinal,GlobalOrdinal,Node> >& Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getBlockRangeMap | ( | ) | const [inline] |
Returns the block-range map.
Definition at line 132 of file Xpetra_TpetraVbrMatrix.hpp.
| const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getPointRowMap | ( | ) | const [inline] |
Returns the point-row map.
Definition at line 135 of file Xpetra_TpetraVbrMatrix.hpp.
| const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getPointColMap | ( | ) | const [inline] |
Returns the point-column map.
Definition at line 138 of file Xpetra_TpetraVbrMatrix.hpp.
| bool Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isFillComplete | ( | ) | const [inline] |
Return true if fillComplete has been called, false otherwise.
Definition at line 141 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::putScalar | ( | Scalar | s | ) | [inline] |
Set the specified scalar throughout the matrix.
This method may be called any time (before or after fillComplete()).
Definition at line 151 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setGlobalBlockEntry | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| const Teuchos::SerialDenseMatrix< GlobalOrdinal, Scalar > & | blockEntry | ||
| ) | [inline] |
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 163 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setLocalBlockEntry | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| const Teuchos::SerialDenseMatrix< LocalOrdinal, Scalar > & | blockEntry | ||
| ) | [inline] |
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 173 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoGlobalBlockEntry | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| const Teuchos::SerialDenseMatrix< GlobalOrdinal, Scalar > & | blockEntry | ||
| ) | [inline] |
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 185 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoLocalBlockEntry | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| const Teuchos::SerialDenseMatrix< LocalOrdinal, Scalar > & | blockEntry | ||
| ) | [inline] |
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 195 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setGlobalBlockEntry | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| LocalOrdinal | blkRowSize, | ||
| LocalOrdinal | blkColSize, | ||
| LocalOrdinal | LDA, | ||
| const Teuchos::ArrayView< const Scalar > & | blockEntry | ||
| ) | [inline] |
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 207 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setLocalBlockEntry | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| LocalOrdinal | blkRowSize, | ||
| LocalOrdinal | blkColSize, | ||
| LocalOrdinal | LDA, | ||
| const Teuchos::ArrayView< const Scalar > & | blockEntry | ||
| ) | [inline] |
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 217 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoGlobalBlockEntry | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| LocalOrdinal | blkRowSize, | ||
| LocalOrdinal | blkColSize, | ||
| LocalOrdinal | LDA, | ||
| const Teuchos::ArrayView< const Scalar > & | blockEntry | ||
| ) | [inline] |
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 229 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoLocalBlockEntry | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| LocalOrdinal | blkRowSize, | ||
| LocalOrdinal | blkColSize, | ||
| LocalOrdinal | LDA, | ||
| const Teuchos::ArrayView< const Scalar > & | blockEntry | ||
| ) | [inline] |
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 239 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fillComplete | ( | const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & | blockDomainMap, |
| const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & | blockRangeMap | ||
| ) | [inline] |
Transition the matrix to the packed, optimized-storage state.
This method also sets the domain and range maps.
Definition at line 250 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fillComplete | ( | ) | [inline] |
Transition the matrix to the packed, optimized-storage state.
This method internally calls fillComplete(getBlockRowMap(),getBlockRowMap()).
Definition at line 256 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalBlockEntryView | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| LocalOrdinal & | numPtRows, | ||
| LocalOrdinal & | numPtCols, | ||
| Teuchos::ArrayRCP< const Scalar > & | blockEntry | ||
| ) | const [inline] |
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 271 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalBlockEntryViewNonConst | ( | GlobalOrdinal | globalBlockRow, |
| GlobalOrdinal | globalBlockCol, | ||
| LocalOrdinal & | numPtRows, | ||
| LocalOrdinal & | numPtCols, | ||
| Teuchos::ArrayRCP< Scalar > & | blockEntry | ||
| ) | [inline] |
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 288 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalBlockEntryView | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| LocalOrdinal & | numPtRows, | ||
| LocalOrdinal & | numPtCols, | ||
| Teuchos::ArrayRCP< const Scalar > & | blockEntry | ||
| ) | const [inline] |
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 305 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalBlockEntryViewNonConst | ( | LocalOrdinal | localBlockRow, |
| LocalOrdinal | localBlockCol, | ||
| LocalOrdinal & | numPtRows, | ||
| LocalOrdinal & | numPtCols, | ||
| Teuchos::ArrayRCP< Scalar > & | blockEntry | ||
| ) | [inline] |
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 328 of file Xpetra_TpetraVbrMatrix.hpp.
| std::string Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::description | ( | ) | const [inline] |
Definition at line 338 of file Xpetra_TpetraVbrMatrix.hpp.
| void Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::describe | ( | Teuchos::FancyOStream & | out, |
| const Teuchos::EVerbosityLevel | verbLevel = Teuchos::Describable::verbLevel_default |
||
| ) | const [inline] |
Print the object with some verbosity level to a FancyOStream object.
Definition at line 342 of file Xpetra_TpetraVbrMatrix.hpp.
| RCP< const Tpetra::VbrMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getTpetra_VbrMatrix | ( | ) | const [inline] |
Definition at line 345 of file Xpetra_TpetraVbrMatrix.hpp.
const RCP< const Tpetra::VbrMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Xpetra::TpetraVbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mtx_ [private] |
Definition at line 349 of file Xpetra_TpetraVbrMatrix.hpp.
1.7.6.1