GEOSX
Protected Types | Protected Member Functions | Protected Attributes | Friends | List of all members
geos::MatrixBase< MATRIX, VECTOR > Class Template Referenceabstract

Common base template for all matrix wrapper types. More...

#include <MatrixBase.hpp>

Inheritance diagram for geos::MatrixBase< MATRIX, VECTOR >:
Inheritance graph
[legend]

Protected Types

using Base = LinearOperator< VECTOR >
 Alias for base type.
 
using Matrix = MATRIX
 Type alias for actual derived matrix class.
 
using Vector = VECTOR
 Type alias for a compatible vector class.
 

Protected Member Functions

virtual localIndex numLocalRows () const=0
 Get the number of local rows. More...
 
virtual localIndex numLocalCols () const=0
 Get the number of local columns. More...
 
virtual globalIndex numGlobalRows () const=0
 Get the number of global rows. More...
 
virtual globalIndex numGlobalCols () const=0
 Get the number of global columns. More...
 
Status query methods
bool closed () const
 Query matrix closed status. More...
 
bool assembled () const
 Query matrix assembled status. More...
 
bool ready () const
 Query matrix ready status. More...
 
bool modifiable () const
 Query matrix status. More...
 
bool insertable () const
 Query matrix status. More...
 
virtual bool created () const =0
 Query matrix creation status. More...
 
DofManager related methods

Some solvers and preconditioners rely on information about degrees-of-freedom of the linear problem represented by a matrix (for example, to decompose the matrix into blocks representing various parts of the physical problem). This lightweight interface allows one to associate a DofManager instance with the matrix, and thus avoid having to pass it to the preconditioner separately. The association is non-owning, the user must ensure the lifetime of DofManager object while any solvers/preconditioners set up with this matrix are in use.

void setDofManager (DofManager const *const dofManager)
 Associate a DofManager with this matrix. More...
 
const DofManagerdofManager () const
 
Create Methods
virtual void createWithLocalSize (localIndex const localSize, localIndex const maxEntriesPerRow, MPI_Comm const &comm)
 Create a square matrix from local number of rows. More...
 
virtual void createWithGlobalSize (globalIndex const globalSize, localIndex const maxEntriesPerRow, MPI_Comm const &comm)
 Create a square matrix from global number of rows. More...
 
virtual void createWithLocalSize (localIndex const localRows, localIndex const localCols, localIndex const maxEntriesPerRow, MPI_Comm const &comm)=0
 Create a rectangular matrix from number of rows/columns. More...
 
virtual void createWithGlobalSize (globalIndex const globalRows, globalIndex const globalCols, localIndex const maxEntriesPerRow, MPI_Comm const &comm)=0
 Create a rectangular matrix from number of rows/columns. More...
 
virtual void create (CRSMatrixView< real64 const, globalIndex const > const &localMatrix, localIndex const numLocalColumns, MPI_Comm const &comm)
 Create parallel matrix from a local CRS matrix. More...
 
Open/close methods
virtual void open ()=0
 Open matrix for adding new entries. More...
 
virtual void close ()=0
 Assemble and compress the matrix. More...
 
virtual void reset ()
 Reset the matrix to default state.
 
Global modification methods
virtual void set (real64 const value)=0
 Set all non-zero elements to a value. More...
 
virtual void zero ()=0
 Set all elements to zero.
 
Add/Set/Insert Methods

insert method can only be used while the matrix is being filled for the first time, e.g. during sparsity pattern construction. add and set methods can only be used after the matrix has been assembled, and both only accept entries that already exist in the sparsity pattern.

Note
Caution: these methods are not thread-safe.
virtual void add (globalIndex const rowIndex, globalIndex const colIndex, real64 const value)=0
 Add to one element. More...
 
virtual void set (globalIndex const rowIndex, globalIndex const colIndex, real64 const value)=0
 Set one element. More...
 
virtual void insert (globalIndex const rowIndex, globalIndex const colIndex, real64 const value)=0
 Insert one element. More...
 
virtual void add (globalIndex const rowIndex, globalIndex const *colIndices, real64 const *values, localIndex const size)=0
 Add elements to one row using c-style arrays. More...
 
virtual void set (globalIndex const rowIndex, globalIndex const *colIndices, real64 const *values, localIndex const size)=0
 Set elements to one row using c-style arrays. More...
 
virtual void insert (globalIndex const rowIndex, globalIndex const *colIndices, real64 const *values, localIndex const size)=0
 Insert elements to one row using c-style arrays. More...
 
virtual void add (globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values)=0
 Add elements to one row using array1d. More...
 
virtual void set (globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values)=0
 Set elements of one row using array1d. More...
 
virtual void insert (globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values)=0
 Insert elements of one row using array1d. More...
 
virtual void add (arraySlice1d< globalIndex const > const &rowIndices, arraySlice1d< globalIndex const > const &colIndices, arraySlice2d< real64 const > const &values)=0
 Add a dense block of values. More...
 
virtual void set (arraySlice1d< globalIndex const > const &rowIndices, arraySlice1d< globalIndex const > const &colIndices, arraySlice2d< real64 const > const &values)=0
 Set a dense block of values. More...
 
virtual void insert (arraySlice1d< globalIndex const > const &rowIndices, arraySlice1d< globalIndex const > const &colIndices, arraySlice2d< real64 const > const &values)=0
 Insert a dense block of values. More...
 
virtual void add (globalIndex const *rowIndices, globalIndex const *colIndices, real64 const *values, localIndex const numRows, localIndex const numCols)=0
 Add a dense block of values. More...
 
virtual void set (globalIndex const *rowIndices, globalIndex const *colIndices, real64 const *values, localIndex const numRows, localIndex const numCols)=0
 Set a dense block of values. More...
 
virtual void insert (globalIndex const *rowIndices, globalIndex const *colIndices, real64 const *values, localIndex const numRows, localIndex const numCols)=0
 Insert dense matrix. More...
 
virtual void insert (arrayView1d< globalIndex const > const &rowIndices, arrayView1d< globalIndex const > const &colIndices, arrayView1d< real64 const > const &values)=0
 Insert values stored in 3 linear vectors. More...
 
Linear Algebra Methods
virtual void apply (Vector const &src, Vector &dst) const override=0
 Apply operator to a vector, dst = this(src). More...
 
virtual void applyTranspose (Vector const &src, Vector &dst) const =0
 Apply transpose of the matrix to a vector. More...
 
virtual void residual (Vector const &x, Vector const &b, Vector &r) const override
 Compute residual r = b - A * x. More...
 
virtual void multiply (Matrix const &src, Matrix &dst) const =0
 Matrix/Matrix multiplication. More...
 
virtual void leftMultiplyTranspose (Matrix const &src, Matrix &dst) const =0
 Matrix/Matrix transpose multiplication. More...
 
virtual void rightMultiplyTranspose (Matrix const &src, Matrix &dst) const =0
 Matrix/Matrix transpose multiplication. More...
 
virtual void multiplyRAP (Matrix const &R, Matrix const &P, Matrix &dst) const
 Compute the triple product dst = R * this * P More...
 
virtual void multiplyPtAP (Matrix const &P, Matrix &dst) const
 Compute the triple product dst = P^T * this * P More...
 
virtual void multiplyRARt (Matrix const &R, Matrix &dst) const
 Compute the triple product dst = R * this * R^T More...
 
virtual void gemv (real64 const alpha, Vector const &x, real64 const beta, Vector &y, bool useTranspose=false) const =0
 Compute gemv y = alpha*A*x + beta*y. More...
 
virtual void separateComponentFilter (Matrix &dst, integer const dofPerPoint) const =0
 Apply a separate component approximation (filter) to this matrix. More...
 
virtual void scale (real64 const scalingFactor)=0
 Multiply all elements by scalingFactor. More...
 
virtual void leftScale (Vector const &vec)=0
 Pre-multiplies (left) with diagonal matrix consisting of the values in vec. More...
 
virtual void rightScale (Vector const &vec)=0
 Post-multiplies (right) with diagonal matrix consisting of the values in vec. More...
 
virtual void leftRightScale (Vector const &vecLeft, Vector const &vecRight)=0
 Post-multiplies (right) with diagonal matrix consisting of the values in vecRight and pre-multiplies (left) with diagonal matrix consisting of the values in vec. More...
 
virtual void rescaleRows (arrayView1d< globalIndex const > const &rowIndices, RowSumType const rowSumType)=0
 Rescales selected rows of matrix using row sum reciprocal as a factor. More...
 
virtual void transpose (Matrix &dst) const =0
 Matrix transposition. More...
 
virtual real64 clearRow (globalIndex const row, bool const keepDiag=false, real64 const diagValue=0.0)=0
 Clear a row, and optionally set diagonal element to diagValue. More...
 
virtual void addEntries (Matrix const &src, MatrixPatternOp const op, real64 const scale)=0
 Add entries of another matrix to this. More...
 
virtual void addDiagonal (Vector const &src, real64 const scale)=0
 Add (scaled) entries of a vector to the diagonal of this matrix. More...
 
virtual void clampEntries (real64 const lo, real64 const hi, bool const excludeDiag)=0
 Clamp each matrix value between values of lo and hi. More...
 
Accessors Methods
virtual localIndex maxRowLength () const =0
 Returns the number of nonzero entries in the longest row of the matrix. More...
 
virtual localIndex rowLength (globalIndex const globalRowIndex) const =0
 Get row length via global row index. More...
 
virtual void getRowLengths (arrayView1d< localIndex > const &lengths) const =0
 Get the row lengths of every local row. More...
 
virtual void getRowCopy (globalIndex const globalRow, arraySlice1d< globalIndex > const &colIndices, arraySlice1d< real64 > const &values) const =0
 Returns a copy of the data in row globalRow. More...
 
virtual void extractDiagonal (Vector &dst) const =0
 Extract diagonal values into a vector. More...
 
virtual void getRowSums (Vector &dst, RowSumType const rowSumType) const =0
 Populate a vector with row sums of this. More...
 
virtual globalIndex ilower () const =0
 Returns the index of the first global row owned by that processor. More...
 
virtual globalIndex iupper () const =0
 Returns index one past the last global row owned by that processor. More...
 
virtual globalIndex jlower () const =0
 Returns the index of the first global col owned by that processor. More...
 
virtual globalIndex jupper () const =0
 Returns index one past the last global col owned by that processor. More...
 
virtual localIndex numLocalNonzeros () const =0
 Returns the number of nonzeros in the local portion of the matrix. More...
 
virtual globalIndex numGlobalNonzeros () const =0
 Returns the total number of nonzeros in the matrix. More...
 
virtual real64 normInf () const =0
 Returns the infinity norm of the matrix. More...
 
virtual real64 norm1 () const =0
 Returns the one norm of the matrix. More...
 
virtual real64 normFrobenius () const =0
 Returns the Frobenius norm of the matrix. More...
 
virtual real64 normMax () const =0
 Returns the max norm of the matrix (the largest absolute element value). More...
 
virtual real64 normMax (arrayView1d< globalIndex const > const &rowIndices) const =0
 Returns the max norm of the matrix on a subset of rows. More...
 
virtual localIndex getLocalRowID (globalIndex const index) const =0
 Map a global row index to local row index. More...
 
virtual globalIndex getGlobalRowID (localIndex const index) const =0
 Map a local row index to global row index. More...
 
I/O Methods
virtual void print (std::ostream &os=std::cout) const =0
 Print the matrix in Trilinos format to a stream. More...
 
virtual void write (string const &filename, LAIOutputFormat const format=LAIOutputFormat::MATRIX_MARKET) const =0
 Write the matrix to filename in a matlab-compatible format. More...
 

Protected Attributes

bool m_closed = true
 Flag indicating whether the matrix is currently open for adding new entries.
 
bool m_assembled = false
 Flag indicating whether the matrix (sparsity pattern) has been assembled.
 
const DofManagerm_dofManager {}
 (optional) DofManager associated with this matrix
 

Friends

std::ostream & operator<< (std::ostream &os, Matrix const &matrix)
 Stream insertion operator for all matrix types. More...
 

Additional Inherited Members

- Public Types inherited from geos::LinearOperator< VECTOR >
using Vector = VECTOR
 Alias for template parameter.
 
- Public Member Functions inherited from geos::LinearOperator< VECTOR >
 LinearOperator ()=default
 Constructor.
 
virtual ~LinearOperator ()=default
 Destructor.
 
virtual globalIndex numGlobalRows () const =0
 Get the number of global rows. More...
 
virtual globalIndex numGlobalCols () const =0
 Get the number of global columns. More...
 
virtual localIndex numLocalRows () const =0
 Get the number of local rows. More...
 
virtual localIndex numLocalCols () const =0
 Get the number of local columns. More...
 
virtual MPI_Comm comm () const =0
 Get the MPI communicator the matrix was created with. More...
 

Detailed Description

template<typename MATRIX, typename VECTOR>
class geos::MatrixBase< MATRIX, VECTOR >

Common base template for all matrix wrapper types.

Template Parameters
MATRIXderived matrix type
VECTORcompatible vector type

This class template provides a common interface for all derived matrix wrapper types. Most methods are pure abstract in order to get the compiler to enforce a common interface in all derived classes; the class should be inherited from privately by implementations to avoid the possibility of accidentally invoking virtual function calls. Some basic facilities are provided related to matrix lifecycle. These include status flags and corresponding status query functions.

As an added benefit, the documentation for matrix interface also lives here and does not need to be duplicated across matrix wrappers. Derived classes should still document specific functions in case a particular LA package's behavior deviates from expectations or has unexpected performance impacts. In that case, @copydoc tag can be used to copy over the documentation.

Definition at line 74 of file MatrixBase.hpp.

Member Function Documentation

◆ add() [1/5]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::add ( arraySlice1d< globalIndex const > const &  rowIndices,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice2d< real64 const > const &  values 
)
protectedpure virtual

Add a dense block of values.

Parameters
rowIndicesGlobal row indices
colIndicesGlobal col indices
valuesDense local matrix of values

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ add() [2/5]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::add ( globalIndex const *  rowIndices,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  numRows,
localIndex const  numCols 
)
protectedpure virtual

Add a dense block of values.

Parameters
rowIndicesGlobal row indices
colIndicesGlobal col indices
valuesDense local matrix of values
numRowsNumber of row indices
numColsNumber of column indices
Note
Row major layout assumed in values

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ add() [3/5]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::add ( globalIndex const  rowIndex,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice1d< real64 const > const &  values 
)
protectedpure virtual

Add elements to one row using array1d.

Parameters
rowIndexGlobal row index
colIndicesGlobal column indices
valuesValues to add to prescribed locations

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ add() [4/5]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::add ( globalIndex const  rowIndex,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  size 
)
protectedpure virtual

Add elements to one row using c-style arrays.

Parameters
rowIndexGlobal row index
colIndicesGlobal column indices
valuesValues to add to prescribed locations
sizeNumber of elements

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ add() [5/5]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::add ( globalIndex const  rowIndex,
globalIndex const  colIndex,
real64 const  value 
)
protectedpure virtual

Add to one element.

Parameters
rowIndexGlobal row index
colIndexGlobal column index
valueValue to add to prescribed location

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ addDiagonal()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::addDiagonal ( Vector const &  src,
real64 const  scale 
)
protectedpure virtual

Add (scaled) entries of a vector to the diagonal of this matrix.

Parameters
srcthe source vector
scaleoptional scaling factor
Note
this must be square and have a (possibly zero) diagonal entry in every row. this and src must have the same parallel row distribution.

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ addEntries()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::addEntries ( Matrix const &  src,
MatrixPatternOp const  op,
real64 const  scale 
)
protectedpure virtual

Add entries of another matrix to this.

Parameters
srcthe source matrix
ophandling of nonzero patterns, see MatrixPatternOp
scalefactor to scale entries of src by
Note
Sparsity pattern of this must be a superset of sparsity of src. this and src must have the same parallel row distribution.

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ apply()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::apply ( Vector const &  src,
Vector dst 
) const
overrideprotectedpure virtual

Apply operator to a vector, dst = this(src).

Parameters
srcinput vector
dstoutput vector
Warning
src and dst cannot alias the same vector (some implementations may allow this).

Implements geos::LinearOperator< VECTOR >.

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ applyTranspose()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::applyTranspose ( Vector const &  src,
Vector dst 
) const
protectedpure virtual

Apply transpose of the matrix to a vector.

Parameters
srcInput vector (x).
dstOutput vector (b).

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ assembled()

template<typename MATRIX , typename VECTOR >
bool geos::MatrixBase< MATRIX, VECTOR >::assembled ( ) const
inlineprotected

Query matrix assembled status.

Returns
true if matrix has been opened and closed since creation; false otherwise

Definition at line 107 of file MatrixBase.hpp.

◆ clampEntries()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::clampEntries ( real64 const  lo,
real64 const  hi,
bool const  excludeDiag 
)
protectedpure virtual

Clamp each matrix value between values of lo and hi.

Parameters
lomin value
himax value
excludeDiagiff true, diagonal values are unchanged

Effectively sets each matrix value v to min(max(v, lo), hi).

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ clearRow()

template<typename MATRIX , typename VECTOR >
virtual real64 geos::MatrixBase< MATRIX, VECTOR >::clearRow ( globalIndex const  row,
bool const  keepDiag = false,
real64 const  diagValue = 0.0 
)
protectedpure virtual

Clear a row, and optionally set diagonal element to diagValue.

Parameters
rowglobalIndex of the row to be cleared.
diagValue(Optional) set diagonal element to desired value.
keepDiagif true, diagValue is ignored and original diagonal is preserved
Returns
original diagonal value if matrix is square; zero otherwise
Note
diagValue and keepDiag are ignored if the matrix is not square

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ close()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::close ( )
protectedpure virtual

Assemble and compress the matrix.

Compresses the matrix to CSR format with contiguous memory on each processor. Prevents from adding new entries in the sparsity pattern but allows for modification of existing entries.

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ closed()

template<typename MATRIX , typename VECTOR >
bool geos::MatrixBase< MATRIX, VECTOR >::closed ( ) const
inlineprotected

Query matrix closed status.

Returns
true if matrix has been opened and has not been closed since; false otherwise

Definition at line 101 of file MatrixBase.hpp.

◆ create()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::create ( CRSMatrixView< real64 const, globalIndex const > const &  localMatrix,
localIndex const  numLocalColumns,
MPI_Comm const &  comm 
)
inlineprotectedvirtual

Create parallel matrix from a local CRS matrix.

Parameters
localMatrixThe input local matrix.
numLocalColumnsnumber of local columns (not available from localMatrix in general)
commThe MPI communicator to use.
Note
Copies values, so that localMatrix does not need to retain its values after the call.
Todo:
Replace generic implementation with more efficient ones in each package.

Reimplemented in geos::HypreMatrix.

Definition at line 247 of file MatrixBase.hpp.

◆ created()

template<typename MATRIX , typename VECTOR >
virtual bool geos::MatrixBase< MATRIX, VECTOR >::created ( ) const
protectedpure virtual

Query matrix creation status.

Returns
true if matrix has been created

◆ createWithGlobalSize() [1/2]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::createWithGlobalSize ( globalIndex const  globalRows,
globalIndex const  globalCols,
localIndex const  maxEntriesPerRow,
MPI_Comm const &  comm 
)
protectedpure virtual

Create a rectangular matrix from number of rows/columns.

Parameters
commMPI communicator.
globalRowsGlobal number of rows.
globalColsGlobal number of columns.
maxEntriesPerRowMaximum number of entries per row (hint).

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ createWithGlobalSize() [2/2]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::createWithGlobalSize ( globalIndex const  globalSize,
localIndex const  maxEntriesPerRow,
MPI_Comm const &  comm 
)
inlineprotectedvirtual

Create a square matrix from global number of rows.

Create a square matrix with an (approximately) even partitioning of rows.

Parameters
globalSizeGlobal dimensions for a square matrix.
maxEntriesPerRowMaximum number of non-zero entries per row.
commMPI communicator.

Definition at line 203 of file MatrixBase.hpp.

◆ createWithLocalSize() [1/2]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::createWithLocalSize ( localIndex const  localRows,
localIndex const  localCols,
localIndex const  maxEntriesPerRow,
MPI_Comm const &  comm 
)
protectedpure virtual

Create a rectangular matrix from number of rows/columns.

Parameters
commMPI communicator.
localRowsLocal number of rows.
localColsLocal number of columns.
maxEntriesPerRowMaximum number of entries per row (hint).

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ createWithLocalSize() [2/2]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::createWithLocalSize ( localIndex const  localSize,
localIndex const  maxEntriesPerRow,
MPI_Comm const &  comm 
)
inlineprotectedvirtual

Create a square matrix from local number of rows.

Parameters
localSizelocal number of rows for square matrix.
maxEntriesPerRowMaximum number of non-zero entries per row.
commMPI communicator.

Definition at line 185 of file MatrixBase.hpp.

◆ dofManager()

template<typename MATRIX , typename VECTOR >
const DofManager* geos::MatrixBase< MATRIX, VECTOR >::dofManager ( ) const
inlineprotected
Returns
the associated DofManager

Definition at line 164 of file MatrixBase.hpp.

◆ extractDiagonal()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::extractDiagonal ( Vector dst) const
protectedpure virtual

Extract diagonal values into a vector.

Parameters
dstthe target vector, must have the same row partitioning as this

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ gemv()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::gemv ( real64 const  alpha,
Vector const &  x,
real64 const  beta,
Vector y,
bool  useTranspose = false 
) const
protectedpure virtual

Compute gemv y = alpha*A*x + beta*y.

Note
The naming convention follows the BLAS library.
Parameters
alphaScalar factor for added matvec product.
xInput vector.
betaScalar factor for right hand side.
yOutput vector.
useTransposeBoolean, set to true to use A^T.
Warning
x and y cannot alias the same vector.

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ getGlobalRowID()

template<typename MATRIX , typename VECTOR >
virtual globalIndex geos::MatrixBase< MATRIX, VECTOR >::getGlobalRowID ( localIndex const  index) const
protectedpure virtual

Map a local row index to global row index.

Parameters
indexthe local row index (between 0 and number of local rows)
Returns
the global row index corresponding to index

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ getLocalRowID()

template<typename MATRIX , typename VECTOR >
virtual localIndex geos::MatrixBase< MATRIX, VECTOR >::getLocalRowID ( globalIndex const  index) const
protectedpure virtual

Map a global row index to local row index.

Parameters
indexthe global row index
Returns
the local row index corresponding to index, or -1 if not a local row

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ getRowCopy()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::getRowCopy ( globalIndex const  globalRow,
arraySlice1d< globalIndex > const &  colIndices,
arraySlice1d< real64 > const &  values 
) const
protectedpure virtual

Returns a copy of the data in row globalRow.

Parameters
[in]globalRowthe index of global row to extract
[out]colIndicesthe output array of global column indices (must have a large enough size)
[out]valuesthe output array of values (must have a large enough size)

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ getRowLengths()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::getRowLengths ( arrayView1d< localIndex > const &  lengths) const
protectedpure virtual

Get the row lengths of every local row.

Parameters
lengthsan array view to be populated with row lengths
Note
The implementation may move the view's buffer to a different memory space.

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ getRowSums()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::getRowSums ( Vector dst,
RowSumType const  rowSumType 
) const
protectedpure virtual

Populate a vector with row sums of this.

Parameters
dstthe target vector, must have the same row partitioning as this
rowSumTypetype of row sum operation to perform

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ ilower()

template<typename MATRIX , typename VECTOR >
virtual globalIndex geos::MatrixBase< MATRIX, VECTOR >::ilower ( ) const
protectedpure virtual

Returns the index of the first global row owned by that processor.

Returns
the index of the first global row owned by that processor

◆ insert() [1/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::insert ( arraySlice1d< globalIndex const > const &  rowIndices,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice2d< real64 const > const &  values 
)
protectedpure virtual

Insert a dense block of values.

Parameters
rowIndicesGlobal row indices
colIndicesGlobal col indices
valuesDense local matrix of values

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ insert() [2/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::insert ( arrayView1d< globalIndex const > const &  rowIndices,
arrayView1d< globalIndex const > const &  colIndices,
arrayView1d< real64 const > const &  values 
)
protectedpure virtual

Insert values stored in 3 linear vectors.

Parameters
rowIndicesArray of global row indices
colIndicesArray of global column indices
valuesArray of values

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ insert() [3/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::insert ( globalIndex const *  rowIndices,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  numRows,
localIndex const  numCols 
)
protectedpure virtual

Insert dense matrix.

Parameters
rowIndicesGlobal row indices
colIndicesGlobal col indices
valuesDense local matrix of values
numRowsNumber of row indices
numColsNumber of column indices
Note
Row major layout assumed in values

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ insert() [4/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::insert ( globalIndex const  rowIndex,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice1d< real64 const > const &  values 
)
protectedpure virtual

Insert elements of one row using array1d.

Parameters
rowIndexGlobal row index
colIndicesGlobal column indices
valuesValues to add to prescribed locations

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ insert() [5/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::insert ( globalIndex const  rowIndex,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  size 
)
protectedpure virtual

Insert elements to one row using c-style arrays.

Parameters
rowIndexGlobal row index
colIndicesGlobal column indices
valuesValues to add to prescribed locations
sizeNumber of elements

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ insert() [6/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::insert ( globalIndex const  rowIndex,
globalIndex const  colIndex,
real64 const  value 
)
protectedpure virtual

Insert one element.

Parameters
rowIndexGlobal row index
colIndexGlobal column index
valueValue to insert at prescribed location

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ insertable()

template<typename MATRIX , typename VECTOR >
bool geos::MatrixBase< MATRIX, VECTOR >::insertable ( ) const
inlineprotected

Query matrix status.

Returns
true if matrix has NOT been assembled yet (not closed since last create() call) and is currently open for insertion of new entries

Definition at line 129 of file MatrixBase.hpp.

◆ iupper()

template<typename MATRIX , typename VECTOR >
virtual globalIndex geos::MatrixBase< MATRIX, VECTOR >::iupper ( ) const
protectedpure virtual

Returns index one past the last global row owned by that processor.

Returns
the next index after last global row owned by that processor
Note
The intention is for [ilower; iupper) to be used as a half-open index range

◆ jlower()

template<typename MATRIX , typename VECTOR >
virtual globalIndex geos::MatrixBase< MATRIX, VECTOR >::jlower ( ) const
protectedpure virtual

Returns the index of the first global col owned by that processor.

Returns
index of the first owned global col
Note
Matrix implementations don't physically "own" column ranges the same way they do row ranges. Instead, the column range refers to the "diagonal" block of columns which would correspond to the local range of entries of a vector created with the same local/global size as the number of matrix columns.

◆ jupper()

template<typename MATRIX , typename VECTOR >
virtual globalIndex geos::MatrixBase< MATRIX, VECTOR >::jupper ( ) const
protectedpure virtual

Returns index one past the last global col owned by that processor.

Returns
index one past the last owned global col
Note
The intention is for [jlower; jupper) to be used as a half-open index range.
Also see note for jlower() about the meaning of "owned" columns.

◆ leftMultiplyTranspose()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::leftMultiplyTranspose ( Matrix const &  src,
Matrix dst 
) const
protectedpure virtual

Matrix/Matrix transpose multiplication.

Compute this^T * B = C.

Parameters
srcInput matrix (B).
dstOutput matrix (C).
Note
The output matrix dst doesn't need to be created beforehand.

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ leftRightScale()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::leftRightScale ( Vector const &  vecLeft,
Vector const &  vecRight 
)
protectedpure virtual

Post-multiplies (right) with diagonal matrix consisting of the values in vecRight and pre-multiplies (left) with diagonal matrix consisting of the values in vec.

Parameters
vecLeftvec to pre-multiply with.
vecRightvec to post-multiply with.

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ leftScale()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::leftScale ( Vector const &  vec)
protectedpure virtual

Pre-multiplies (left) with diagonal matrix consisting of the values in vec.

Parameters
vecVector to pre-multiply with.

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ maxRowLength()

template<typename MATRIX , typename VECTOR >
virtual localIndex geos::MatrixBase< MATRIX, VECTOR >::maxRowLength ( ) const
protectedpure virtual

Returns the number of nonzero entries in the longest row of the matrix.

Returns
the max length of a row

Collective.

◆ modifiable()

template<typename MATRIX , typename VECTOR >
bool geos::MatrixBase< MATRIX, VECTOR >::modifiable ( ) const
inlineprotected

Query matrix status.

Returns
true if matrix has been assembled and is currently open; this implies individual entries within existing sparsity pattern can be altered via set()/add() methods.

Definition at line 122 of file MatrixBase.hpp.

◆ multiply()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::multiply ( Matrix const &  src,
Matrix dst 
) const
protectedpure virtual

Matrix/Matrix multiplication.

Compute this * B = C.

Parameters
srcInput matrix (B).
dstOutput matrix (C).
Note
The output matrix dst doesn't need to be created beforehand.

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ multiplyPtAP()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::multiplyPtAP ( Matrix const &  P,
Matrix dst 
) const
inlineprotectedvirtual

Compute the triple product dst = P^T * this * P

Parameters
Pthe "prolongation" matrix
dstthe resulting product matrix (will be re-created as needed)

Reimplemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

Definition at line 624 of file MatrixBase.hpp.

◆ multiplyRAP()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::multiplyRAP ( Matrix const &  R,
Matrix const &  P,
Matrix dst 
) const
inlineprotectedvirtual

Compute the triple product dst = R * this * P

Parameters
Rthe "restriction" matrix
Pthe "prolongation" matrix
dstthe resulting product matrix (will be re-created as needed)

Reimplemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

Definition at line 604 of file MatrixBase.hpp.

◆ multiplyRARt()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::multiplyRARt ( Matrix const &  R,
Matrix dst 
) const
inlineprotectedvirtual

Compute the triple product dst = R * this * R^T

Parameters
Rthe "restriction" matrix
dstthe resulting product matrix (will be re-created as needed)

Definition at line 642 of file MatrixBase.hpp.

◆ norm1()

template<typename MATRIX , typename VECTOR >
virtual real64 geos::MatrixBase< MATRIX, VECTOR >::norm1 ( ) const
protectedpure virtual

Returns the one norm of the matrix.

Returns
the value of 1-norm

◆ normFrobenius()

template<typename MATRIX , typename VECTOR >
virtual real64 geos::MatrixBase< MATRIX, VECTOR >::normFrobenius ( ) const
protectedpure virtual

Returns the Frobenius norm of the matrix.

Returns
the value of Frobenius norm

◆ normInf()

template<typename MATRIX , typename VECTOR >
virtual real64 geos::MatrixBase< MATRIX, VECTOR >::normInf ( ) const
protectedpure virtual

Returns the infinity norm of the matrix.

Returns
the value of infinity norm

◆ normMax() [1/2]

template<typename MATRIX , typename VECTOR >
virtual real64 geos::MatrixBase< MATRIX, VECTOR >::normMax ( ) const
protectedpure virtual

Returns the max norm of the matrix (the largest absolute element value).

Returns
the value of max norm

◆ normMax() [2/2]

template<typename MATRIX , typename VECTOR >
virtual real64 geos::MatrixBase< MATRIX, VECTOR >::normMax ( arrayView1d< globalIndex const > const &  rowIndices) const
protectedpure virtual

Returns the max norm of the matrix on a subset of rows.

Parameters
rowIndicesglobal indices of rows to compute norm over
Returns
the value of max norm

Implemented in geos::HypreMatrix, geos::EpetraMatrix, and geos::PetscMatrix.

◆ numGlobalCols()

template<typename MATRIX , typename VECTOR >
virtual globalIndex geos::LinearOperator< VECTOR >::numGlobalCols
protected

Get the number of global columns.

Returns
Number of global columns in the operator.

◆ numGlobalNonzeros()

template<typename MATRIX , typename VECTOR >
virtual globalIndex geos::MatrixBase< MATRIX, VECTOR >::numGlobalNonzeros ( ) const
protectedpure virtual

Returns the total number of nonzeros in the matrix.

Returns
the total number of nonzeros in the matrix

◆ numGlobalRows()

template<typename MATRIX , typename VECTOR >
virtual globalIndex geos::LinearOperator< VECTOR >::numGlobalRows
protected

Get the number of global rows.

Returns
Number of global rows in the operator.

◆ numLocalCols()

template<typename MATRIX , typename VECTOR >
virtual localIndex geos::LinearOperator< VECTOR >::numLocalCols
protected

Get the number of local columns.

Returns
Number of local columns in the operator.
Note
The use of term "local columns" refers not to physical partitioning of columns across ranks (as e.g. matrices are partitioned by rows and typically physically store all column entries), but to the partitioning of a compatible vector object that this operator can be applied to.

◆ numLocalNonzeros()

template<typename MATRIX , typename VECTOR >
virtual localIndex geos::MatrixBase< MATRIX, VECTOR >::numLocalNonzeros ( ) const
protectedpure virtual

Returns the number of nonzeros in the local portion of the matrix.

Returns
the number of nonzeros in the local portion of the matrix

◆ numLocalRows()

template<typename MATRIX , typename VECTOR >
virtual localIndex geos::LinearOperator< VECTOR >::numLocalRows
protected

Get the number of local rows.

Returns
Number of local rows in the operator.

◆ open()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::open ( )
protectedpure virtual

Open matrix for adding new entries.

Note
Adding entries that result in modifications of sparsity pattern may not be allowed by most implementations. An error will be raised in that case.

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ print()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::print ( std::ostream &  os = std::cout) const
protectedpure virtual

Print the matrix in Trilinos format to a stream.

Parameters
osthe output stream

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ ready()

template<typename MATRIX , typename VECTOR >
bool geos::MatrixBase< MATRIX, VECTOR >::ready ( ) const
inlineprotected

Query matrix ready status.

Returns
true if matrix has been assembled and is currently closed; this implies it's ready to be used or re-opened for adding/setting values

Definition at line 114 of file MatrixBase.hpp.

◆ rescaleRows()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::rescaleRows ( arrayView1d< globalIndex const > const &  rowIndices,
RowSumType const  rowSumType 
)
protectedpure virtual

Rescales selected rows of matrix using row sum reciprocal as a factor.

Parameters
rowIndicesglobal indicies of rows to scale (all must be locally owned)
rowSumTypetype of row sums to use as scaling factors

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ residual()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::residual ( Vector const &  x,
Vector const &  b,
Vector r 
) const
inlineoverrideprotectedvirtual

Compute residual r = b - A * x.

Overrides LinearOperator::residual().

Parameters
xInput solution.
bInput right hand side.
rOutput residual.
Warning
x and r cannot alias the same vector, but b and r can.

Reimplemented from geos::LinearOperator< VECTOR >.

Definition at line 550 of file MatrixBase.hpp.

◆ rightMultiplyTranspose()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::rightMultiplyTranspose ( Matrix const &  src,
Matrix dst 
) const
protectedpure virtual

Matrix/Matrix transpose multiplication.

Compute B * this^T = C.

Parameters
srcInput matrix (B).
dstOutput matrix (C).
Note
The output matrix dst doesn't need to be created beforehand.

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ rightScale()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::rightScale ( Vector const &  vec)
protectedpure virtual

Post-multiplies (right) with diagonal matrix consisting of the values in vec.

Parameters
vecVector to post-multiply with.

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ rowLength()

template<typename MATRIX , typename VECTOR >
virtual localIndex geos::MatrixBase< MATRIX, VECTOR >::rowLength ( globalIndex const  globalRowIndex) const
protectedpure virtual

Get row length via global row index.

Parameters
[in]globalRowIndexthe global row index
Returns
the number of nonzero entries in the row

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ scale()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::scale ( real64 const  scalingFactor)
protectedpure virtual

Multiply all elements by scalingFactor.

Parameters
scalingFactorScaling factor.

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ separateComponentFilter()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::separateComponentFilter ( Matrix dst,
integer const  dofPerPoint 
) const
protectedpure virtual

Apply a separate component approximation (filter) to this matrix.

Parameters
dstthe target (filtered) matrix
dofPerPointnumber of degrees-of-freedom per node

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ set() [1/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::set ( arraySlice1d< globalIndex const > const &  rowIndices,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice2d< real64 const > const &  values 
)
protectedpure virtual

Set a dense block of values.

Parameters
rowIndicesGlobal row indices
colIndicesGlobal col indices
valuesDense local matrix of values

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ set() [2/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::set ( globalIndex const *  rowIndices,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  numRows,
localIndex const  numCols 
)
protectedpure virtual

Set a dense block of values.

Parameters
rowIndicesGlobal row indices
colIndicesGlobal col indices
valuesDense local matrix of values
numRowsNumber of row indices
numColsNumber of column indices
Note
Row major layout assumed in values

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ set() [3/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::set ( globalIndex const  rowIndex,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice1d< real64 const > const &  values 
)
protectedpure virtual

Set elements of one row using array1d.

Parameters
rowIndexGlobal row index
colIndicesGlobal column indices
valuesValues to add to prescribed locations

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ set() [4/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::set ( globalIndex const  rowIndex,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  size 
)
protectedpure virtual

Set elements to one row using c-style arrays.

Parameters
rowIndexGlobal row index
colIndicesGlobal column indices
valuesValues to add to prescribed locations
sizeNumber of elements

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ set() [5/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::set ( globalIndex const  rowIndex,
globalIndex const  colIndex,
real64 const  value 
)
protectedpure virtual

Set one element.

Parameters
rowIndexGlobal row index
colIndexGlobal column index
valueValue to set at prescribed location

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ set() [6/6]

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::set ( real64 const  value)
protectedpure virtual

Set all non-zero elements to a value.

Parameters
valuethe value to set all elements to

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

◆ setDofManager()

template<typename MATRIX , typename VECTOR >
void geos::MatrixBase< MATRIX, VECTOR >::setDofManager ( DofManager const *const  dofManager)
inlineprotected

Associate a DofManager with this matrix.

Parameters
dofManagerthe DofManager containing the relevant degrees of freedom

Definition at line 156 of file MatrixBase.hpp.

◆ transpose()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::transpose ( Matrix dst) const
protectedpure virtual

Matrix transposition.

Compute B = this^T.

Parameters
dstOutput matrix (B).

Implemented in geos::PetscMatrix, geos::HypreMatrix, and geos::EpetraMatrix.

◆ write()

template<typename MATRIX , typename VECTOR >
virtual void geos::MatrixBase< MATRIX, VECTOR >::write ( string const &  filename,
LAIOutputFormat const  format = LAIOutputFormat::MATRIX_MARKET 
) const
protectedpure virtual

Write the matrix to filename in a matlab-compatible format.

Parameters
filenamename of the output file
formatoutput format

Within octave / matlab:

load filename M = spconvert(filename_root)

Implemented in geos::HypreMatrix, geos::PetscMatrix, and geos::EpetraMatrix.

Friends And Related Function Documentation

◆ operator<<

template<typename MATRIX , typename VECTOR >
std::ostream& operator<< ( std::ostream &  os,
Matrix const &  matrix 
)
friend

Stream insertion operator for all matrix types.

Parameters
osthe output stream
matrixthe matrix to be printed
Returns
reference to the output stream

Definition at line 949 of file MatrixBase.hpp.


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