GEOSX
Public Types | Public Member Functions | List of all members
geos::EpetraMatrix Class Referencefinal

Wrapper class for Epetra's CrsMatrix. More...

#include <EpetraMatrix.hpp>

Inheritance diagram for geos::EpetraMatrix:
Inheritance graph
[legend]

Public Types

using Vector = EpetraVector
 Compatible vector type.
 
using Export = EpetraExport
 Associated exporter type.
 
- Public Types inherited from geos::LinearOperator< EpetraVector >
using Vector = EpetraVector
 Alias for template parameter.
 

Public Member Functions

const Epetra_FECrsMatrix & unwrapped () const
 Returns a const pointer to the underlying matrix. More...
 
Epetra_FECrsMatrix & unwrapped ()
 Returns a non-const pointer to the underlying matrix. More...
 
Constructor/Destructor methods
 EpetraMatrix ()
 Empty matrix constructor. More...
 
 EpetraMatrix (EpetraMatrix const &src)
 Copy constructor. More...
 
 EpetraMatrix (EpetraMatrix &&src) noexcept
 Move constructor. More...
 
EpetraMatrixoperator= (EpetraMatrix const &src)
 Copy assignment. More...
 
EpetraMatrixoperator= (EpetraMatrix &&src) noexcept
 Move assignment. More...
 
virtual ~EpetraMatrix () override
 Virtual destructor.
 
MatrixBase interface
virtual void createWithLocalSize (localIndex const localRows, localIndex const localCols, localIndex const maxEntriesPerRow, MPI_Comm const &comm) override
 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) override
 Create a rectangular matrix from number of rows/columns. More...
 
virtual void open () override
 Open matrix for adding new entries. More...
 
virtual void close () override
 Assemble and compress the matrix. More...
 
virtual bool created () const override
 Query matrix creation status. More...
 
virtual void reset () override
 Reset the matrix to default state.
 
virtual void set (real64 const value) override
 Set all non-zero elements to a value. More...
 
virtual void zero () override
 Set all elements to zero.
 
virtual void add (globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
 Add to one element. More...
 
virtual void set (globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
 Set one element. More...
 
virtual void insert (globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
 Insert one element. More...
 
virtual void add (globalIndex const rowIndex, globalIndex const *colIndices, real64 const *values, localIndex const size) override
 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) override
 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) override
 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) override
 Add elements to one row using array1d. More...
 
virtual void set (globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values) override
 Set elements of one row using array1d. More...
 
virtual void insert (globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values) override
 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) override
 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) override
 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) override
 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) override
 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) override
 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) override
 Insert dense matrix. More...
 
virtual void insert (arrayView1d< globalIndex const > const &rowIndices, arrayView1d< globalIndex const > const &colIndices, arrayView1d< real64 const > const &values) override
 Insert values stored in 3 linear vectors. More...
 
virtual void apply (EpetraVector const &src, EpetraVector &dst) const override
 Apply operator to a vector, dst = this(src). More...
 
virtual void applyTranspose (EpetraVector const &src, EpetraVector &dst) const override
 Apply transpose of the matrix to a vector. More...
 
virtual void multiply (EpetraMatrix const &src, EpetraMatrix &dst) const override
 Matrix/Matrix multiplication. More...
 
virtual void leftMultiplyTranspose (EpetraMatrix const &src, EpetraMatrix &dst) const override
 Matrix/Matrix transpose multiplication. More...
 
virtual void rightMultiplyTranspose (EpetraMatrix const &src, EpetraMatrix &dst) const override
 Matrix/Matrix transpose multiplication. More...
 
virtual void multiplyRAP (EpetraMatrix const &R, EpetraMatrix const &P, EpetraMatrix &dst) const override
 Compute the triple product dst = R * this * P More...
 
virtual void multiplyPtAP (EpetraMatrix const &P, EpetraMatrix &dst) const override
 Compute the triple product dst = P^T * this * P More...
 
virtual void gemv (real64 const alpha, EpetraVector const &x, real64 const beta, EpetraVector &y, bool useTranspose=false) const override
 Compute gemv y = alpha*A*x + beta*y. More...
 
virtual void scale (real64 const scalingFactor) override
 Multiply all elements by scalingFactor. More...
 
virtual void leftScale (EpetraVector const &vec) override
 Pre-multiplies (left) with diagonal matrix consisting of the values in vec. More...
 
virtual void rightScale (EpetraVector const &vec) override
 Post-multiplies (right) with diagonal matrix consisting of the values in vec. More...
 
virtual void leftRightScale (EpetraVector const &vecLeft, EpetraVector const &vecRight) override
 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) override
 Rescales selected rows of matrix using row sum reciprocal as a factor. More...
 
virtual void transpose (EpetraMatrix &dst) const override
 Matrix transposition. More...
 
virtual void separateComponentFilter (EpetraMatrix &dst, integer const dofsPerNode) const override
 Apply a separate component approximation (filter) to this matrix. More...
 
virtual real64 clearRow (globalIndex const row, bool const keepDiag=false, real64 const diagValue=0.0) override
 Clear a row, and optionally set diagonal element to diagValue. More...
 
virtual void addEntries (EpetraMatrix const &src, MatrixPatternOp const op, real64 const scale) override
 Add entries of another matrix to this. More...
 
virtual void addDiagonal (EpetraVector const &src, real64 const scale) override
 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) override
 Clamp each matrix value between values of lo and hi. More...
 
virtual localIndex maxRowLength () const override
 Returns the number of nonzero entries in the longest row of the matrix. More...
 
virtual localIndex rowLength (globalIndex const globalRowIndex) const override
 Get row length via global row index. More...
 
virtual void getRowLengths (arrayView1d< localIndex > const &lengths) const override
 Get the row lengths of every local row. More...
 
virtual void getRowCopy (globalIndex globalRow, arraySlice1d< globalIndex > const &colIndices, arraySlice1d< real64 > const &values) const override
 Returns a copy of the data in row globalRow. More...
 
virtual void extractDiagonal (EpetraVector &dst) const override
 Extract diagonal values into a vector. More...
 
virtual void getRowSums (EpetraVector &dst, RowSumType const rowSumType) const override
 Populate a vector with row sums of this. More...
 
virtual globalIndex numGlobalRows () const override
 Get the number of global rows. More...
 
virtual globalIndex numGlobalCols () const override
 Get the number of global columns. More...
 
virtual localIndex numLocalRows () const override
 Get the number of local rows. More...
 
virtual localIndex numLocalCols () const override
 Get the number of local columns. More...
 
virtual globalIndex ilower () const override
 Returns the index of the first global row owned by that processor. More...
 
virtual globalIndex iupper () const override
 Returns index one past the last global row owned by that processor. More...
 
virtual globalIndex jlower () const override
 Returns the index of the first global col owned by that processor. More...
 
virtual globalIndex jupper () const override
 Returns index one past the last global col owned by that processor. More...
 
virtual localIndex numLocalNonzeros () const override
 Returns the number of nonzeros in the local portion of the matrix. More...
 
virtual globalIndex numGlobalNonzeros () const override
 Returns the total number of nonzeros in the matrix. More...
 
virtual real64 normInf () const override
 Returns the infinity norm of the matrix. More...
 
virtual real64 norm1 () const override
 Returns the one norm of the matrix. More...
 
virtual real64 normFrobenius () const override
 Returns the Frobenius norm of the matrix. More...
 
virtual real64 normMax () const override
 Returns the max norm of the matrix (the largest absolute element value). More...
 
virtual real64 normMax (arrayView1d< globalIndex const > const &rowIndices) const override
 Returns the max norm of the matrix on a subset of rows. More...
 
virtual localIndex getLocalRowID (globalIndex const index) const override
 Map a global row index to local row index. More...
 
virtual globalIndex getGlobalRowID (localIndex const index) const override
 Map a local row index to global row index. More...
 
virtual MPI_Comm comm () const override
 Get the MPI communicator the matrix was created with. More...
 
virtual void print (std::ostream &os=std::cout) const override
 Print the matrix in Trilinos format to a stream. More...
 
virtual void write (string const &filename, LAIOutputFormat const format=LAIOutputFormat::MATRIX_MARKET) const override
 Write the matrix to filename in a matlab-compatible format. More...
 
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 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 globalSize, localIndex const maxEntriesPerRow, MPI_Comm const &comm)
 Create a square matrix from global number of rows. 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...
 
bool closed () const
 Query matrix closed status. More...
 
bool assembled () const
 Query matrix assembled status. More...
 
bool insertable () const
 Query matrix status. More...
 
bool modifiable () const
 Query matrix status. More...
 
bool ready () const
 Query matrix ready status. More...
 
virtual void residual (Vector const &x, Vector const &b, Vector &r) const override
 Compute residual r = b - A * x. More...
 
void setDofManager (DofManager const *const dofManager)
 Associate a DofManager with this matrix. More...
 
const DofManagerdofManager () const
 
- Public Member Functions inherited from geos::LinearOperator< EpetraVector >
 LinearOperator ()=default
 Constructor.
 
virtual ~LinearOperator ()=default
 Destructor.
 
virtual MPI_Comm comm () const=0
 Get the MPI communicator the matrix was created with. More...
 

Detailed Description

Wrapper class for Epetra's CrsMatrix.

Definition at line 38 of file EpetraMatrix.hpp.

Constructor & Destructor Documentation

◆ EpetraMatrix() [1/3]

geos::EpetraMatrix::EpetraMatrix ( )

Empty matrix constructor.

Create an empty (distributed) matrix.

◆ EpetraMatrix() [2/3]

geos::EpetraMatrix::EpetraMatrix ( EpetraMatrix const &  src)

Copy constructor.

Parameters
[in]srcthe matrix to be copied

◆ EpetraMatrix() [3/3]

geos::EpetraMatrix::EpetraMatrix ( EpetraMatrix &&  src)
noexcept

Move constructor.

Parameters
[in]srcthe matrix to be copied

Member Function Documentation

◆ add() [1/5]

virtual void geos::EpetraMatrix::add ( arraySlice1d< globalIndex const > const &  rowIndices,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice2d< real64 const > const &  values 
)
overridevirtual

Add a dense block of values.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ add() [2/5]

virtual void geos::EpetraMatrix::add ( globalIndex const *  rowIndices,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  numRows,
localIndex const  numCols 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ add() [3/5]

virtual void geos::EpetraMatrix::add ( globalIndex const  rowIndex,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice1d< real64 const > const &  values 
)
overridevirtual

Add elements to one row using array1d.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ add() [4/5]

virtual void geos::EpetraMatrix::add ( globalIndex const  rowIndex,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  size 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ add() [5/5]

virtual void geos::EpetraMatrix::add ( globalIndex const  rowIndex,
globalIndex const  colIndex,
real64 const  value 
)
overridevirtual

Add to one element.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ addDiagonal()

virtual void geos::EpetraMatrix::addDiagonal ( EpetraVector const &  src,
real64 const  scale 
)
overridevirtual

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.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ addEntries()

virtual void geos::EpetraMatrix::addEntries ( EpetraMatrix const &  src,
MatrixPatternOp const  op,
real64 const  scale 
)
overridevirtual

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.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ apply()

virtual void geos::EpetraMatrix::apply ( EpetraVector const &  src,
EpetraVector dst 
) const
overridevirtual

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< EpetraVector >.

◆ applyTranspose()

virtual void geos::EpetraMatrix::applyTranspose ( EpetraVector const &  src,
EpetraVector dst 
) const
overridevirtual

Apply transpose of the matrix to a vector.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ assembled()

bool geos::MatrixBase< MATRIX, VECTOR >::assembled
inline

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

virtual void geos::EpetraMatrix::clampEntries ( real64 const  lo,
real64 const  hi,
bool const  excludeDiag 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ clearRow()

virtual real64 geos::EpetraMatrix::clearRow ( globalIndex const  row,
bool const  keepDiag = false,
real64 const  diagValue = 0.0 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ close()

virtual void geos::EpetraMatrix::close ( )
overridevirtual

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.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ closed()

bool geos::MatrixBase< MATRIX, VECTOR >::closed
inline

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.

◆ comm()

virtual MPI_Comm geos::EpetraMatrix::comm ( ) const
overridevirtual

Get the MPI communicator the matrix was created with.

Returns
MPI communicator passed in create...()
Note
when build without MPI, may return anything (MPI_Comm will be a mock type defined in MpiWrapper)

◆ create()

virtual void geos::MatrixBase< MATRIX, VECTOR >::create
inline

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.

Definition at line 247 of file MatrixBase.hpp.

◆ created()

virtual bool geos::EpetraMatrix::created ( ) const
overridevirtual

Query matrix creation status.

Returns
true if matrix has been created

◆ createWithGlobalSize() [1/3]

virtual void geos::EpetraMatrix::createWithGlobalSize ( globalIndex const  globalRows,
globalIndex const  globalCols,
localIndex const  maxEntriesPerRow,
MPI_Comm const &  comm 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ createWithGlobalSize() [2/3]

virtual void geos::MatrixBase< MATRIX, VECTOR >::createWithGlobalSize

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

◆ createWithGlobalSize() [3/3]

virtual void geos::MatrixBase< MATRIX, VECTOR >::createWithGlobalSize
inline

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/3]

virtual void geos::EpetraMatrix::createWithLocalSize ( localIndex const  localRows,
localIndex const  localCols,
localIndex const  maxEntriesPerRow,
MPI_Comm const &  comm 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ createWithLocalSize() [2/3]

virtual void geos::MatrixBase< MATRIX, VECTOR >::createWithLocalSize

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

◆ createWithLocalSize() [3/3]

virtual void geos::MatrixBase< MATRIX, VECTOR >::createWithLocalSize
inline

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

const DofManager* geos::MatrixBase< MATRIX, VECTOR >::dofManager
inline
Returns
the associated DofManager

Definition at line 164 of file MatrixBase.hpp.

◆ extractDiagonal()

virtual void geos::EpetraMatrix::extractDiagonal ( EpetraVector dst) const
overridevirtual

Extract diagonal values into a vector.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ gemv()

virtual void geos::EpetraMatrix::gemv ( real64 const  alpha,
EpetraVector const &  x,
real64 const  beta,
EpetraVector y,
bool  useTranspose = false 
) const
overridevirtual

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.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ getGlobalRowID()

virtual globalIndex geos::EpetraMatrix::getGlobalRowID ( localIndex const  index) const
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ getLocalRowID()

virtual localIndex geos::EpetraMatrix::getLocalRowID ( globalIndex const  index) const
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ getRowCopy()

virtual void geos::EpetraMatrix::getRowCopy ( globalIndex  globalRow,
arraySlice1d< globalIndex > const &  colIndices,
arraySlice1d< real64 > const &  values 
) const
overridevirtual

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)

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ getRowLengths()

virtual void geos::EpetraMatrix::getRowLengths ( arrayView1d< localIndex > const &  lengths) const
overridevirtual

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.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ getRowSums()

virtual void geos::EpetraMatrix::getRowSums ( EpetraVector dst,
RowSumType const  rowSumType 
) const
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ ilower()

virtual globalIndex geos::EpetraMatrix::ilower ( ) const
overridevirtual

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]

virtual void geos::EpetraMatrix::insert ( arraySlice1d< globalIndex const > const &  rowIndices,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice2d< real64 const > const &  values 
)
overridevirtual

Insert a dense block of values.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ insert() [2/6]

virtual void geos::EpetraMatrix::insert ( arrayView1d< globalIndex const > const &  rowIndices,
arrayView1d< globalIndex const > const &  colIndices,
arrayView1d< real64 const > const &  values 
)
overridevirtual

Insert values stored in 3 linear vectors.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ insert() [3/6]

virtual void geos::EpetraMatrix::insert ( globalIndex const *  rowIndices,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  numRows,
localIndex const  numCols 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ insert() [4/6]

virtual void geos::EpetraMatrix::insert ( globalIndex const  rowIndex,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice1d< real64 const > const &  values 
)
overridevirtual

Insert elements of one row using array1d.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ insert() [5/6]

virtual void geos::EpetraMatrix::insert ( globalIndex const  rowIndex,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  size 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ insert() [6/6]

virtual void geos::EpetraMatrix::insert ( globalIndex const  rowIndex,
globalIndex const  colIndex,
real64 const  value 
)
overridevirtual

Insert one element.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ insertable()

bool geos::MatrixBase< MATRIX, VECTOR >::insertable
inline

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

virtual globalIndex geos::EpetraMatrix::iupper ( ) const
overridevirtual

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

virtual globalIndex geos::EpetraMatrix::jlower ( ) const
overridevirtual

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

virtual globalIndex geos::EpetraMatrix::jupper ( ) const
overridevirtual

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

virtual void geos::EpetraMatrix::leftMultiplyTranspose ( EpetraMatrix const &  src,
EpetraMatrix dst 
) const
overridevirtual

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.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ leftRightScale()

virtual void geos::EpetraMatrix::leftRightScale ( EpetraVector const &  vecLeft,
EpetraVector const &  vecRight 
)
overridevirtual

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.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ leftScale()

virtual void geos::EpetraMatrix::leftScale ( EpetraVector const &  vec)
overridevirtual

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

Parameters
vecVector to pre-multiply with.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ maxRowLength()

virtual localIndex geos::EpetraMatrix::maxRowLength ( ) const
overridevirtual

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

Returns
the max length of a row

Collective.

◆ modifiable()

bool geos::MatrixBase< MATRIX, VECTOR >::modifiable
inline

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

virtual void geos::EpetraMatrix::multiply ( EpetraMatrix const &  src,
EpetraMatrix dst 
) const
overridevirtual

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.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ multiplyPtAP()

virtual void geos::EpetraMatrix::multiplyPtAP ( EpetraMatrix const &  P,
EpetraMatrix dst 
) const
overridevirtual

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

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

Reimplemented from geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ multiplyRAP()

virtual void geos::EpetraMatrix::multiplyRAP ( EpetraMatrix const &  R,
EpetraMatrix const &  P,
EpetraMatrix dst 
) const
overridevirtual

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 from geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ norm1()

virtual real64 geos::EpetraMatrix::norm1 ( ) const
overridevirtual

Returns the one norm of the matrix.

Returns
the value of 1-norm

◆ normFrobenius()

virtual real64 geos::EpetraMatrix::normFrobenius ( ) const
overridevirtual

Returns the Frobenius norm of the matrix.

Returns
the value of Frobenius norm

◆ normInf()

virtual real64 geos::EpetraMatrix::normInf ( ) const
overridevirtual

Returns the infinity norm of the matrix.

Returns
the value of infinity norm

◆ normMax() [1/2]

virtual real64 geos::EpetraMatrix::normMax ( ) const
overridevirtual

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

Returns
the value of max norm

◆ normMax() [2/2]

virtual real64 geos::EpetraMatrix::normMax ( arrayView1d< globalIndex const > const &  rowIndices) const
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ numGlobalCols()

virtual globalIndex geos::EpetraMatrix::numGlobalCols ( ) const
overridevirtual

Get the number of global columns.

Returns
Number of global columns in the operator.

◆ numGlobalNonzeros()

virtual globalIndex geos::EpetraMatrix::numGlobalNonzeros ( ) const
overridevirtual

Returns the total number of nonzeros in the matrix.

Returns
the total number of nonzeros in the matrix

◆ numGlobalRows()

virtual globalIndex geos::EpetraMatrix::numGlobalRows ( ) const
overridevirtual

Get the number of global rows.

Returns
Number of global rows in the operator.

◆ numLocalCols()

virtual localIndex geos::EpetraMatrix::numLocalCols ( ) const
overridevirtual

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

virtual localIndex geos::EpetraMatrix::numLocalNonzeros ( ) const
overridevirtual

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

virtual localIndex geos::EpetraMatrix::numLocalRows ( ) const
overridevirtual

Get the number of local rows.

Returns
Number of local rows in the operator.

◆ open()

virtual void geos::EpetraMatrix::open ( )
overridevirtual

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.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ operator=() [1/2]

EpetraMatrix& geos::EpetraMatrix::operator= ( EpetraMatrix &&  src)
noexcept

Move assignment.

Parameters
srcmatrix to be moved from.
Returns
the new matrix.

◆ operator=() [2/2]

EpetraMatrix& geos::EpetraMatrix::operator= ( EpetraMatrix const &  src)

Copy assignment.

Parameters
srcmatrix to be copied.
Returns
the new vector.

◆ print()

virtual void geos::EpetraMatrix::print ( std::ostream &  os = std::cout) const
overridevirtual

Print the matrix in Trilinos format to a stream.

Parameters
osthe output stream

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ ready()

bool geos::MatrixBase< MATRIX, VECTOR >::ready
inline

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

virtual void geos::EpetraMatrix::rescaleRows ( arrayView1d< globalIndex const > const &  rowIndices,
RowSumType const  rowSumType 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ residual()

virtual void geos::MatrixBase< MATRIX, VECTOR >::residual
inlineoverride

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.

Definition at line 550 of file MatrixBase.hpp.

◆ rightMultiplyTranspose()

virtual void geos::EpetraMatrix::rightMultiplyTranspose ( EpetraMatrix const &  src,
EpetraMatrix dst 
) const
overridevirtual

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.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ rightScale()

virtual void geos::EpetraMatrix::rightScale ( EpetraVector const &  vec)
overridevirtual

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

Parameters
vecVector to post-multiply with.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ rowLength()

virtual localIndex geos::EpetraMatrix::rowLength ( globalIndex const  globalRowIndex) const
overridevirtual

Get row length via global row index.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ scale()

virtual void geos::EpetraMatrix::scale ( real64 const  scalingFactor)
overridevirtual

Multiply all elements by scalingFactor.

Parameters
scalingFactorScaling factor.

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ separateComponentFilter()

virtual void geos::EpetraMatrix::separateComponentFilter ( EpetraMatrix dst,
integer const  dofPerPoint 
) const
overridevirtual

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

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ set() [1/6]

virtual void geos::EpetraMatrix::set ( arraySlice1d< globalIndex const > const &  rowIndices,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice2d< real64 const > const &  values 
)
overridevirtual

Set a dense block of values.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ set() [2/6]

virtual void geos::EpetraMatrix::set ( globalIndex const *  rowIndices,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  numRows,
localIndex const  numCols 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ set() [3/6]

virtual void geos::EpetraMatrix::set ( globalIndex const  rowIndex,
arraySlice1d< globalIndex const > const &  colIndices,
arraySlice1d< real64 const > const &  values 
)
overridevirtual

Set elements of one row using array1d.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ set() [4/6]

virtual void geos::EpetraMatrix::set ( globalIndex const  rowIndex,
globalIndex const *  colIndices,
real64 const *  values,
localIndex const  size 
)
overridevirtual

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ set() [5/6]

virtual void geos::EpetraMatrix::set ( globalIndex const  rowIndex,
globalIndex const  colIndex,
real64 const  value 
)
overridevirtual

Set one element.

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

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ set() [6/6]

virtual void geos::EpetraMatrix::set ( real64 const  value)
overridevirtual

Set all non-zero elements to a value.

Parameters
valuethe value to set all elements to

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ setDofManager()

void geos::MatrixBase< MATRIX, VECTOR >::setDofManager
inline

Associate a DofManager with this matrix.

Parameters
dofManagerthe DofManager containing the relevant degrees of freedom

Definition at line 156 of file MatrixBase.hpp.

◆ transpose()

virtual void geos::EpetraMatrix::transpose ( EpetraMatrix dst) const
overridevirtual

Matrix transposition.

Compute B = this^T.

Parameters
dstOutput matrix (B).

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.

◆ unwrapped() [1/2]

Epetra_FECrsMatrix& geos::EpetraMatrix::unwrapped ( )

Returns a non-const pointer to the underlying matrix.

Returns
non-const pointer to the underlying matrix

◆ unwrapped() [2/2]

const Epetra_FECrsMatrix& geos::EpetraMatrix::unwrapped ( ) const

Returns a const pointer to the underlying matrix.

Returns
const pointer to the underlying matrix

◆ write()

virtual void geos::EpetraMatrix::write ( string const &  filename,
LAIOutputFormat const  format = LAIOutputFormat::MATRIX_MARKET 
) const
overridevirtual

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)

Implements geos::MatrixBase< EpetraMatrix, EpetraVector >.


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