Xpetra
Version of the Day
|
#include <Xpetra_BlockedCrsMatrix_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::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > | |
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 | |
global_size_t | getGlobalNumRows () const |
Returns the number of global rows. More... | |
global_size_t | getGlobalNumCols () const |
Returns the number of global columns in the matrix. More... | |
size_t | getLocalNumRows () const |
Returns the number of matrix rows owned on the calling node. More... | |
global_size_t | getGlobalNumEntries () const |
Returns the global number of entries in this matrix. More... | |
size_t | getLocalNumEntries () const |
Returns the local number of entries in this matrix. More... | |
size_t | getNumEntriesInLocalRow (LocalOrdinal localRow) const |
Returns the current number of entries on this node in the specified local row. More... | |
size_t | getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const |
Returns the current number of entries in the specified (locally owned) global row. More... | |
size_t | getGlobalMaxNumRowEntries () const |
Returns the maximum number of entries across all rows/columns on all nodes. More... | |
size_t | getLocalMaxNumRowEntries () const |
Returns the maximum number of entries across all rows/columns on this node. More... | |
bool | isLocallyIndexed () const |
If matrix indices of all matrix blocks are in the local range, this function returns true. Otherwise, this function returns false. More... | |
bool | isGloballyIndexed () const |
If matrix indices are in the global range, this function returns true. Otherwise, this function returns false. More... | |
bool | isFillComplete () const |
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 |
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calling routine. More... | |
void | getGlobalRowView (GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const |
Extract a const, non-persisting view of global indices in a specified row of the matrix. More... | |
void | getLocalRowView (LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const |
Extract a const, non-persisting view of local indices in a specified row of the matrix. More... | |
void | getLocalDiagCopy (Vector &diag) const |
Get a copy of the diagonal entries owned by this node, with local row indices. More... | |
void | leftScale (const Vector &x) |
Left scale matrix using the given vector entries. More... | |
void | rightScale (const Vector &x) |
Right scale matrix using the given vector entries. More... | |
virtual ScalarTraits< Scalar > ::magnitudeType | getFrobeniusNorm () const |
Get Frobenius norm of the matrix. More... | |
virtual bool | haveGlobalConstants () const |
Returns true if globalConstants have been computed; false otherwise. More... | |
virtual void | bgs_apply (const MultiVector &X, MultiVector &Y, size_t row, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const |
Special multiplication routine (for BGS/Jacobi smoother) More... | |
const Teuchos::RCP< const Map > | getMap () const |
Implements DistObject interface. More... | |
void | doImport (const Matrix &source, const Import &importer, CombineMode CM) |
Import. More... | |
void | doExport (const Matrix &dest, const Import &importer, CombineMode CM) |
Export. More... | |
void | doImport (const Matrix &source, const Export &exporter, CombineMode CM) |
Import (using an Exporter). More... | |
void | doExport (const Matrix &dest, const Export &exporter, CombineMode CM) |
Export (using an Importer). More... | |
bool | hasCrsGraph () const |
Supports the getCrsGraph() call. More... | |
RCP< const CrsGraph > | getCrsGraph () const |
Returns the CrsGraph associated with this matrix. More... | |
local_matrix_type | getLocalMatrixDevice () const |
Access the underlying local Kokkos::CrsMatrix object. More... | |
local_matrix_type::HostMirror | getLocalMatrixHost () const |
Access the underlying local Kokkos::CrsMatrix object. More... | |
LocalOrdinal | GetStorageBlockSize () const |
Returns the block size of the storage mechanism. More... | |
void | residual (const MultiVector &X, const MultiVector &B, MultiVector &R) const |
Compute a residual R = B - (*this) * X. More... | |
Public Member Functions inherited from Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > | |
void | SetFixedBlockSize (LocalOrdinal blksize, GlobalOrdinal offset=0) |
LocalOrdinal | GetFixedBlockSize () const |
bool | IsFixedBlockSizeSet () const |
Returns true, if SetFixedBlockSize has been called before. More... | |
virtual void | SetMaxEigenvalueEstimate (Scalar const &sigma) |
virtual Scalar | GetMaxEigenvalueEstimate () const |
Matrix () | |
virtual | ~Matrix () |
Destructor. More... | |
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 |
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... | |
Public Member Functions inherited from Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > | |
virtual | ~Operator () |
virtual void | removeEmptyProcessesInPlace (const RCP< const map_type > &) |
virtual bool | hasTransposeApply () const |
Whether this operator supports applying the transpose or conjugate transpose. More... | |
Private Member Functions | |
void | CreateDefaultView () |
Private Attributes | |
bool | is_diagonal_ |
If we're diagonal, a bunch of the extraction stuff should work. More... | |
Teuchos::RCP< const MapExtractor > | domainmaps_ |
full domain map together with all partial domain maps More... | |
Teuchos::RCP< const MapExtractor > | rangemaps_ |
full range map together with all partial domain maps More... | |
std::vector< Teuchos::RCP < Matrix > > | blocks_ |
row major matrix block storage More... | |
bool | bRangeThyraMode_ |
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub blocks, i.e. all GIDs of submaps are contiguous and start from 0. More... | |
bool | bDomainThyraMode_ |
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub blocks, i.e. all GIDs of submaps are contiguous and start from 0. More... | |
Constructor/Destructor Methods | |
BlockedCrsMatrix (const Teuchos::RCP< const BlockedMap > &rangeMaps, const Teuchos::RCP< const BlockedMap > &domainMaps, size_t numEntriesPerRow) | |
Constructor. More... | |
BlockedCrsMatrix (Teuchos::RCP< const MapExtractor > &rangeMapExtractor, Teuchos::RCP< const MapExtractor > &domainMapExtractor, size_t numEntriesPerRow) | |
Constructor. More... | |
virtual | ~BlockedCrsMatrix () |
Destructor. More... | |
Insertion/Removal Methods | |
void | insertGlobalValues (GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) |
Insert matrix entries, using global IDs. More... | |
void | insertLocalValues (LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals) |
Insert matrix entries, using local IDs. More... | |
void | removeEmptyProcessesInPlace (const Teuchos::RCP< const Map > &newMap) |
void | replaceGlobalValues (GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) |
Replace matrix entries, using global IDs. More... | |
void | replaceLocalValues (LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals) |
Replace matrix entries, using local IDs. More... | |
virtual void | setAllToScalar (const Scalar &alpha) |
Set all matrix entries equal to scalar. More... | |
void | scale (const Scalar &alpha) |
Scale the current values of a matrix, this = alpha*this. More... | |
Transformational Methods | |
void | resumeFill (const RCP< ParameterList > ¶ms=null) |
void | fillComplete (const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > ¶ms=null) |
Signal that data entry is complete. More... | |
void | fillComplete (const RCP< ParameterList > ¶ms=null) |
Signal that data entry is complete. More... | |
Methods implementing Matrix | |
Multiplies this matrix by a MultiVector.
Both are required to have constant stride, and they are not permitted to ocupy overlapping space. No runtime checking will be performed in a non-debug build. This method is templated on the scalar type of MultiVector objects, allowing this method to be applied to MultiVector objects of arbitrary type. However, it is recommended that multiply() not be called directly; instead, use the CrsMatrixMultiplyOp, as it will handle the import/exprt operations required to apply a matrix with non-trivial communication needs. If | |
virtual void | apply (const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP< Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > ®ionInterfaceImporter, const Teuchos::ArrayRCP< LocalOrdinal > ®ionInterfaceLIDs) const |
sparse matrix-multivector multiplication for the region layout matrices (currently no blocked implementation) More... | |
virtual void | apply (const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const |
Computes the sparse matrix-multivector multiplication. More... | |
RCP< const Map > | getFullDomainMap () const |
Returns the Map associated with the full domain of this operator. More... | |
RCP< const BlockedMap > | getBlockedDomainMap () const |
Returns the BlockedMap associated with the domain of this operator. More... | |
const RCP< const Map > | getDomainMap () const |
Returns the Map associated with the domain of this operator. More... | |
RCP< const Map > | getDomainMap (size_t i) const |
Returns the Map associated with the i'th block domain of this operator. More... | |
RCP< const Map > | getDomainMap (size_t i, bool bThyraMode) const |
Returns the Map associated with the i'th block domain of this operator. More... | |
RCP< const Map > | getFullRangeMap () const |
Returns the Map associated with the full range of this operator. More... | |
RCP< const BlockedMap > | getBlockedRangeMap () const |
Returns the BlockedMap associated with the range of this operator. More... | |
const RCP< const Map > | getRangeMap () const |
Returns the Map associated with the range of this operator. More... | |
RCP< const Map > | getRangeMap (size_t i) const |
Returns the Map associated with the i'th block range of this operator. More... | |
RCP< const Map > | getRangeMap (size_t i, bool bThyraMode) const |
Returns the Map associated with the i'th block range of this operator. More... | |
RCP< const MapExtractor > | getRangeMapExtractor () const |
Returns map extractor class for range map. More... | |
RCP< const MapExtractor > | getDomainMapExtractor () const |
Returns map extractor for domain map. More... | |
Overridden from Teuchos::Describable | |
std::string | description () const |
Return a simple one-line description of this object. More... | |
void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const |
Print the object with some verbosity level to an FancyOStream object. More... | |
Overridden from Teuchos::LabeledObject | |
void | setObjectLabel (const std::string &objectLabel) |
Block matrix access | |
virtual bool | isDiagonal () const |
virtual size_t | Rows () const |
number of row blocks More... | |
virtual size_t | Cols () const |
number of column blocks More... | |
Teuchos::RCP< Matrix > | getCrsMatrix () const |
return unwrap 1x1 blocked operators More... | |
Teuchos::RCP< Matrix > | getInnermostCrsMatrix () |
helper routine recursively returns the first inner-most non-null matrix block from a (nested) blocked operator More... | |
Teuchos::RCP< Matrix > | getMatrix (size_t r, size_t c) const |
return block (r,c) More... | |
void | setMatrix (size_t r, size_t c, Teuchos::RCP< Matrix > mat) |
set matrix block More... | |
Teuchos::RCP< Matrix > | Merge () const |
merge BlockedCrsMatrix blocks in a CrsMatrix More... | |
helper functions | |
void | Add (const Matrix &A, const Scalar scalarA, Matrix &B, const Scalar scalarB) const |
Add a Xpetra::CrsMatrix to another: B = B*scalarB + A*scalarA. More... | |
Additional Inherited Members | |
Protected Attributes inherited from Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > | |
Teuchos::Hashtable < viewLabel_t, RCP< MatrixView > > | operatorViewTable_ |
viewLabel_t | defaultViewLabel_ |
viewLabel_t | currentViewLabel_ |
Definition at line 65 of file Xpetra_BlockedCrsMatrix_decl.hpp.
typedef Scalar Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scalar_type |
Definition at line 67 of file Xpetra_BlockedCrsMatrix_decl.hpp.
typedef LocalOrdinal Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type |
Definition at line 68 of file Xpetra_BlockedCrsMatrix_decl.hpp.
typedef GlobalOrdinal Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type |
Definition at line 69 of file Xpetra_BlockedCrsMatrix_decl.hpp.
typedef Node Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type |
Definition at line 70 of file Xpetra_BlockedCrsMatrix_decl.hpp.
typedef CrsMatrix::local_matrix_type Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type |
Definition at line 533 of file Xpetra_BlockedCrsMatrix_decl.hpp.
Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BlockedCrsMatrix | ( | const Teuchos::RCP< const BlockedMap > & | rangeMaps, |
const Teuchos::RCP< const BlockedMap > & | domainMaps, | ||
size_t | numEntriesPerRow | ||
) |
Constructor.
rangeMaps | range maps for all blocks |
domainMaps | domain maps for all blocks |
numEntriesPerRow | estimated number of entries per row in each block(!) |
Definition at line 56 of file Xpetra_BlockedCrsMatrix_def.hpp.
Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BlockedCrsMatrix | ( | Teuchos::RCP< const MapExtractor > & | rangeMapExtractor, |
Teuchos::RCP< const MapExtractor > & | domainMapExtractor, | ||
size_t | numEntriesPerRow | ||
) |
Constructor.
rangeMapExtractor | range map extractor for all blocks |
domainMapExtractor | domain map extractor for all blocks |
numEntriesPerRow | estimated number of entries per row in each block(!) |
Definition at line 77 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtualdefault |
Destructor.
|
virtual |
Insert matrix entries, using global IDs.
Note: this routine throws for Rows() > 1 and/or Cols() > 1
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. Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 225 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Insert matrix entries, using local IDs.
Note: this routine throws if Rows() > 1 and/or Cols() > 1
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
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 235 of file Xpetra_BlockedCrsMatrix_def.hpp.
void Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::removeEmptyProcessesInPlace | ( | const Teuchos::RCP< const Map > & | newMap | ) |
Definition at line 245 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
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.Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 255 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
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.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 267 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Set all matrix entries equal to scalar.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 280 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Scale the current values of a matrix, this = alpha*this.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 293 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Resume fill operations. After calling fillComplete(), resumeFill() must be called before initiating any changes to the matrix.
For BlockedCrsMatrix objects we call the routine iteratively for all sub-blocks.
resumeFill() may be called repeatedly.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 305 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Signal that data entry is complete.
Note: for blocked operators the specified domain and range maps have no meaning. We just call fillComplete for all underlying blocks
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 317 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
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
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 327 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the number of global rows.
Undefined if isFillActive().
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 383 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the number of global columns in the matrix.
Undefined if isFillActive().
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 398 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the number of matrix rows owned on the calling node.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 413 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the global number of entries in this matrix.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 428 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the local number of entries in this matrix.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 441 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
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.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 454 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the current number of entries in the specified (locally owned) global row.
Returns OrdinalTraits<size_t>::invalid() if the specified local row is not valid for this matrix.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 467 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the maximum number of entries across all rows/columns on all nodes.
Undefined if isFillActive().
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 478 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the maximum number of entries across all rows/columns on this node.
Undefined if isFillActive().
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 496 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
If matrix indices of all matrix blocks are in the local range, this function returns true. Otherwise, this function returns false.
if false, then this does not automatically mean that all blocks are globally indexed. The user has to make sure, that all matrix blocks are indexed in the same way (locally or globally). Otherwise the block matrix is not valid...
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 514 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
If matrix indices are in the global range, this function returns true. Otherwise, this function returns false.
if false, then this does not automatically mean that all blocks are locally indexed. The user has to make sure, that all matrix blocks are indexed in the same way (locally or globally). Otherwise the block matrix is not valid...
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 523 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns true
if fillComplete() has been called and the matrix is in compute mode.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 532 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
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
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 541 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
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.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 554 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
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) Global 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.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 564 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
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 the zero and non-zero diagonals owned by this node.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 579 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Left scale matrix using the given vector entries.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 614 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Right scale matrix using the given vector entries.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 652 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Get Frobenius norm of the matrix.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 690 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns true if globalConstants have been computed; false otherwise.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 705 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
sparse matrix-multivector multiplication for the region layout matrices (currently no blocked implementation)
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
virtual |
Computes the sparse matrix-multivector multiplication.
Performs \(Y = \alpha A^{\textrm{mode}} X + \beta Y\), with one special exceptions:
beta == 0
, apply() overwrites Y
, so that any values in Y
(including NaNs) are ignored. Implements Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Reimplemented in Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 713 of file Xpetra_BlockedCrsMatrix_def.hpp.
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getFullDomainMap | ( | ) | const |
Returns the Map associated with the full domain of this operator.
Definition at line 807 of file Xpetra_BlockedCrsMatrix_def.hpp.
RCP< const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getBlockedDomainMap | ( | ) | const |
Returns the BlockedMap associated with the domain of this operator.
Definition at line 813 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the Map associated with the domain of this operator.
Implements Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 819 of file Xpetra_BlockedCrsMatrix_def.hpp.
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getDomainMap | ( | size_t | i | ) | const |
Returns the Map associated with the i'th block domain of this operator.
Definition at line 825 of file Xpetra_BlockedCrsMatrix_def.hpp.
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getDomainMap | ( | size_t | i, |
bool | bThyraMode | ||
) | const |
Returns the Map associated with the i'th block domain of this operator.
Definition at line 831 of file Xpetra_BlockedCrsMatrix_def.hpp.
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getFullRangeMap | ( | ) | const |
Returns the Map associated with the full range of this operator.
Definition at line 837 of file Xpetra_BlockedCrsMatrix_def.hpp.
RCP< const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getBlockedRangeMap | ( | ) | const |
Returns the BlockedMap associated with the range of this operator.
Definition at line 843 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the Map associated with the range of this operator.
Implements Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 849 of file Xpetra_BlockedCrsMatrix_def.hpp.
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRangeMap | ( | size_t | i | ) | const |
Returns the Map associated with the i'th block range of this operator.
Definition at line 855 of file Xpetra_BlockedCrsMatrix_def.hpp.
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRangeMap | ( | size_t | i, |
bool | bThyraMode | ||
) | const |
Returns the Map associated with the i'th block range of this operator.
Definition at line 861 of file Xpetra_BlockedCrsMatrix_def.hpp.
RCP< const Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRangeMapExtractor | ( | ) | const |
Returns map extractor class for range map.
Definition at line 867 of file Xpetra_BlockedCrsMatrix_def.hpp.
RCP< const Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getDomainMapExtractor | ( | ) | const |
Returns map extractor for domain map.
Definition at line 873 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Special multiplication routine (for BGS/Jacobi smoother)
Computes the sparse matrix-multivector multiplication (plus linear combination with input/result vector)
Performs \(Y = \alpha A^{\textrm{mode}} X + \beta Y\), with one special exception:
beta == 0
, apply() overwrites Y
, so that any values in Y
(including NaNs) are ignored.X | Vector to be multiplied by matrix (input) |
Y | result vector |
row | Index of block row to be treated |
mode | Transpose mode |
alpha | scaling factor for result of matrix-vector product |
beta | scaling factor for linear combination with result vector |
Definition at line 879 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Implements DistObject interface.
Access function for the Tpetra::Map this DistObject was constructed with.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 944 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 953 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 963 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Import (using an Exporter).
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 973 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Export (using an Importer).
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 984 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Return a simple one-line description of this object.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Reimplemented in Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 994 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Print the object with some verbosity level to an FancyOStream object.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Reimplemented in Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 997 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1028 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Supports the getCrsGraph() call.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1041 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the CrsGraph associated with this matrix.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1049 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Definition at line 1058 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
number of row blocks
Definition at line 1061 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
number of column blocks
Definition at line 1067 of file Xpetra_BlockedCrsMatrix_def.hpp.
Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getCrsMatrix | ( | ) | const |
return unwrap 1x1 blocked operators
Definition at line 1073 of file Xpetra_BlockedCrsMatrix_def.hpp.
Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getInnermostCrsMatrix | ( | ) |
helper routine recursively returns the first inner-most non-null matrix block from a (nested) blocked operator
Definition at line 1085 of file Xpetra_BlockedCrsMatrix_def.hpp.
Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getMatrix | ( | size_t | r, |
size_t | c | ||
) | const |
return block (r,c)
Definition at line 1103 of file Xpetra_BlockedCrsMatrix_def.hpp.
void Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setMatrix | ( | size_t | r, |
size_t | c, | ||
Teuchos::RCP< Matrix > | mat | ||
) |
set matrix block
Definition at line 1116 of file Xpetra_BlockedCrsMatrix_def.hpp.
Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Merge | ( | ) | const |
merge BlockedCrsMatrix blocks in a CrsMatrix
Definition at line 1128 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Access the underlying local Kokkos::CrsMatrix object.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1247 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Access the underlying local Kokkos::CrsMatrix object.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1255 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Returns the block size of the storage mechanism.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1277 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
virtual |
Compute a residual R = B - (*this) * X.
Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1280 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
private |
Add a Xpetra::CrsMatrix to another: B = B*scalarB + A*scalarA.
Note, that this routine works only correctly if A only has entries which are empty (zero) in B. We use the insertGlobalValues routine for inserting the new values from A in B. The sumIntoGlobalValues routine is not implemented in Xpetra (and would not extend the graph of B for new entries). Here we need something to catch the exceptions of a future implementation of sumIntoGlobalValues that then adds the remaining new entries with insertGlobal Values.
This routine is private and used only by Merge. Since the blocks in BlockedCrsMatrix are seperated, this routine works for merging a BlockedCrsMatrix.
Definition at line 1289 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
private |
Definition at line 1333 of file Xpetra_BlockedCrsMatrix_def.hpp.
|
private |
If we're diagonal, a bunch of the extraction stuff should work.
Definition at line 574 of file Xpetra_BlockedCrsMatrix_decl.hpp.
|
private |
full domain map together with all partial domain maps
Definition at line 575 of file Xpetra_BlockedCrsMatrix_decl.hpp.
|
private |
full range map together with all partial domain maps
Definition at line 576 of file Xpetra_BlockedCrsMatrix_decl.hpp.
|
private |
row major matrix block storage
Definition at line 578 of file Xpetra_BlockedCrsMatrix_decl.hpp.
|
private |
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub blocks, i.e. all GIDs of submaps are contiguous and start from 0.
Definition at line 582 of file Xpetra_BlockedCrsMatrix_decl.hpp.
|
private |
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub blocks, i.e. all GIDs of submaps are contiguous and start from 0.
Definition at line 583 of file Xpetra_BlockedCrsMatrix_decl.hpp.