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

Concrete implementation of Xpetra::Matrix. More...

#include <Xpetra_CrsMatrixWrap_fwd.hpp>

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

Public Member Functions

global_size_t getGlobalNumRows () const
 Returns the number of global rows in this matrix. More...
 
global_size_t getGlobalNumCols () const
 Returns the number of global columns in the matrix. More...
 
size_t getNodeNumRows () 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 getNodeNumEntries () 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 global row. More...
 
size_t getGlobalMaxNumRowEntries () const
 Returns the maximum number of entries across all rows/columns on all nodes. More...
 
size_t getNodeMaxNumRowEntries () const
 Returns the maximum number of entries across all rows/columns on this node. More...
 
bool isLocallyIndexed () const
 If matrix indices 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...
 
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 (Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
 Get a copy of the diagonal entries owned by this node, with local row idices. More...
 
void getLocalDiagOffsets (Teuchos::ArrayRCP< size_t > &offsets) const
 Get offsets of the diagonal entries in the matrix. More...
 
void getLocalDiagCopy (Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const
 Get a copy of the diagonal entries owned by this node, with local row indices, using row offsets. More...
 
ScalarTraits< Scalar >
::magnitudeType 
getFrobeniusNorm () const
 Get Frobenius norm of the matrix. More...
 
void leftScale (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
 Left scale matrix using the given vector entries. More...
 
void rightScale (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
 Right scale matrix using the given vector entries. More...
 
bool haveGlobalConstants () const
 Returns true if globalConstants have been computed; false otherwise. More...
 
const Teuchos::RCP< const
Xpetra::Map< LocalOrdinal,
GlobalOrdinal, Node > > 
getMap () const
 Implements DistObject interface. More...
 
void doImport (const Matrix &source, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
 Import. More...
 
void doExport (const Matrix &dest, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
 Export. More...
 
void doImport (const Matrix &source, const Xpetra::Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
 Import (using an Exporter). More...
 
void doExport (const Matrix &dest, const Xpetra::Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
 Export (using an Importer). More...
 
bool hasCrsGraph () const
 Supports the getCrsGraph() call. More...
 
RCP< const CrsGraphgetCrsGraph () const
 Returns the CrsGraph associated with this matrix. More...
 
RCP< CrsMatrixgetCrsMatrix () const
 
void residual (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &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
 
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_tGetDefaultViewLabel () const
 
const viewLabel_tGetCurrentViewLabel () const
 
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 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...
 
- Public Member Functions inherited from Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >
virtual ~Operator ()
 
virtual void removeEmptyProcessesInPlace (const RCP< const Map > &)
 
virtual void residual (const MultiVector &X, const MultiVector &B, MultiVector &R) const =0
 Compute a residual R = B - (*this) * X. More...
 
virtual bool hasTransposeApply () const
 Whether this operator supports applying the transpose or conjugate transpose. More...
 

Private Types

typedef Xpetra::Map
< LocalOrdinal, GlobalOrdinal,
Node > 
Map
 
typedef Xpetra::CrsMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > 
CrsMatrix
 
typedef Xpetra::Matrix< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > 
Matrix
 
typedef Xpetra::CrsGraph
< LocalOrdinal, GlobalOrdinal,
Node > 
CrsGraph
 
typedef
Xpetra::TpetraCrsMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > 
TpetraCrsMatrix
 
typedef
Xpetra::CrsMatrixFactory
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > 
CrsMatrixFactory
 
typedef Xpetra::MatrixView
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > 
MatrixView
 

Private Member Functions

void CreateDefaultView ()
 
void updateDefaultView () const
 

Private Attributes

bool finalDefaultView_
 
RCP< CrsMatrixmatrixData_
 

Constructor/Destructor Methods

 CrsMatrixWrap (const RCP< const Map > &rowMap)
 Constructor for a dynamic profile matrix (Epetra only) More...
 
 CrsMatrixWrap (const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow)
 Constructor specifying fixed number of entries for each row. More...
 
 CrsMatrixWrap (const RCP< const Map > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc)
 Constructor specifying (possibly different) number of entries in each row. More...
 
 CrsMatrixWrap (const RCP< const Map > &rowMap, const RCP< const Map > &colMap, size_t maxNumEntriesPerRow)
 Constructor specifying fixed number of entries for each row and column map. More...
 
 CrsMatrixWrap (const RCP< const Map > &rowMap, const RCP< const Map > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc)
 Constructor specifying fixed number of entries for each row and column map. More...
 
 CrsMatrixWrap (RCP< CrsMatrix > matrix)
 
 CrsMatrixWrap (const RCP< const CrsGraph > &graph, const RCP< ParameterList > &paramList=Teuchos::null)
 
virtual ~CrsMatrixWrap ()
 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 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 > &params=null)
 
void fillComplete (const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< Teuchos::ParameterList > &params=null)
 Signal that data entry is complete, specifying domain and range maps. More...
 
void fillComplete (const RCP< ParameterList > &params=null)
 Signal that data entry is complete. More...
 

Methods implementing Matrix

Multiplies this matrix by a MultiVector.

\c X is required to be post-imported, i.e., described by the column map of the matrix. \c Y is required to be pre-exported, i.e., described by the row map of the matrix.

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 beta is equal to zero, the operation will enjoy overwrite semantics (Y will be overwritten with the result of the multiplication). Otherwise, the result of the multiplication will be accumulated into Y.

virtual void apply (const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &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 Xpetra::Map
< LocalOrdinal, GlobalOrdinal,
Node > > 
getDomainMap () const
 Returns the Map associated with the domain of this operator. This will be null until fillComplete() is called. More...
 
RCP< const Xpetra::Map
< LocalOrdinal, GlobalOrdinal,
Node > > 
getRangeMap () const
 
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...
 
const RCP< const Map > & getColMap (viewLabel_t viewLabel) const
 Returns the Map that describes the column distribution in this matrix. More...
 
void removeEmptyProcessesInPlace (const Teuchos::RCP< const Map > &newMap)
 

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)
 

Additional Inherited Members

- 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
 
- 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...
 
- Protected Attributes inherited from Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >
Teuchos::Hashtable
< viewLabel_t, RCP< MatrixView > > 
operatorViewTable_
 
viewLabel_t defaultViewLabel_
 
viewLabel_t currentViewLabel_
 

Detailed Description

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

Concrete implementation of Xpetra::Matrix.

Definition at line 51 of file Xpetra_CrsMatrixWrap_fwd.hpp.

Member Typedef Documentation

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

Definition at line 86 of file Xpetra_CrsMatrixWrap_decl.hpp.

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

Definition at line 87 of file Xpetra_CrsMatrixWrap_decl.hpp.

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

Definition at line 88 of file Xpetra_CrsMatrixWrap_decl.hpp.

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

Definition at line 89 of file Xpetra_CrsMatrixWrap_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::TpetraCrsMatrix
private

Definition at line 91 of file Xpetra_CrsMatrixWrap_decl.hpp.

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

Definition at line 93 of file Xpetra_CrsMatrixWrap_decl.hpp.

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

Definition at line 94 of file Xpetra_CrsMatrixWrap_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrixWrap ( const RCP< const Map > &  rowMap)

Constructor for a dynamic profile matrix (Epetra only)

Definition at line 71 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrixWrap ( const RCP< const Map > &  rowMap,
size_t  maxNumEntriesPerRow 
)

Constructor specifying fixed number of entries for each row.

Definition at line 79 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrixWrap ( const RCP< const Map > &  rowMap,
const ArrayRCP< const size_t > &  NumEntriesPerRowToAlloc 
)

Constructor specifying (possibly different) number of entries in each row.

Definition at line 88 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrixWrap ( const RCP< const Map > &  rowMap,
const RCP< const Map > &  colMap,
size_t  maxNumEntriesPerRow 
)

Constructor specifying fixed number of entries for each row and column map.

Definition at line 97 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrixWrap ( const RCP< const Map > &  rowMap,
const RCP< const Map > &  colMap,
const ArrayRCP< const size_t > &  NumEntriesPerRowToAlloc 
)

Constructor specifying fixed number of entries for each row and column map.

Definition at line 108 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrixWrap ( RCP< CrsMatrix matrix)

Definition at line 151 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrixWrap ( const RCP< const CrsGraph > &  graph,
const RCP< ParameterList > &  paramList = Teuchos::null 
)

Definition at line 162 of file Xpetra_CrsMatrixWrap_def.hpp.

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

Destructor.

Definition at line 173 of file Xpetra_CrsMatrixWrap_def.hpp.

Member Function Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::insertGlobalValues ( GlobalOrdinal  globalRow,
const ArrayView< const GlobalOrdinal > &  cols,
const ArrayView< const Scalar > &  vals 
)
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.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 176 of file Xpetra_CrsMatrixWrap_def.hpp.

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

Insert matrix entries, using local IDs.

All index values must be in the local space.

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

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 181 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceGlobalValues ( GlobalOrdinal  globalRow,
const ArrayView< const GlobalOrdinal > &  cols,
const ArrayView< const Scalar > &  vals 
)
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.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 186 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceLocalValues ( LocalOrdinal  localRow,
const ArrayView< const LocalOrdinal > &  cols,
const ArrayView< const Scalar > &  vals 
)
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 191 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setAllToScalar ( const Scalar &  alpha)
virtual

Set all matrix entries equal to scalar.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 196 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scale ( const Scalar &  alpha)
virtual

Scale the current values of a matrix, this = alpha*this.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 199 of file Xpetra_CrsMatrixWrap_def.hpp.

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

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 204 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fillComplete ( const RCP< const Map > &  domainMap,
const RCP< const Map > &  rangeMap,
const RCP< Teuchos::ParameterList > &  params = null 
)

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

Definition at line 209 of file Xpetra_CrsMatrixWrap_def.hpp.

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

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 217 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumRows ( ) const
virtual

Returns the number of global rows in this matrix.

Undefined if isFillActive().

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 225 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumCols ( ) const
virtual

Returns the number of global columns in the matrix.

Undefined if isFillActive().

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 230 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNodeNumRows ( ) const
virtual

Returns the number of matrix rows owned on the calling node.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 235 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumEntries ( ) const
virtual

Returns the global number of entries in this matrix.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 240 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNodeNumEntries ( ) const
virtual

Returns the local number of entries in this matrix.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 245 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInLocalRow ( LocalOrdinal  localRow) const
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 250 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInGlobalRow ( GlobalOrdinal  globalRow) const
virtual

Returns the current number of entries in the specified global row.

Returns OrdinalTraits<size_t>::invalid() if the row is not owned by this process.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 255 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalMaxNumRowEntries ( ) const
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 260 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNodeMaxNumRowEntries ( ) const
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 265 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isLocallyIndexed ( ) const
virtual

If matrix indices are in the local range, this function returns true. Otherwise, this function returns false. */.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 270 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isGloballyIndexed ( ) const
virtual

If matrix indices are in the global range, this function returns true. Otherwise, this function returns false. */.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 275 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isFillComplete ( ) const
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 280 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowCopy ( LocalOrdinal  LocalRow,
const ArrayView< LocalOrdinal > &  Indices,
const ArrayView< Scalar > &  Values,
size_t &  NumEntries 
) const
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

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 285 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowView ( GlobalOrdinal  GlobalRow,
ArrayView< const GlobalOrdinal > &  indices,
ArrayView< const Scalar > &  values 
) const
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.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 294 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowView ( LocalOrdinal  LocalRow,
ArrayView< const LocalOrdinal > &  indices,
ArrayView< const Scalar > &  values 
) const
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) Global 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.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 299 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagCopy ( Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  diag) const
virtual

Get a copy of the diagonal entries owned by this node, with local row idices.

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 304 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagOffsets ( Teuchos::ArrayRCP< size_t > &  offsets) const

Get offsets of the diagonal entries in the matrix.

Definition at line 309 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagCopy ( Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  diag,
const Teuchos::ArrayView< const size_t > &  offsets 
) const

Get a copy of the diagonal entries owned by this node, with local row indices, using row offsets.

Definition at line 314 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
ScalarTraits< Scalar >::magnitudeType Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getFrobeniusNorm ( ) const
virtual

Get Frobenius norm of the matrix.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 319 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::leftScale ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  x)
virtual

Left scale matrix using the given vector entries.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 324 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::rightScale ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  x)
virtual

Right scale matrix using the given vector entries.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 329 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::haveGlobalConstants ( ) const
virtual

Returns true if globalConstants have been computed; false otherwise.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 334 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::apply ( const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  X,
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  Y,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
Scalar  alpha = ScalarTraits<Scalar>::one(),
Scalar  beta = ScalarTraits<Scalar>::zero() 
) const
virtual

Computes the sparse matrix-multivector multiplication.

Performs \(Y = \alpha A^{\textrm{mode}} X + \beta Y\), with one special exceptions:

  • if beta == 0, apply() overwrites Y, so that any values in Y (including NaNs) are ignored.

Implements Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 339 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getDomainMap ( ) const
virtual

Returns the Map associated with the domain of this operator. This will be null until fillComplete() is called.

Implements Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 349 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRangeMap ( ) const
virtual

Returns the Map associated with the domain of this operator. This will be null until fillComplete() is called.

Implements Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 354 of file Xpetra_CrsMatrixWrap_def.hpp.

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

Returns the Map that describes the column distribution in this matrix. This might be null until fillComplete() is called.

Reimplemented from Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 359 of file Xpetra_CrsMatrixWrap_def.hpp.

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

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

Reimplemented from Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 362 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::removeEmptyProcessesInPlace ( const Teuchos::RCP< const Map > &  newMap)

Definition at line 369 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getMap ( ) const
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 376 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::doImport ( const Matrix source,
const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
CombineMode  CM 
)
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::doExport ( const Matrix dest,
const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
CombineMode  CM 
)
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::doImport ( const Matrix source,
const Xpetra::Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter,
CombineMode  CM 
)
virtual

Import (using an Exporter).

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 395 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::doExport ( const Matrix dest,
const Xpetra::Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter,
CombineMode  CM 
)
virtual

Export (using an Importer).

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 402 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
std::string Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::description ( ) const
virtual

Return a simple one-line description of this object.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 409 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
virtual

Print the object with some verbosity level to an FancyOStream object.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 414 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setObjectLabel ( const std::string &  objectLabel)
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::hasCrsGraph ( ) const
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< const Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getCrsGraph ( ) const
virtual

Returns the CrsGraph associated with this matrix.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 453 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getCrsMatrix ( ) const

Definition at line 457 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< 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
virtual

Compute a residual R = B - (*this) * X.

Implements Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 483 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CreateDefaultView ( )
private

Definition at line 462 of file Xpetra_CrsMatrixWrap_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::updateDefaultView ( ) const
private

Definition at line 473 of file Xpetra_CrsMatrixWrap_def.hpp.

Member Data Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
bool Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::finalDefaultView_
mutableprivate

Definition at line 493 of file Xpetra_CrsMatrixWrap_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP<CrsMatrix> Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node >::matrixData_
private

Definition at line 496 of file Xpetra_CrsMatrixWrap_decl.hpp.


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