Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Class Template Referenceabstract

Xpetra-specific matrix class. More...

#include <Xpetra_Matrix_fwd.hpp>

Inheritance diagram for Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >:
Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node > Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >

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 CrsGraphgetCrsGraph () 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 > > &regionInterfaceImporter, const Teuchos::ArrayRCP< LocalOrdinal > &regionInterfaceLIDs) 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_tGetDefaultViewLabel () const
 
const viewLabel_tGetCurrentViewLabel () 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 > &params=null)=0
 
virtual void fillComplete (const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)=0
 Signal that data entry is complete, specifying domain and range maps. More...
 
virtual void fillComplete (const RCP< ParameterList > &params=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
 

Detailed Description

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
class Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >

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

RCP<Xpetra::CrsMatrix> crsA;
RCP<Xpetra::Matrix> A = rcp(new CrsMatrixWrap(crsA));

Definition at line 51 of file Xpetra_Matrix_fwd.hpp.

Member Typedef Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Map
private

Definition at line 98 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix
private

Definition at line 99 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef Xpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node> Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsGraph
private

Definition at line 100 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef Xpetra::CrsMatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node> Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrixFactory
private

Definition at line 101 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef Xpetra::MatrixView<Scalar, LocalOrdinal, GlobalOrdinal, Node> Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::MatrixView
private

Definition at line 102 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef Scalar Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scalar_type

Definition at line 105 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef LocalOrdinal Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type

Definition at line 106 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef GlobalOrdinal Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type

Definition at line 107 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef Node Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type

Definition at line 108 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef CrsMatrix::local_matrix_type Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type

Definition at line 111 of file Xpetra_Matrix.hpp.

Constructor & Destructor Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Matrix ( )
inline

Definition at line 117 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::~Matrix ( )
inlinevirtual

Destructor.

Definition at line 120 of file Xpetra_Matrix.hpp.

Member Function Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CreateView ( viewLabel_t  viewLabel,
const RCP< const Map > &  rowMap,
const RCP< const Map > &  colMap 
)
inline

Definition at line 126 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
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 
)
inline

Definition at line 133 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::PrintViews ( Teuchos::FancyOStream &  out) const
inline

Print all of the views associated with the Matrix.

Definition at line 176 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::RemoveView ( const viewLabel_t  viewLabel)
inline

Definition at line 189 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
const viewLabel_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SwitchToView ( const viewLabel_t  viewLabel)
inline

Definition at line 195 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
bool Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::IsView ( const viewLabel_t  viewLabel) const
inline

Definition at line 202 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
const viewLabel_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SwitchToDefaultView ( )
inline

Definition at line 206 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
const viewLabel_t& Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetDefaultViewLabel ( ) const
inline

Definition at line 208 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
const viewLabel_t& Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetCurrentViewLabel ( ) const
inline

Definition at line 210 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::insertGlobalValues ( GlobalOrdinal  globalRow,
const ArrayView< const GlobalOrdinal > &  cols,
const ArrayView< const Scalar > &  vals 
)
pure virtual

Insert matrix entries, using global IDs.

All index values must be in the global space.

Precondition
globalRow exists as an ID in the global row map
isLocallyIndexed() == false
isStorageOptimized() == false
Postcondition
isGloballyIndexed() == true
Note
If 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.
If the matrix row already contains values at the indices corresponding to values in cols, then the new values will be summed with the old values; this may happen at insertion or during the next call to fillComplete().
If 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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::insertLocalValues ( LocalOrdinal  localRow,
const ArrayView< const LocalOrdinal > &  cols,
const ArrayView< const Scalar > &  vals 
)
pure virtual

Insert matrix entries, using local IDs.

All index values must be in the local space.

Precondition
localRow exists as an ID in the local row map
isGloballyIndexed() == false
isStorageOptimized() == false
Postcondition
isLocallyIndexed() == true

Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceGlobalValues ( GlobalOrdinal  globalRow,
const ArrayView< const GlobalOrdinal > &  cols,
const ArrayView< const Scalar > &  vals 
)
pure virtual

Replace matrix entries, using global IDs.

All index values must be in the global space.

Precondition
globalRow is a global row belonging to the matrix on this node.
Note
If (globalRow,cols[i]) corresponds to an entry that is duplicated in this matrix row (likely because it was inserted more than once and fillComplete() has not been called in the interim), the behavior of this function is not defined.

Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceLocalValues ( LocalOrdinal  localRow,
const ArrayView< const LocalOrdinal > &  cols,
const ArrayView< const Scalar > &  vals 
)
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setAllToScalar ( const Scalar &  alpha)
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scale ( const Scalar &  alpha)
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::resumeFill ( const RCP< ParameterList > &  params = null)
pure virtual

Resume fill operations. After calling fillComplete(), resumeFill() must be called before initiating any changes to the matrix.

resumeFill() may be called repeatedly.

Postcondition
isFillActive() == true
isFillComplete() == false

Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fillComplete ( const RCP< const Map > &  domainMap,
const RCP< const Map > &  rangeMap,
const RCP< ParameterList > &  params = null 
)
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.

Precondition
isFillActive() == true
isFillComplete()() == false
Postcondition
isFillActive() == false
isFillComplete() == true
if os == DoOptimizeStorage, then isStorageOptimized() == true

Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fillComplete ( const RCP< ParameterList > &  params = null)
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.

Note
This method calls fillComplete( getRowMap(), getRowMap(), os ).
Precondition
isFillActive() == true
isFillComplete()() == false
Postcondition
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual const RCP<const Map>& Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRowMap ( ) const
inlinevirtual

Returns the Map that describes the row distribution in this matrix.

Definition at line 315 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual const RCP<const Map>& Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRowMap ( viewLabel_t  viewLabel) const
inlinevirtual

Returns the Map that describes the row distribution in this matrix.

Definition at line 318 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual const RCP<const Map>& Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getColMap ( ) const
inlinevirtual

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 325 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual const RCP<const Map>& Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getColMap ( viewLabel_t  viewLabel) const
inlinevirtual

Returns the Map that describes the column distribution in this matrix.

Reimplemented in Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 328 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual global_size_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumRows ( ) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual global_size_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumCols ( ) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual size_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalNumRows ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual global_size_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumEntries ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual size_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalNumEntries ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual size_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInLocalRow ( LocalOrdinal  localRow) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual size_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInGlobalRow ( GlobalOrdinal  globalRow) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual size_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalMaxNumRowEntries ( ) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual size_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalMaxNumRowEntries ( ) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual bool Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isLocallyIndexed ( ) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual bool Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isGloballyIndexed ( ) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual bool Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isFillComplete ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowCopy ( LocalOrdinal  LocalRow,
const ArrayView< LocalOrdinal > &  Indices,
const ArrayView< Scalar > &  Values,
size_t &  NumEntries 
) const
pure virtual

Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calling routine.

Parameters
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().

Precondition
isLocallyIndexed()==true or hasColMap() == true

Implemented in Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowView ( GlobalOrdinal  GlobalRow,
ArrayView< const GlobalOrdinal > &  indices,
ArrayView< const Scalar > &  values 
) const
pure virtual

Extract a const, non-persisting view of global indices in a specified row of the matrix.

Parameters
GlobalRow- (In) Global row number for which indices are desired.
Indices- (Out) Global column indices corresponding to values.
Values- (Out) Row values
Precondition
isLocallyIndexed() == false
Postcondition
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowView ( LocalOrdinal  LocalRow,
ArrayView< const LocalOrdinal > &  indices,
ArrayView< const Scalar > &  values 
) const
pure virtual

Extract a const, non-persisting view of local indices in a specified row of the matrix.

Parameters
LocalRow- (In) Local row number for which indices are desired.
Indices- (Out) Local column indices corresponding to values.
Values- (Out) Row values
Precondition
isGloballyIndexed() == false
Postcondition
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagCopy ( Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  diag) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual ScalarTraits<Scalar>::magnitudeType Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getFrobeniusNorm ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::leftScale ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  x)
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::rightScale ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  x)
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual bool Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::haveGlobalConstants ( ) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getMap ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::doImport ( const Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
CombineMode  CM 
)
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::doExport ( const Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  dest,
const Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
CombineMode  CM 
)
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::doImport ( const Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter,
CombineMode  CM 
)
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::doExport ( const Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  dest,
const Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter,
CombineMode  CM 
)
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual std::string Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::description ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setObjectLabel ( const std::string &  objectLabel)
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual bool Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::hasCrsGraph ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual RCP<const CrsGraph> Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getCrsGraph ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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 > > &  regionInterfaceImporter,
const Teuchos::ArrayRCP< LocalOrdinal > &  regionInterfaceLIDs 
) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetFixedBlockSize ( LocalOrdinal  blksize,
GlobalOrdinal  offset = 0 
)
inline

Set fixed block size of operator (e.g., 3 for 3 DOFs per node).

Parameters
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 513 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
LocalOrdinal Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetFixedBlockSize ( ) const
inline

Definition at line 537 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
bool Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::IsFixedBlockSizeSet ( ) const
inline

Returns true, if SetFixedBlockSize has been called before.

Definition at line 551 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual LocalOrdinal Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetStorageBlockSize ( ) const
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 >.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetMaxEigenvalueEstimate ( Scalar const &  sigma)
inlinevirtual

Definition at line 560 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual Scalar Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetMaxEigenvalueEstimate ( ) const
inlinevirtual

Definition at line 566 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual local_matrix_type Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalMatrixDevice ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual local_matrix_type::HostMirror Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalMatrixHost ( ) const
pure virtual
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::residual ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  X,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B,
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  R 
) const
pure virtual

Member Data Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::Hashtable<viewLabel_t, RCP<MatrixView> > Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::operatorViewTable_
protected

Definition at line 587 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
viewLabel_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::defaultViewLabel_
protected

Definition at line 589 of file Xpetra_Matrix.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
viewLabel_t Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::currentViewLabel_
protected

Definition at line 590 of file Xpetra_Matrix.hpp.


The documentation for this class was generated from the following files: