GEOSX
PetscMatrix.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_PETSCSPARSEMATRIX_HPP_
20 #define GEOSX_LINEARALGEBRA_INTERFACES_PETSCSPARSEMATRIX_HPP_
21 
22 #include "common/DataTypes.hpp"
23 #include "PetscVector.hpp"
26 
33 
36 extern "C" struct _p_Mat;
37 
39 
40 namespace geosx
41 {
42 
47 class PetscMatrix final : public virtual LinearOperator< PetscVector >,
48  private MatrixBase< PetscMatrix, PetscVector >
49 {
50 public:
51 
54 
56  using Mat = struct _p_Mat *;
57 
61 
66  PetscMatrix();
67 
72  PetscMatrix( PetscMatrix const & src );
73 
77  ~PetscMatrix() override;
78 
80 
84 
88  using MatrixBase::create;
89  using MatrixBase::closed;
93  using MatrixBase::ready;
95 
96  virtual void createWithLocalSize( localIndex const localRows,
97  localIndex const localCols,
98  localIndex const maxEntriesPerRow,
99  MPI_Comm const & comm ) override;
100 
101  virtual void createWithGlobalSize( globalIndex const globalRows,
102  globalIndex const globalCols,
103  localIndex const maxEntriesPerRow,
104  MPI_Comm const & comm ) override;
105 
109  virtual bool created() const override;
110 
111  virtual void reset() override;
112 
113  virtual void set( real64 const value ) override;
114 
115  virtual void zero() override;
116 
117  virtual void open() override;
118 
119  virtual void close() override;
120 
121  virtual void add( globalIndex const rowIndex,
122  globalIndex const colIndex,
123  real64 const value ) override;
124 
125  virtual void set( globalIndex const rowIndex,
126  globalIndex const colIndex,
127  real64 const value ) override;
128 
129  virtual void insert( globalIndex const rowIndex,
130  globalIndex const colIndex,
131  real64 const value ) override;
132 
133  virtual void add( globalIndex const rowIndex,
134  globalIndex const * colIndices,
135  real64 const * values,
136  localIndex const size ) override;
137 
138  virtual void set( globalIndex const rowIndex,
139  globalIndex const * colIndices,
140  real64 const * values,
141  localIndex const size ) override;
142 
143  virtual void insert( globalIndex const rowIndex,
144  globalIndex const * colIndices,
145  real64 const * values,
146  localIndex const size ) override;
147 
148  virtual void add( globalIndex const rowIndex,
149  arraySlice1d< globalIndex const > const & colIndices,
150  arraySlice1d< real64 const > const & values ) override;
151 
152  virtual void set( globalIndex const rowIndex,
153  arraySlice1d< globalIndex const > const & colIndices,
154  arraySlice1d< real64 const > const & values ) override;
155 
156  virtual void insert( globalIndex const rowIndex,
157  arraySlice1d< globalIndex const > const & colIndices,
158  arraySlice1d< real64 const > const & values ) override;
159 
160  virtual void add( arraySlice1d< globalIndex const > const & rowIndices,
161  arraySlice1d< globalIndex const > const & colIndices,
163 
164  virtual void set( arraySlice1d< globalIndex const > const & rowIndices,
165  arraySlice1d< globalIndex const > const & colIndices,
167 
168  virtual void insert( arraySlice1d< globalIndex const > const & rowIndices,
169  arraySlice1d< globalIndex const > const & colIndices,
171 
172  virtual void add( arraySlice1d< globalIndex const > const & rowIndices,
173  arraySlice1d< globalIndex const > const & colIndices,
175 
176  virtual void set( arraySlice1d< globalIndex const > const & rowIndices,
177  arraySlice1d< globalIndex const > const & colIndices,
179 
180  virtual void insert( arraySlice1d< globalIndex const > const & rowIndices,
181  arraySlice1d< globalIndex const > const & colIndices,
183 
184  virtual void add( globalIndex const * rowIndices,
185  globalIndex const * colIndices,
186  real64 const * values,
187  localIndex const numRows,
188  localIndex const numCols ) override;
189 
190  virtual void set( globalIndex const * rowIndices,
191  globalIndex const * colIndices,
192  real64 const * values,
193  localIndex const numRows,
194  localIndex const numCols ) override;
195 
196  virtual void insert( globalIndex const * rowIndices,
197  globalIndex const * colIndices,
198  real64 const * values,
199  localIndex const numRows,
200  localIndex const numCols ) override;
201 
202  virtual void apply( PetscVector const & src,
203  PetscVector & dst ) const override;
204 
205  virtual void multiply( PetscMatrix const & src,
206  PetscMatrix & dst ) const override;
207 
208  virtual void applyTranspose( Vector const & src,
209  Vector & dst ) const override;
210 
211  virtual void leftMultiplyTranspose( PetscMatrix const & src,
212  PetscMatrix & dst ) const override;
213 
214  virtual void rightMultiplyTranspose( PetscMatrix const & src,
215  PetscMatrix & dst ) const override;
216 
217  virtual void multiplyRAP( PetscMatrix const & R,
218  PetscMatrix const & P,
219  PetscMatrix & dst ) const override;
220 
221  virtual void multiplyPtAP( PetscMatrix const & P,
222  PetscMatrix & dst ) const override;
223 
224  virtual void gemv( real64 const alpha,
225  PetscVector const & x,
226  real64 const beta,
227  PetscVector & y,
228  bool useTranspose = false ) const override;
229 
230  virtual void scale( real64 const scalingFactor ) override;
231 
232  virtual void leftScale( PetscVector const & vec ) override;
233 
234  virtual void rightScale( PetscVector const & vec ) override;
235 
236  virtual void leftRightScale( PetscVector const & vecLeft,
237  PetscVector const & vecRight ) override;
238 
239  virtual void transpose( PetscMatrix & dst ) const override;
240 
241  virtual real64 clearRow( globalIndex const row,
242  bool const keepDiag = false,
243  real64 const diagValue = 0.0 ) override;
244 
245  virtual void addEntries( PetscMatrix const & src,
246  real64 const scale = 1.0,
247  bool const samePattern = true ) override;
248 
249  virtual void addDiagonal( PetscVector const & src ) override;
250 
254  virtual localIndex maxRowLength() const override;
255 
256  virtual localIndex localRowLength( localIndex localRowIndex ) const override;
257 
258  virtual localIndex globalRowLength( globalIndex globalRowIndex ) const override;
259 
260  virtual real64 getDiagValue( globalIndex globalRow ) const override;
261 
262  virtual void extractDiagonal( PetscVector & dst ) const override;
263 
264  virtual void getRowCopy( globalIndex globalRow,
265  arraySlice1d< globalIndex > const & colIndices,
266  arraySlice1d< real64 > const & values ) const override;
267 
271  virtual globalIndex numGlobalRows() const override;
272 
276  virtual globalIndex numGlobalCols() const override;
277 
281  virtual localIndex numLocalRows() const override;
282 
286  virtual localIndex numLocalCols() const override;
287 
291  virtual globalIndex ilower() const override;
292 
296  virtual globalIndex iupper() const override;
297 
301  virtual globalIndex jlower() const override;
302 
306  virtual globalIndex jupper() const override;
307 
311  virtual localIndex numLocalNonzeros() const override;
312 
316  virtual globalIndex numGlobalNonzeros() const override;
317 
321  virtual real64 normInf() const override;
322 
326  virtual real64 norm1() const override;
327 
331  virtual real64 normFrobenius() const override;
332 
333  virtual localIndex getLocalRowID( globalIndex const index ) const override;
334 
335  virtual globalIndex getGlobalRowID( localIndex const index ) const override;
336 
340  virtual MPI_Comm getComm() const override;
341 
342  virtual void print( std::ostream & os = std::cout ) const override;
343 
344  virtual void write( string const & filename,
345  LAIOutputFormat const format = LAIOutputFormat::MATRIX_MARKET ) const override;
346 
348 
353  const Mat & unwrapped() const;
354 
359  Mat & unwrapped();
360 
361 private:
362 
364  Mat m_mat;
365 
367  array1d< globalIndex > m_rowsToClear;
368 
370  array1d< real64 > m_diagValues;
371 
372 };
373 
374 } // namespace geosx
375 
376 #endif /*GEOSX_LINEARALGEBRA_INTERFACES_PETSCSPARSEMATRIX_HPP_*/
virtual void add(globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
Add to one element.
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 multiplyRAP(PetscMatrix const &R, PetscMatrix const &P, PetscMatrix &dst) const override
Compute the triple product dst = R * this * P
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 globalIndex numGlobalNonzeros() const override
Returns the total number of nonzeros in the matrix.
virtual localIndex globalRowLength(globalIndex globalRowIndex) const override
Get row length via global row index.
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 localIndex numLocalRows() const override
Return the local number of columns on each processor.
virtual globalIndex iupper() const override
Returns index one past the last global row owned by that processor.
long long int globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:128
virtual void leftMultiplyTranspose(PetscMatrix const &src, PetscMatrix &dst) const override
Matrix/Matrix transpose multiplication.
virtual void addDiagonal(PetscVector const &src) override
Add entries of a vector to the diagonal of this matrix.
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.
This class serves to provide a sliced multidimensional interface to the family of LvArray classes...
Definition: ArraySlice.hpp:89
virtual localIndex numLocalCols() const override
Return the local number of columns on each processor.
virtual void extractDiagonal(PetscVector &dst) const override
Extract diagonal values into a vector.
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 print(std::ostream &os=std::cout) const override
Print the matrix in Trilinos format to a stream.
virtual void scale(real64 const scalingFactor) override
Multiply all elements by scalingFactor.
virtual localIndex maxRowLength() const override
Returns the number of nonzero entries in the longest row of the matrix.
bool ready() const
Query matrix ready status.
Definition: MatrixBase.hpp:125
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 real64 getDiagValue(globalIndex globalRow) const override
get diagonal element value on a given row
virtual void zero() override
Set all elements to zero.
virtual globalIndex jlower() const override
Returns the index of the first global col owned by that processor.
~PetscMatrix() override
Destructor.
struct _p_Mat * Mat
Alias for PETSc matrix struct pointer.
Definition: PetscMatrix.hpp:56
virtual void addEntries(PetscMatrix const &src, real64 const scale=1.0, bool const samePattern=true) override
Add entries of another matrix to this.
virtual void leftRightScale(PetscVector const &vecLeft, PetscVector const &vecRight) override
Post-multiplies (right) with diagonal matrix consisting of the values in vecRight and pre-multiplies ...
const Mat & unwrapped() const
Returns a const pointer to the underlying matrix.
virtual void applyTranspose(Vector const &src, Vector &dst) const override
Apply transpose of the matrix to a vector.
double real64
64-bit floating point type.
Definition: DataTypes.hpp:136
This class creates and provides basic support for Vec vector object type used in PETSc.
Definition: PetscVector.hpp:44
virtual localIndex localRowLength(localIndex localRowIndex) const override
Get row length via local row index.
virtual void leftScale(PetscVector const &vec) override
Pre-multiplies (left) with diagonal matrix consisting of the values in vec.
virtual void multiplyPtAP(PetscMatrix const &P, PetscMatrix &dst) const override
Compute the triple product dst = P^T * this * P
virtual globalIndex numGlobalCols() const override
Returns the number of global columns.
virtual bool created() const override
Returns the number of global rows.
virtual real64 normFrobenius() const override
Returns the Frobenius norm of the matrix.
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 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 globalIndex ilower() const override
Returns the index of the first global row owned by that processor.
PetscMatrix()
Empty matrix constructor.
virtual localIndex getLocalRowID(globalIndex const index) const override
Map a global row index to local row index.
virtual localIndex numLocalNonzeros() const override
Returns the number of nonzeros in the local portion of the matrix.
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 rightMultiplyTranspose(PetscMatrix const &src, PetscMatrix &dst) const override
Matrix/Matrix transpose multiplication.
LAIOutputFormat
Definition: common.hpp:139
virtual void apply(PetscVector const &src, PetscVector &dst) const override
Apply operator to a vector.
This class creates and provides basic support for the Mat matrix object type used in PETSc...
Definition: PetscMatrix.hpp:47
bool closed() const
Query matrix closed status.
Definition: MatrixBase.hpp:112
virtual void close() override
Assemble and compress the matrix.
bool modifiable() const
Query matrix status.
Definition: MatrixBase.hpp:133
virtual void rightScale(PetscVector const &vec) override
Post-multiplies (right) with diagonal matrix consisting of the values in vec.
virtual void transpose(PetscMatrix &dst) const override
Matrix transposition.
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 multiply(PetscMatrix const &src, PetscMatrix &dst) const override
Matrix/Matrix multiplication.
virtual void gemv(real64 const alpha, PetscVector const &x, real64 const beta, PetscVector &y, bool useTranspose=false) const override
Compute gemv y = alpha*A*x + beta*y.
virtual void residual(Vector const &x, Vector const &b, Vector &r) const override
Compute residual r = Ax - b.
Definition: MatrixBase.hpp:546
virtual real64 norm1() const override
Returns the one norm of the matrix.
Abstract base class for linear operators.
virtual real64 normInf() const override
Returns the infinity norm of the matrix.
virtual void insert(globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
Insert one element.
virtual globalIndex jupper() const override
Returns index one past the last global col owned by that processor.
This class provides a fixed dimensional resizeable array interface in addition to an interface simila...
Definition: Array.hpp:55
virtual void open() override
Open matrix for adding new entries.
virtual MPI_Comm getComm() const override
Get the MPI communicator the matrix was created with.
Common base template for all matrix wrapper types.
Definition: MatrixBase.hpp:49
virtual void reset() override
Reset the matrix to default state.
virtual globalIndex getGlobalRowID(localIndex const index) const override
Map a local row index to global row index.
virtual globalIndex numGlobalRows() const override
Returns the number of global rows.