19 #ifndef GEOSX_LINEARALGEBRA_INTERFACES_HYPREMATRIX_HPP_ 20 #define GEOSX_LINEARALGEBRA_INTERFACES_HYPREMATRIX_HPP_ 36 extern "C" struct hypre_IJMatrix_struct;
42 extern "C" struct hypre_ParCSRMatrix_struct;
110 MPI_Comm
const & comm )
override;
115 MPI_Comm
const & comm )
override;
117 virtual void open()
override;
119 virtual void close()
override;
124 virtual bool created()
const override;
126 virtual void reset()
override;
128 virtual void set(
real64 const value )
override;
130 virtual void zero()
override;
134 real64 const value )
override;
138 real64 const value )
override;
142 real64 const value )
override;
217 Vector & dst )
const override;
239 bool useTranspose =
false )
const override;
241 virtual void scale(
real64 const scalingFactor )
override;
253 bool const keepDiag =
false,
254 real64 const diagValue = 0.0 )
override;
258 bool const samePattern =
true )
override;
351 virtual MPI_Comm
getComm()
const override;
353 virtual void print( std::ostream & os = std::cout )
const override;
355 virtual void write(
string const & filename,
HypreMatrix()
Empty matrix constructor.
Wrapper class for hypre's ParVector.
virtual real64 normFrobenius() const override
Returns the Frobenius norm of the matrix.
virtual void insert(globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
Insert one element.
virtual void createWithLocalSize(localIndex const localSize, localIndex const maxEntriesPerRow, MPI_Comm const &comm)
Create a square matrix from local number of rows.
bool assembled() const
Query matrix assembled status.
virtual void addDiagonal(HypreVector const &src) override
Add entries of a vector to the diagonal of this matrix.
virtual void open() override
Open matrix for adding new entries.
virtual globalIndex jlower() const override
Returns the index of the first global col owned by that processor.
Wrapper class for hypre's ParCSRMatrix.
virtual localIndex numLocalNonzeros() const override
Returns the number of nonzeros in the local portion of the matrix.
virtual localIndex numLocalRows() const override
Return the local number of columns on each processor.
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 apply(HypreVector const &src, HypreVector &dst) const override
Apply operator to a vector.
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.
long long int globalIndex
Global index type (for indexing objects across MPI partitions).
This class serves to provide a sliced multidimensional interface to the family of LvArray classes...
std::string format(int NDIM, INDEX_TYPE const *const dims)
This function returns a string that may be used as the "type" in a call to TV_ttf_add_row(). This will either be a single value or an array.
virtual void leftScale(HypreVector const &vec) override
Pre-multiplies (left) with diagonal matrix consisting of the values in vec.
virtual void add(globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
Add to one element.
virtual void addEntries(HypreMatrix const &src, real64 const scale=1.0, bool const samePattern=true) override
Add entries of another matrix to this.
virtual real64 normInf() const override
Returns the infinity norm of the matrix.
virtual void leftMultiplyTranspose(HypreMatrix const &src, HypreMatrix &dst) const override
Matrix/Matrix transpose multiplication.
bool ready() const
Query matrix ready status.
virtual localIndex globalRowLength(globalIndex globalRowIndex) const override
Get row length via global row index.
~HypreMatrix() override
Virtual destructor.
virtual void extractDiagonal(HypreVector &dst) const override
Extract diagonal values into a vector.
virtual MPI_Comm getComm() const override
Get the MPI communicator the matrix was created with.
double real64
64-bit floating point type.
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 globalIndex getGlobalRowID(localIndex const index) const override
Map a local row index to global row index.
virtual localIndex localRowLength(localIndex localRowIndex) const override
Get row length via local row index.
hypre_IJMatrix_struct * HYPRE_IJMatrix
IJMatrix pointer alias.
virtual void close() override
Assemble and compress the matrix.
virtual void transpose(HypreMatrix &dst) const override
Matrix transposition.
virtual globalIndex iupper() const override
Returns index one past the last global row owned by that processor.
virtual void print(std::ostream &os=std::cout) const override
Print the matrix in Trilinos format to a stream.
virtual globalIndex jupper() const override
Returns index one past the last global col owned by that processor.
virtual void create(CRSMatrixView< real64 const, globalIndex const > const &localMatrix, MPI_Comm const &comm)
Create parallel matrix from a local CRS 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 globalIndex numGlobalRows() const override
Returns the number of global rows.
virtual void createWithGlobalSize(globalIndex const globalSize, localIndex const maxEntriesPerRow, MPI_Comm const &comm)
Create a square matrix from global number of rows.
virtual void multiplyPtAP(HypreMatrix const &P, HypreMatrix &dst) const override
Compute the triple product dst = P^T * this * P
virtual void rightMultiplyTranspose(HypreMatrix const &src, HypreMatrix &dst) const override
Matrix/Matrix transpose multiplication.
virtual localIndex getLocalRowID(globalIndex const index) const override
Map a global row index to local row index.
virtual void zero() override
Set all elements to zero.
virtual void leftRightScale(HypreVector const &vecLeft, HypreVector const &vecRight) override
Post-multiplies (right) with diagonal matrix consisting of the values in vecRight and pre-multiplies ...
bool closed() const
Query matrix closed status.
virtual void reset() override
Reset the matrix to default state.
hypre_ParCSRMatrix_struct * HYPRE_ParCSRMatrix
ParCSRMatrix pointer alias.
virtual real64 norm1() const override
Returns the one norm of the matrix.
virtual void write(string const &filename, LAIOutputFormat const format=LAIOutputFormat::MATRIX_MARKET) const override
Write the matrix to filename in a matlab-compatible format.
bool modifiable() const
Query matrix status.
virtual real64 getDiagValue(globalIndex globalRow) const override
get diagonal element value on a given row
virtual localIndex numLocalCols() const override
Return the local number of columns on each processor.
bool insertable() const
Query matrix status.
std::ptrdiff_t localIndex
Local index type (for indexing objects within an MPI partition).
virtual void multiplyRAP(HypreMatrix const &R, HypreMatrix const &P, HypreMatrix &dst) const override
Compute the triple product dst = R * this * P
virtual bool created() const override
Query matrix creation status.
virtual globalIndex numGlobalCols() const override
Returns the number of global columns.
virtual void residual(Vector const &x, Vector const &b, Vector &r) const override
Compute residual r = Ax - b.
virtual void multiply(HypreMatrix const &src, HypreMatrix &dst) const override
Matrix/Matrix multiplication.
virtual localIndex maxRowLength() const override
Returns the number of nonzero entries in the longest row of the matrix.
virtual void gemv(real64 const alpha, HypreVector const &x, real64 const beta, HypreVector &y, bool useTranspose=false) const override
Compute gemv y = alpha*A*x + beta*y.
HYPRE_ParCSRMatrix const & unwrapped() const
Returns a pointer to implementation.
virtual void rightScale(HypreVector const &vec) override
Post-multiplies (right) with diagonal matrix consisting of the values in vec.
Abstract base class for linear operators.
virtual void scale(real64 const scalingFactor) override
Multiply all elements by scalingFactor.
HYPRE_IJMatrix const & unwrappedIJ() const
Returns a pointer to implementation.
virtual void applyTranspose(Vector const &src, Vector &dst) const override
Apply transpose of the matrix to a vector.
virtual globalIndex numGlobalNonzeros() const override
Returns the total number of nonzeros in the matrix.
Common base template for all matrix wrapper types.
virtual globalIndex ilower() const override
Returns the index of the first global row owned by that processor.