20 #ifndef GEOS_LINEARALGEBRA_INTERFACES_EPETRAMATRIX_HPP_
21 #define GEOS_LINEARALGEBRA_INTERFACES_EPETRAMATRIX_HPP_
30 class Epetra_CrsMatrix;
31 class Epetra_FECrsMatrix;
40 private MatrixBase< EpetraMatrix, EpetraVector >
119 MPI_Comm
const &
comm )
override;
124 MPI_Comm
const &
comm )
override;
143 real64 const value )
override;
147 real64 const value )
override;
151 real64 const value )
override;
240 bool useTranspose =
false )
const override;
259 integer const dofsPerNode )
const override;
262 bool const keepDiag =
false,
263 real64 const diagValue = 0.0 )
override;
274 bool const excludeDiag )
override;
386 virtual MPI_Comm
comm()
const override;
388 virtual void print( std::ostream & os = std::cout )
const override;
390 virtual void write(
string const & filename,
391 LAIOutputFormat const format = LAIOutputFormat::MATRIX_MARKET )
const override;
421 void create( Epetra_CrsMatrix
const & src );
424 std::unique_ptr< Epetra_FECrsMatrix > m_matrix;
427 std::unique_ptr< Epetra_Map > m_src_map;
430 std::unique_ptr< Epetra_Map > m_dst_map;
Facilitates exporting Epetra matrix and associated vector objects (either in parallel or serial).
Wrapper class for Epetra's CrsMatrix.
virtual void rightMultiplyTranspose(EpetraMatrix const &src, EpetraMatrix &dst) const override
Matrix/Matrix transpose multiplication.
virtual void write(string const &filename, LAIOutputFormat const format=LAIOutputFormat::MATRIX_MARKET) const override
Write the matrix to filename in a matlab-compatible format.
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.
virtual globalIndex numGlobalCols() const override
virtual void open() override
Open matrix for adding new entries.
virtual void rescaleRows(arrayView1d< globalIndex const > const &rowIndices, RowSumType const rowSumType) override
Rescales selected rows of matrix using row sum reciprocal as a factor.
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.
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.
virtual void extractLocal(CRSMatrixView< real64, localIndex > const &localMat) const override
Extract block of the matrix corresponding to locally owned columns using previously allocated storage...
virtual globalIndex iupper() const override
Returns index one past the last global row owned by that processor.
virtual void rightScale(EpetraVector const &vec) override
Post-multiplies (right) with diagonal matrix consisting of the values in vec.
virtual localIndex rowLength(globalIndex const globalRowIndex) const override
Get row length via global row index.
virtual localIndex maxRowLength() const override
Returns the number of nonzero entries in the longest row of the matrix.
virtual void apply(EpetraVector const &src, EpetraVector &dst) const override
Apply operator to a vector, dst = this(src).
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.
virtual void insert(globalIndex const *rowIndices, globalIndex const *colIndices, real64 const *values, localIndex const numRows, localIndex const numCols) override
Insert dense matrix.
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.
virtual globalIndex numGlobalNonzeros() const override
Epetra_FECrsMatrix & unwrapped()
Returns a non-const pointer to the underlying matrix.
EpetraMatrix & operator=(EpetraMatrix &&src) noexcept
Move assignment.
virtual ~EpetraMatrix() override
Virtual destructor.
virtual void multiply(EpetraMatrix const &src, EpetraMatrix &dst) const override
Matrix/Matrix multiplication.
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.
virtual void extract(CRSMatrixView< real64, globalIndex const > const &localMat) const override
Extract local part of the matrix using previously allocated storage and structure.
virtual void clampEntries(real64 const lo, real64 const hi, bool const excludeDiag) override
Clamp each matrix value between values of lo and hi.
virtual void getRowSums(EpetraVector &dst, RowSumType const rowSumType) const override
Populate a vector with row sums of this.
virtual void add(globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
Add to one element.
EpetraMatrix & operator=(EpetraMatrix const &src)
Copy assignment.
virtual void getRowCopy(globalIndex globalRow, arraySlice1d< globalIndex > const &colIndices, arraySlice1d< real64 > const &values) const override
Returns a copy of the data in row globalRow.
virtual void transpose(EpetraMatrix &dst) const override
Matrix transposition.
virtual void computeScalingVector(EpetraVector &scaling) const override
Compute left and right scaling vectors for diagonal scaling.
virtual void insert(globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values) override
Insert elements of one row using array1d.
virtual void insert(globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
Insert one element.
virtual void separateComponentFilter(EpetraMatrix &dst, integer const dofsPerNode) const override
Apply a separate component approximation (filter) to this matrix.
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.
virtual void extract(CRSMatrixView< real64, globalIndex > const &localMat) const override
Extract local part of the matrix using previously allocated storage.
virtual void print(std::ostream &os=std::cout) const override
Print the matrix in Trilinos format to a stream.
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.
virtual void leftMultiplyTranspose(EpetraMatrix const &src, EpetraMatrix &dst) const override
Matrix/Matrix transpose multiplication.
virtual globalIndex jlower() const override
Returns the index of the first global col owned by that processor.
virtual globalIndex ilower() const override
Returns the index of the first global row owned by that processor.
virtual globalIndex numGlobalRows() const override
virtual real64 normFrobenius() const override
Returns the Frobenius norm of the matrix.
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 ...
virtual localIndex getLocalRowID(globalIndex const index) const override
Map a global row index to local row index.
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.
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.
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.
virtual void zero() override
Set all elements to zero.
virtual globalIndex getGlobalRowID(localIndex const index) const override
Map a local row index to global row index.
EpetraMatrix(EpetraMatrix &&src) noexcept
Move constructor.
virtual real64 normMax() const override
Returns the max norm of the matrix (the largest absolute element value).
virtual void reset() override
Reset the matrix to default state.
virtual real64 normInf() const override
Returns the infinity norm of the matrix.
virtual void set(real64 const value) override
Set all non-zero elements to a value.
virtual MPI_Comm comm() const override
Get the MPI communicator the matrix was created with.
virtual localIndex maxRowLengthLocal() const override
Returns the number of nonzero entries in the longest local row of the matrix.
virtual void add(globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values) override
Add elements to one row using array1d.
virtual real64 normMax(arrayView1d< globalIndex const > const &rowIndices) const override
Returns the max norm of the matrix on a subset of rows.
virtual void extractDiagonal(EpetraVector &dst) const override
Extract diagonal values into a vector.
virtual void getRowLocalLengths(arrayView1d< localIndex > const &lengths) const override
Get the local row lengths of every local row (number of nonzeros in local columns)
virtual localIndex numLocalRows() const override
virtual void addEntries(EpetraMatrix const &src, MatrixPatternOp const op, real64 const scale) override
Add entries of another matrix to this.
virtual localIndex numLocalNonzeros() const override
EpetraMatrix(EpetraMatrix const &src)
Copy constructor.
virtual bool created() const override
Query matrix creation status.
Epetra_FECrsMatrix const & unwrapped() const
Returns a const pointer to the underlying matrix.
virtual void close() override
Assemble and compress the matrix.
virtual void addDiagonal(EpetraVector const &src, real64 const scale) override
Add (scaled) entries of a vector to the diagonal of this matrix.
virtual void leftScale(EpetraVector const &vec) override
Pre-multiplies (left) with diagonal matrix consisting of the values in vec.
virtual void getRowLengths(arrayView1d< localIndex > const &lengths) const override
Get the row lengths of every local row.
virtual void scale(real64 const scalingFactor) override
Multiply all elements by scalingFactor.
virtual void multiplyPtAP(EpetraMatrix const &P, EpetraMatrix &dst) const override
Compute the triple product dst = P^T * this * P
virtual globalIndex jupper() const override
Returns index one past the last global col owned by that processor.
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.
virtual void multiplyRAP(EpetraMatrix const &R, EpetraMatrix const &P, EpetraMatrix &dst) const override
Compute the triple product dst = R * this * P
virtual real64 norm1() const override
Returns the one norm of the matrix.
virtual void applyTranspose(EpetraVector const &src, EpetraVector &dst) const override
Apply transpose of the matrix to a vector.
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.
virtual localIndex numLocalCols() const override
EpetraMatrix()
Empty matrix constructor.
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.
virtual void set(globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
Set one element.
virtual void set(globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values) override
Set elements of one row using array1d.
Wrapper around Trilinos' Epetra_Vector object.
Abstract base class for linear operators.
Common base template for all matrix wrapper types.
virtual void multiplyPtAP(Matrix const &P, Matrix &dst) const
Compute the triple product dst = P^T * this * P
virtual void residual(Vector const &x, Vector const &b, Vector &r) const override
Compute residual r = b - A * x.
void setDofManager(DofManager const *const dofManager)
Associate a DofManager with this matrix.
bool assembled() const
Query matrix assembled status.
bool closed() const
Query matrix closed status.
virtual void createWithGlobalSize(globalIndex const globalSize, localIndex const maxEntriesPerRow, MPI_Comm const &comm)
Create a square matrix from global number of rows.
bool insertable() const
Query matrix status.
virtual void createWithLocalSize(localIndex const localSize, localIndex const maxEntriesPerRow, MPI_Comm const &comm)
Create a square matrix from local number of rows.
CRSMatrix< real64, globalIndex > extract() const
Extract local part of the matrix.
DofManager const * dofManager() const
bool ready() const
Query matrix ready status.
bool modifiable() const
Query matrix status.
virtual void multiplyP1tAP2(Matrix const &P1, Matrix const &P2, Matrix &dst) const
Compute the triple product dst = P1^T * this * P2
CRSMatrix< real64, localIndex > extractLocal() const
Extract block of the matrix corresponding to locally owned columns.
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.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
RowSumType
Type of row sum to compute.
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
ArraySlice< T, 2, USD > arraySlice2d
Alias for 2D array slice.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
ArraySlice< T, 1, USD > arraySlice1d
Alias for 1D array slice.
LvArray::CRSMatrixView< T, COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
int integer
Signed integer type.
MatrixPatternOp
Describes relationship between and treatment of nonzero patterns of arguments in matrix functions lik...