GEOSX
HypreMatrix.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2018-2020 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2020 The Board of Trustees of the Leland Stanford Junior University
7  * Copyright (c) 2018-2020 Total, S.A
8  * Copyright (c) 2019- GEOSX Contributors
9  * All rights reserved
10  *
11  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
12  * ------------------------------------------------------------------------------------------------------------
13  */
14 
19 #ifndef GEOSX_LINEARALGEBRA_INTERFACES_HYPREMATRIX_HPP_
20 #define GEOSX_LINEARALGEBRA_INTERFACES_HYPREMATRIX_HPP_
21 
22 #include "common/DataTypes.hpp"
23 #include "HypreVector.hpp"
26 
33 
36 extern "C" struct hypre_IJMatrix_struct;
37 
39 using HYPRE_IJMatrix = hypre_IJMatrix_struct *;
40 
42 extern "C" struct hypre_ParCSRMatrix_struct;
43 
45 using HYPRE_ParCSRMatrix = hypre_ParCSRMatrix_struct *;
46 
48 
49 namespace geosx
50 {
51 
58 class HypreMatrix final : public virtual LinearOperator< HypreVector >,
59  private MatrixBase< HypreMatrix, HypreVector >
60 {
61 public:
62 
65 
69 
77  HypreMatrix();
78 
83  HypreMatrix( HypreMatrix const & src );
84 
88  ~HypreMatrix() override;
89 
91 
95 
99  using MatrixBase::create;
100  using MatrixBase::closed;
101  using MatrixBase::assembled;
104  using MatrixBase::ready;
105  using MatrixBase::residual;
106 
107  virtual void createWithLocalSize( localIndex const localRows,
108  localIndex const localCols,
109  localIndex const maxEntriesPerRow,
110  MPI_Comm const & comm ) override;
111 
112  virtual void createWithGlobalSize( globalIndex const globalRows,
113  globalIndex const globalCols,
114  localIndex const maxEntriesPerRow,
115  MPI_Comm const & comm ) override;
116 
117  virtual void open() override;
118 
119  virtual void close() override;
120 
124  virtual bool created() const override;
125 
126  virtual void reset() override;
127 
128  virtual void set( real64 const value ) override;
129 
130  virtual void zero() override;
131 
132  virtual void add( globalIndex const rowIndex,
133  globalIndex const colIndex,
134  real64 const value ) override;
135 
136  virtual void set( globalIndex const rowIndex,
137  globalIndex const colIndex,
138  real64 const value ) override;
139 
140  virtual void insert( globalIndex const rowIndex,
141  globalIndex const colIndex,
142  real64 const value ) override;
143 
144  virtual void add( globalIndex const rowIndex,
145  globalIndex const * colIndices,
146  real64 const * values,
147  localIndex const size ) override;
148 
149  virtual void set( globalIndex const rowIndex,
150  globalIndex const * colIndices,
151  real64 const * values,
152  localIndex const size ) override;
153 
154  virtual void insert( globalIndex const rowIndex,
155  globalIndex const * colIndices,
156  real64 const * values,
157  localIndex const size ) override;
158 
159  virtual void add( globalIndex const rowIndex,
160  arraySlice1d< globalIndex const > const & colIndices,
161  arraySlice1d< real64 const > const & values ) override;
162 
163  virtual void set( globalIndex const rowIndex,
164  arraySlice1d< globalIndex const > const & colIndices,
165  arraySlice1d< real64 const > const & values ) override;
166 
167  virtual void insert( globalIndex const rowIndex,
168  arraySlice1d< globalIndex const > const & colIndices,
169  arraySlice1d< real64 const > const & values ) override;
170 
171  virtual void add( arraySlice1d< globalIndex const > const & rowIndices,
172  arraySlice1d< globalIndex const > const & colIndices,
174 
175  virtual void set( arraySlice1d< globalIndex const > const & rowIndices,
176  arraySlice1d< globalIndex const > const & colIndices,
178 
179  virtual void insert( arraySlice1d< globalIndex const > const & rowIndices,
180  arraySlice1d< globalIndex const > const & colIndices,
182 
183  virtual void add( arraySlice1d< globalIndex const > const & rowIndices,
184  arraySlice1d< globalIndex const > const & colIndices,
186 
187  virtual void set( arraySlice1d< globalIndex const > const & rowIndices,
188  arraySlice1d< globalIndex const > const & colIndices,
190 
191  virtual void insert( arraySlice1d< globalIndex const > const & rowIndices,
192  arraySlice1d< globalIndex const > const & colIndices,
194 
195  virtual void add( globalIndex const * rowIndices,
196  globalIndex const * colIndices,
197  real64 const * values,
198  localIndex const numRows,
199  localIndex const numCols ) override;
200 
201  virtual void set( globalIndex const * rowIndices,
202  globalIndex const * colIndices,
203  real64 const * values,
204  localIndex const numRows,
205  localIndex const numCols ) override;
206 
207  virtual void insert( globalIndex const * rowIndices,
208  globalIndex const * colIndices,
209  real64 const * values,
210  localIndex const numRows,
211  localIndex const numCols ) override;
212 
213  virtual void apply( HypreVector const & src,
214  HypreVector & dst ) const override;
215 
216  virtual void applyTranspose( Vector const & src,
217  Vector & dst ) const override;
218 
219  virtual void multiply( HypreMatrix const & src,
220  HypreMatrix & dst ) const override;
221 
222  virtual void leftMultiplyTranspose( HypreMatrix const & src,
223  HypreMatrix & dst ) const override;
224 
225  virtual void rightMultiplyTranspose( HypreMatrix const & src,
226  HypreMatrix & dst ) const override;
227 
228  virtual void multiplyRAP( HypreMatrix const & R,
229  HypreMatrix const & P,
230  HypreMatrix & dst ) const override;
231 
232  virtual void multiplyPtAP( HypreMatrix const & P,
233  HypreMatrix & dst ) const override;
234 
235  virtual void gemv( real64 const alpha,
236  HypreVector const & x,
237  real64 const beta,
238  HypreVector & y,
239  bool useTranspose = false ) const override;
240 
241  virtual void scale( real64 const scalingFactor ) override;
242 
243  virtual void leftScale( HypreVector const & vec ) override;
244 
245  virtual void rightScale( HypreVector const & vec ) override;
246 
247  virtual void leftRightScale( HypreVector const & vecLeft,
248  HypreVector const & vecRight ) override;
249 
250  virtual void transpose( HypreMatrix & dst ) const override;
251 
252  virtual real64 clearRow( globalIndex const row,
253  bool const keepDiag = false,
254  real64 const diagValue = 0.0 ) override;
255 
256  virtual void addEntries( HypreMatrix const & src,
257  real64 const scale = 1.0,
258  bool const samePattern = true ) override;
259 
260  virtual void addDiagonal( HypreVector const & src ) override;
261 
265  virtual localIndex maxRowLength() const override;
266 
267  virtual localIndex localRowLength( localIndex localRowIndex ) const override;
268 
269  virtual localIndex globalRowLength( globalIndex globalRowIndex ) const override;
270 
271  virtual void getRowCopy( globalIndex globalRow,
272  arraySlice1d< globalIndex > const & colIndices,
273  arraySlice1d< real64 > const & values ) const override;
274 
275  virtual real64 getDiagValue( globalIndex globalRow ) const override;
276 
277  virtual void extractDiagonal( HypreVector & dst ) const override;
278 
282  virtual globalIndex numGlobalRows() const override;
283 
287  virtual globalIndex numGlobalCols() const override;
288 
292  virtual localIndex numLocalRows() const override;
293 
297  virtual localIndex numLocalCols() const override;
298 
302  virtual globalIndex ilower() const override;
303 
307  virtual globalIndex iupper() const override;
308 
312  virtual globalIndex jlower() const override;
313 
317  virtual globalIndex jupper() const override;
318 
322  virtual localIndex numLocalNonzeros() const override;
323 
327  virtual globalIndex numGlobalNonzeros() const override;
328 
332  virtual real64 normInf() const override;
333 
337  virtual real64 norm1() const override;
338 
342  virtual real64 normFrobenius() const override;
343 
344  virtual localIndex getLocalRowID( globalIndex const index ) const override;
345 
346  virtual globalIndex getGlobalRowID( localIndex const index ) const override;
347 
351  virtual MPI_Comm getComm() const override;
352 
353  virtual void print( std::ostream & os = std::cout ) const override;
354 
355  virtual void write( string const & filename,
356  LAIOutputFormat const format = LAIOutputFormat::MATRIX_MARKET ) const override;
357 
359 
364  HYPRE_ParCSRMatrix const & unwrapped() const;
365 
370  HYPRE_IJMatrix const & unwrappedIJ() const;
371 
372 private:
373 
377  void parCSRtoIJ( HYPRE_ParCSRMatrix const & parCSRMatrix );
378 
382  HYPRE_IJMatrix m_ij_mat = nullptr;
383 
387  HYPRE_ParCSRMatrix m_parcsr_mat = nullptr;
388 
389 };
390 
391 } // namespace geosx
392 
393 #endif /*GEOSX_LINEARALGEBRA_INTERFACES_HYPREMATRIX_HPP_*/
HypreMatrix()
Empty matrix constructor.
Wrapper class for hypre&#39;s ParVector.
Definition: HypreVector.hpp:55
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.
Definition: MatrixBase.hpp:164
bool assembled() const
Query matrix assembled status.
Definition: MatrixBase.hpp:118
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&#39;s ParCSRMatrix.
Definition: HypreMatrix.hpp:58
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).
Definition: DataTypes.hpp:128
This class serves to provide a sliced multidimensional interface to the family of LvArray classes...
Definition: ArraySlice.hpp:89
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.
Definition: tv_helpers.hpp:37
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.
Definition: MatrixBase.hpp:125
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.
Definition: DataTypes.hpp:136
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.
Definition: HypreMatrix.hpp:39
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.
Definition: MatrixBase.hpp:225
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.
Definition: MatrixBase.hpp:182
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.
LAIOutputFormat
Definition: common.hpp:139
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.
Definition: MatrixBase.hpp:112
virtual void reset() override
Reset the matrix to default state.
hypre_ParCSRMatrix_struct * HYPRE_ParCSRMatrix
ParCSRMatrix pointer alias.
Definition: HypreMatrix.hpp:45
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.
Definition: MatrixBase.hpp:133
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.
Definition: MatrixBase.hpp:140
std::ptrdiff_t localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:125
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.
Definition: MatrixBase.hpp:546
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.
Definition: MatrixBase.hpp:49
virtual globalIndex ilower() const override
Returns the index of the first global row owned by that processor.