Xpetra
Version of the Day
|
Xpetra-specific matrix class. More...
#include <Xpetra_Matrix_decl.hpp>
Public Types | |
typedef Scalar | scalar_type |
typedef LocalOrdinal | local_ordinal_type |
typedef GlobalOrdinal | global_ordinal_type |
typedef Node | node_type |
typedef CrsMatrix::local_matrix_type | local_matrix_type |
Public Types inherited from Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > | |
typedef Scalar | scalar_type |
The type of the entries of the input and output multivectors. More... | |
typedef LocalOrdinal | local_ordinal_type |
The local index type. More... | |
typedef GlobalOrdinal | global_ordinal_type |
The global index type. More... | |
typedef Node | node_type |
The Kokkos Node type. More... | |
Public Member Functions | |
virtual const Teuchos::RCP < const Xpetra::Map < LocalOrdinal, GlobalOrdinal, Node > > | getMap () const =0 |
Implements DistObject interface. More... | |
virtual void | doImport (const Matrix &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)=0 |
Import. More... | |
virtual void | doExport (const Matrix &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)=0 |
Export. More... | |
virtual void | doImport (const Matrix &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)=0 |
Import (using an Exporter). More... | |
virtual void | doExport (const Matrix &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)=0 |
Export (using an Importer). More... | |
virtual bool | hasCrsGraph () const =0 |
Supports the getCrsGraph() call. More... | |
virtual RCP< const CrsGraph > | getCrsGraph () const =0 |
Returns the CrsGraph associated with this matrix. More... | |
virtual void | apply (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP< Import< LocalOrdinal, GlobalOrdinal, Node > > ®ionInterfaceImporter, const Teuchos::ArrayRCP< LocalOrdinal > ®ionInterfaceLIDs) const =0 |
Computes the matrix-multivector multiplication for region layout matrices. More... | |
void | SetFixedBlockSize (LocalOrdinal blksize, GlobalOrdinal offset=0) |
LocalOrdinal | GetFixedBlockSize () const |
bool | IsFixedBlockSizeSet () const |
Returns true, if SetFixedBlockSize has been called before. More... | |
virtual LocalOrdinal | GetStorageBlockSize () const =0 |
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatrix. More... | |
virtual void | SetMaxEigenvalueEstimate (Scalar const &sigma) |
virtual Scalar | GetMaxEigenvalueEstimate () const |
virtual local_matrix_type | getLocalMatrixDevice () const =0 |
virtual local_matrix_type::HostMirror | getLocalMatrixHost () const =0 |
virtual void | residual (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const =0 |
Compute a residual R = B - (*this) * X. More... | |
Public Member Functions inherited from Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > | |
virtual | ~Operator () |
virtual void | removeEmptyProcessesInPlace (const RCP< const map_type > &) |
virtual const Teuchos::RCP < const map_type > | getDomainMap () const =0 |
The Map associated with the domain of this operator, which must be compatible with X.getMap(). More... | |
virtual const Teuchos::RCP < const map_type > | getRangeMap () const =0 |
The Map associated with the range of this operator, which must be compatible with Y.getMap(). More... | |
virtual void | apply (const mv_type &X, mv_type &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const =0 |
Computes the operator-multivector application. More... | |
virtual bool | hasTransposeApply () const |
Whether this operator supports applying the transpose or conjugate transpose. More... | |
Protected Attributes | |
Teuchos::Hashtable < viewLabel_t, RCP< MatrixView > > | operatorViewTable_ |
viewLabel_t | defaultViewLabel_ |
viewLabel_t | currentViewLabel_ |
Private Types | |
typedef Xpetra::Map < LocalOrdinal, GlobalOrdinal, Node > | Map |
typedef Xpetra::CrsMatrix < Scalar, LocalOrdinal, GlobalOrdinal, Node > | CrsMatrix |
typedef Xpetra::CrsGraph < LocalOrdinal, GlobalOrdinal, Node > | CrsGraph |
typedef Xpetra::CrsMatrixFactory < Scalar, LocalOrdinal, GlobalOrdinal, Node > | CrsMatrixFactory |
typedef Xpetra::MatrixView < Scalar, LocalOrdinal, GlobalOrdinal, Node > | MatrixView |
Constructor/Destructor Methods | |
Matrix () | |
virtual | ~Matrix () |
Destructor. More... | |
View management methods | |
void | CreateView (viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap) |
void | CreateView (const viewLabel_t viewLabel, const RCP< const Matrix > &A, bool transposeA=false, const RCP< const Matrix > &B=Teuchos::null, bool transposeB=false) |
void | PrintViews (Teuchos::FancyOStream &out) const |
Print all of the views associated with the Matrix. More... | |
void | RemoveView (const viewLabel_t viewLabel) |
const viewLabel_t | SwitchToView (const viewLabel_t viewLabel) |
bool | IsView (const viewLabel_t viewLabel) const |
const viewLabel_t | SwitchToDefaultView () |
const viewLabel_t & | GetDefaultViewLabel () const |
const viewLabel_t & | GetCurrentViewLabel () const |
Insertion/Removal Methods | |
virtual void | insertGlobalValues (GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0 |
Insert matrix entries, using global IDs. More... | |
virtual void | insertLocalValues (LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0 |
Insert matrix entries, using local IDs. More... | |
virtual void | replaceGlobalValues (GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0 |
Replace matrix entries, using global IDs. More... | |
virtual void | replaceLocalValues (LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0 |
Replace matrix entries, using local IDs. More... | |
virtual void | setAllToScalar (const Scalar &alpha)=0 |
Set all matrix entries equal to scalar. More... | |
virtual void | scale (const Scalar &alpha)=0 |
Scale the current values of a matrix, this = alpha*this. More... | |
Transformational Methods | |
virtual void | resumeFill (const RCP< ParameterList > ¶ms=null)=0 |
virtual void | fillComplete (const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > ¶ms=null)=0 |
Signal that data entry is complete, specifying domain and range maps. More... | |
virtual void | fillComplete (const RCP< ParameterList > ¶ms=null)=0 |
Signal that data entry is complete. More... | |
Methods implementing RowMatrix | |
virtual const RCP< const Map > & | getRowMap () const |
Returns the Map that describes the row distribution in this matrix. More... | |
virtual const RCP< const Map > & | getRowMap (viewLabel_t viewLabel) const |
Returns the Map that describes the row distribution in this matrix. More... | |
virtual const RCP< const Map > & | getColMap () const |
Returns the Map that describes the column distribution in this matrix. This might be null until fillComplete() is called. More... | |
virtual const RCP< const Map > & | getColMap (viewLabel_t viewLabel) const |
Returns the Map that describes the column distribution in this matrix. More... | |
virtual global_size_t | getGlobalNumRows () const =0 |
Returns the number of global rows in this matrix. More... | |
virtual global_size_t | getGlobalNumCols () const =0 |
Returns the number of global columns in the matrix. More... | |
virtual size_t | getLocalNumRows () const =0 |
Returns the number of matrix rows owned on the calling node. More... | |
virtual global_size_t | getGlobalNumEntries () const =0 |
Returns the global number of entries in this matrix. More... | |
virtual size_t | getLocalNumEntries () const =0 |
Returns the local number of entries in this matrix. More... | |
virtual size_t | getNumEntriesInLocalRow (LocalOrdinal localRow) const =0 |
Returns the current number of entries on this node in the specified local row. More... | |
virtual size_t | getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const =0 |
Returns the current number of entries in the specified global row. More... | |
virtual size_t | getGlobalMaxNumRowEntries () const =0 |
Returns the maximum number of entries across all rows/columns on all nodes. More... | |
virtual size_t | getLocalMaxNumRowEntries () const =0 |
Returns the maximum number of entries across all rows/columns on this node. More... | |
virtual bool | isLocallyIndexed () const =0 |
If matrix indices are in the local range, this function returns true. Otherwise, this function returns false. */. More... | |
virtual bool | isGloballyIndexed () const =0 |
If matrix indices are in the global range, this function returns true. Otherwise, this function returns false. */. More... | |
virtual bool | isFillComplete () const =0 |
Returns true if fillComplete() has been called and the matrix is in compute mode. More... | |
virtual void | getLocalRowCopy (LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const =0 |
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calling routine. More... | |
virtual void | getGlobalRowView (GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const =0 |
Extract a const, non-persisting view of global indices in a specified row of the matrix. More... | |
virtual void | getLocalRowView (LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const =0 |
Extract a const, non-persisting view of local indices in a specified row of the matrix. More... | |
virtual void | getLocalDiagCopy (Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const =0 |
Get a copy of the diagonal entries owned by this node, with local row indices. More... | |
virtual ScalarTraits< Scalar > ::magnitudeType | getFrobeniusNorm () const =0 |
Get Frobenius norm of the matrix. More... | |
virtual void | leftScale (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0 |
Left scale matrix using the given vector entries. More... | |
virtual void | rightScale (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0 |
Right scale matrix using the given vector entries. More... | |
virtual bool | haveGlobalConstants () const =0 |
Returns true if globalConstants have been computed; false otherwise. More... | |
Overridden from Teuchos::Describable | |
virtual std::string | description () const =0 |
Return a simple one-line description of this object. More... | |
virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const =0 |
Print the object with some verbosity level to an FancyOStream object. More... | |
Overridden from Teuchos::LabeledObject | |
virtual void | setObjectLabel (const std::string &objectLabel)=0 |
Xpetra-specific matrix class.
This class is specific to Xpetra and has no analogue in Epetra or Tpetra. The main motivation for this class is to be able to access matrix data in a manner different than how it is stored. For example, it might be more convenient to treat ("view") a matrix stored in compressed row storage as if it were a block matrix. The Xpetra::Matrix class is intended to manage these "views".
How to create a Matrix from an existing CrsMatrix
Definition at line 60 of file Xpetra_Matrix_decl.hpp.
|
private |
Definition at line 61 of file Xpetra_Matrix_decl.hpp.
|
private |
Definition at line 62 of file Xpetra_Matrix_decl.hpp.
|
private |
Definition at line 63 of file Xpetra_Matrix_decl.hpp.
|
private |
Definition at line 64 of file Xpetra_Matrix_decl.hpp.
|
private |
Definition at line 65 of file Xpetra_Matrix_decl.hpp.
typedef Scalar Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scalar_type |
Definition at line 68 of file Xpetra_Matrix_decl.hpp.
typedef LocalOrdinal Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type |
Definition at line 69 of file Xpetra_Matrix_decl.hpp.
typedef GlobalOrdinal Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type |
Definition at line 70 of file Xpetra_Matrix_decl.hpp.
typedef Node Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type |
Definition at line 71 of file Xpetra_Matrix_decl.hpp.
typedef CrsMatrix::local_matrix_type Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type |
Definition at line 74 of file Xpetra_Matrix_decl.hpp.
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Matrix | ( | ) |
Definition at line 35 of file Xpetra_Matrix_def.hpp.
|
virtual |
Destructor.
Definition at line 38 of file Xpetra_Matrix_def.hpp.
void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CreateView | ( | viewLabel_t | viewLabel, |
const RCP< const Map > & | rowMap, | ||
const RCP< const Map > & | colMap | ||
) |
Definition at line 41 of file Xpetra_Matrix_def.hpp.
void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CreateView | ( | const viewLabel_t | viewLabel, |
const RCP< const Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > & | A, | ||
bool | transposeA = false , |
||
const RCP< const Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > & | B = Teuchos::null , |
||
bool | transposeB = false |
||
) |
Definition at line 48 of file Xpetra_Matrix_def.hpp.
void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::PrintViews | ( | Teuchos::FancyOStream & | out | ) | const |
Print all of the views associated with the Matrix.
Definition at line 91 of file Xpetra_Matrix_def.hpp.
void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::RemoveView | ( | const viewLabel_t | viewLabel | ) |
Definition at line 105 of file Xpetra_Matrix_def.hpp.
const std::string Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SwitchToView | ( | const viewLabel_t | viewLabel | ) |
Definition at line 112 of file Xpetra_Matrix_def.hpp.
bool Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::IsView | ( | const viewLabel_t | viewLabel | ) | const |
Definition at line 120 of file Xpetra_Matrix_def.hpp.
const std::string Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SwitchToDefaultView | ( | ) |
Definition at line 125 of file Xpetra_Matrix_def.hpp.
const std::string & Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetDefaultViewLabel | ( | ) | const |
Definition at line 128 of file Xpetra_Matrix_def.hpp.
const std::string & Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetCurrentViewLabel | ( | ) | const |
Definition at line 131 of file Xpetra_Matrix_def.hpp.
|
pure virtual |
Insert matrix entries, using global IDs.
All index values must be in the global space.
globalRow
exists as an ID in the global row map isLocallyIndexed() == false
isStorageOptimized() == false
isGloballyIndexed() == true
globalRow
does not belong to the matrix on this node, then it will be communicated to the appropriate node when globalAssemble() is called (which will, at the latest, occur during the next call to fillComplete().) Otherwise, the entries will be inserted in the local matrix. cols
, then the new values will be summed with the old values; this may happen at insertion or during the next call to fillComplete(). hasColMap() == true
, only (cols[i],vals[i]) where cols[i] belongs to the column map on this node will be inserted into the matrix. Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Insert matrix entries, using local IDs.
All index values must be in the local space.
localRow
exists as an ID in the local row map isGloballyIndexed() == false
isStorageOptimized() == false
isLocallyIndexed() == true
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Replace matrix entries, using global IDs.
All index values must be in the global space.
globalRow
is a global row belonging to the matrix on this node.Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Replace matrix entries, using local IDs.
All index values must be in the local space. Note that if a value is not already present for the specified location in the matrix, the input value will be ignored silently.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Set all matrix entries equal to scalar.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Scale the current values of a matrix, this = alpha*this.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Resume fill operations. After calling fillComplete(), resumeFill() must be called before initiating any changes to the matrix.
resumeFill() may be called repeatedly.
isFillActive() == true
isFillComplete() == false
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Signal that data entry is complete, specifying domain and range maps.
Off-node indices are distributed (via globalAssemble()), indices are sorted, redundant indices are eliminated, and global indices are transformed to local indices.
isFillActive() == true
isFillComplete()() == false
isFillActive() == false
isFillComplete() == true
if os == DoOptimizeStorage, then isStorageOptimized() == true
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Signal that data entry is complete.
Off-node entries are distributed (via globalAssemble()), repeated entries are summed, and global indices are transformed to local indices.
isFillActive() == true
isFillComplete()() == false
isFillActive() == false
isFillComplete() == true
if os == DoOptimizeStorage, then isStorageOptimized() == true
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
virtual |
Returns the Map that describes the row distribution in this matrix.
Definition at line 134 of file Xpetra_Matrix_def.hpp.
|
virtual |
Returns the Map that describes the row distribution in this matrix.
Definition at line 137 of file Xpetra_Matrix_def.hpp.
|
virtual |
Returns the Map that describes the column distribution in this matrix. This might be null
until fillComplete() is called.
Reimplemented in Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 143 of file Xpetra_Matrix_def.hpp.
|
virtual |
Returns the Map that describes the column distribution in this matrix.
Reimplemented in Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 146 of file Xpetra_Matrix_def.hpp.
|
pure virtual |
Returns the number of global rows in this matrix.
Undefined if isFillActive().
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns the number of global columns in the matrix.
Undefined if isFillActive().
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns the number of matrix rows owned on the calling node.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns the global number of entries in this matrix.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns the local number of entries in this matrix.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns the current number of entries on this node in the specified local row.
Returns OrdinalTraits<size_t>::invalid() if the specified local row is not valid for this matrix.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns the current number of entries in the specified global row.
Returns OrdinalTraits<size_t>::invalid() if the specified global row is not owned by this process.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns the maximum number of entries across all rows/columns on all nodes.
Undefined if isFillActive().
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns the maximum number of entries across all rows/columns on this node.
Undefined if isFillActive().
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
If matrix indices are in the local range, this function returns true. Otherwise, this function returns false. */.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
If matrix indices are in the global range, this function returns true. Otherwise, this function returns false. */.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns true
if fillComplete() has been called and the matrix is in compute mode.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calling routine.
LocalRow | - (In) Local row number for which indices are desired. |
Indices | - (Out) Local column indices corresponding to values. |
Values | - (Out) Matrix values. |
NumIndices | - (Out) Number of indices. |
Note: A std::runtime_error exception is thrown if either Indices
or Values
is not large enough to hold the data associated with row LocalRow
. If LocalRow
is not valid for this node, then Indices
and Values
are unchanged and NumIndices
is returned as OrdinalTraits<size_t>::invalid().
isLocallyIndexed()==true
or hasColMap() == true
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Extract a const, non-persisting view of global indices in a specified row of the matrix.
GlobalRow | - (In) Global row number for which indices are desired. |
Indices | - (Out) Global column indices corresponding to values. |
Values | - (Out) Row values |
isLocallyIndexed() == false
indices.size() == getNumEntriesInGlobalRow(GlobalRow)
Note: If GlobalRow
does not belong to this node, then indices
is set to null.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Extract a const, non-persisting view of local indices in a specified row of the matrix.
LocalRow | - (In) Local row number for which indices are desired. |
Indices | - (Out) Local column indices corresponding to values. |
Values | - (Out) Row values |
isGloballyIndexed() == false
indices.size() == getNumEntriesInLocalRow(LocalRow)
Note: If LocalRow
does not belong to this node, then indices
is set to null.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Get a copy of the diagonal entries owned by this node, with local row indices.
Returns a distributed Vector object partitioned according to this matrix's row map, containing the zero and non-zero diagonals owned by this node.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Get Frobenius norm of the matrix.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Left scale matrix using the given vector entries.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Right scale matrix using the given vector entries.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns true if globalConstants have been computed; false otherwise.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Implements DistObject interface.
Access function for the Tpetra::Map this DistObject was constructed with.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
|
pure virtual |
|
pure virtual |
Import (using an Exporter).
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Export (using an Importer).
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Return a simple one-line description of this object.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Print the object with some verbosity level to an FancyOStream object.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
|
pure virtual |
Supports the getCrsGraph() call.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Returns the CrsGraph associated with this matrix.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Computes the matrix-multivector multiplication for region layout matrices.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetFixedBlockSize | ( | LocalOrdinal | blksize, |
GlobalOrdinal | offset = 0 |
||
) |
Set fixed block size of operator (e.g., 3 for 3 DOFs per node).
blksize,: | block size denoting how many DOFs per node are used (LocalOrdinal) |
offset,: | global offset allows to define operators with global indices starting from a given value "offset" instead of 0. (GlobalOrdinal, default = 0) |
Definition at line 152 of file Xpetra_Matrix_def.hpp.
LocalOrdinal Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetFixedBlockSize | ( | ) | const |
Definition at line 174 of file Xpetra_Matrix_def.hpp.
bool Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::IsFixedBlockSizeSet | ( | ) | const |
Returns true, if SetFixedBlockSize
has been called before.
Definition at line 188 of file Xpetra_Matrix_def.hpp.
|
pure virtual |
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatrix.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
virtual |
Definition at line 193 of file Xpetra_Matrix_def.hpp.
|
virtual |
Definition at line 198 of file Xpetra_Matrix_def.hpp.
|
pure virtual |
|
pure virtual |
|
pure virtual |
Compute a residual R = B - (*this) * X.
Implements Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
protected |
Definition at line 441 of file Xpetra_Matrix_decl.hpp.
|
protected |
Definition at line 443 of file Xpetra_Matrix_decl.hpp.
|
protected |
Definition at line 444 of file Xpetra_Matrix_decl.hpp.