GEOS
EpetraMatrix.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 TotalEnergies
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
20 #ifndef GEOS_LINEARALGEBRA_INTERFACES_EPETRAMATRIX_HPP_
21 #define GEOS_LINEARALGEBRA_INTERFACES_EPETRAMATRIX_HPP_
22 
23 #include "common/DataTypes.hpp"
28 
29 class Epetra_Map;
30 class Epetra_CrsMatrix;
31 class Epetra_FECrsMatrix;
32 
33 namespace geos
34 {
35 
39 class EpetraMatrix final : public virtual LinearOperator< EpetraVector >,
40  private MatrixBase< EpetraMatrix, EpetraVector >
41 {
42 public:
43 
46 
49 
54 
61 
66  EpetraMatrix( EpetraMatrix const & src );
67 
72  EpetraMatrix( EpetraMatrix && src ) noexcept;
73 
80 
86  EpetraMatrix & operator=( EpetraMatrix && src ) noexcept;
87 
91  virtual ~EpetraMatrix() override;
92 
94 
99 
102  using MatrixBase::create;
103  using MatrixBase::closed;
104  using MatrixBase::assembled;
107  using MatrixBase::ready;
108  using MatrixBase::residual;
111  using MatrixBase::extract;
115 
116  virtual void createWithLocalSize( localIndex const localRows,
117  localIndex const localCols,
118  localIndex const maxEntriesPerRow,
119  MPI_Comm const & comm ) override;
120 
121  virtual void createWithGlobalSize( globalIndex const globalRows,
122  globalIndex const globalCols,
123  localIndex const maxEntriesPerRow,
124  MPI_Comm const & comm ) override;
125 
126  virtual void open() override;
127 
128  virtual void close() override;
129 
133  virtual bool created() const override;
134 
135  virtual void reset() override;
136 
137  virtual void set( real64 const value ) override;
138 
139  virtual void zero() override;
140 
141  virtual void add( globalIndex const rowIndex,
142  globalIndex const colIndex,
143  real64 const value ) override;
144 
145  virtual void set( globalIndex const rowIndex,
146  globalIndex const colIndex,
147  real64 const value ) override;
148 
149  virtual void insert( globalIndex const rowIndex,
150  globalIndex const colIndex,
151  real64 const value ) override;
152 
153  virtual void add( globalIndex const rowIndex,
154  globalIndex const * colIndices,
155  real64 const * values,
156  localIndex const size ) override;
157 
158  virtual void set( globalIndex const rowIndex,
159  globalIndex const * colIndices,
160  real64 const * values,
161  localIndex const size ) override;
162 
163  virtual void insert( globalIndex const rowIndex,
164  globalIndex const * colIndices,
165  real64 const * values,
166  localIndex const size ) override;
167 
168  virtual void add( globalIndex const rowIndex,
169  arraySlice1d< globalIndex const > const & colIndices,
170  arraySlice1d< real64 const > const & values ) override;
171 
172  virtual void set( globalIndex const rowIndex,
173  arraySlice1d< globalIndex const > const & colIndices,
174  arraySlice1d< real64 const > const & values ) override;
175 
176  virtual void insert( globalIndex const rowIndex,
177  arraySlice1d< globalIndex const > const & colIndices,
178  arraySlice1d< real64 const > const & values ) override;
179 
180  virtual void add( arraySlice1d< globalIndex const > const & rowIndices,
181  arraySlice1d< globalIndex const > const & colIndices,
182  arraySlice2d< real64 const > const & values ) override;
183 
184  virtual void set( arraySlice1d< globalIndex const > const & rowIndices,
185  arraySlice1d< globalIndex const > const & colIndices,
186  arraySlice2d< real64 const > const & values ) override;
187 
188  virtual void insert( arraySlice1d< globalIndex const > const & rowIndices,
189  arraySlice1d< globalIndex const > const & colIndices,
190  arraySlice2d< real64 const > const & values ) override;
191 
192  virtual void add( globalIndex const * rowIndices,
193  globalIndex const * colIndices,
194  real64 const * values,
195  localIndex const numRows,
196  localIndex const numCols ) override;
197 
198  virtual void set( globalIndex const * rowIndices,
199  globalIndex const * colIndices,
200  real64 const * values,
201  localIndex const numRows,
202  localIndex const numCols ) override;
203 
204  virtual void insert( globalIndex const * rowIndices,
205  globalIndex const * colIndices,
206  real64 const * values,
207  localIndex const numRows,
208  localIndex const numCols ) override;
209 
210  virtual void insert( arrayView1d< globalIndex const > const & rowIndices,
211  arrayView1d< globalIndex const > const & colIndices,
212  arrayView1d< real64 const > const & values ) override;
213 
214  virtual void apply( EpetraVector const & src,
215  EpetraVector & dst ) const override;
216 
217  virtual void applyTranspose( EpetraVector const & src,
218  EpetraVector & dst ) const override;
219 
220  virtual void multiply( EpetraMatrix const & src,
221  EpetraMatrix & dst ) const override;
222 
223  virtual void leftMultiplyTranspose( EpetraMatrix const & src,
224  EpetraMatrix & dst ) const override;
225 
226  virtual void rightMultiplyTranspose( EpetraMatrix const & src,
227  EpetraMatrix & dst ) const override;
228 
229  virtual void multiplyRAP( EpetraMatrix const & R,
230  EpetraMatrix const & P,
231  EpetraMatrix & dst ) const override;
232 
233  virtual void multiplyPtAP( EpetraMatrix const & P,
234  EpetraMatrix & dst ) const override;
235 
236  virtual void gemv( real64 const alpha,
237  EpetraVector const & x,
238  real64 const beta,
239  EpetraVector & y,
240  bool useTranspose = false ) const override;
241 
242  virtual void scale( real64 const scalingFactor ) override;
243 
244  virtual void leftScale( EpetraVector const & vec ) override;
245 
246  virtual void rightScale( EpetraVector const & vec ) override;
247 
248  virtual void leftRightScale( EpetraVector const & vecLeft,
249  EpetraVector const & vecRight ) override;
250 
251  virtual void computeScalingVector( EpetraVector & scaling ) const override;
252 
253  virtual void rescaleRows( arrayView1d< globalIndex const > const & rowIndices,
254  RowSumType const rowSumType ) override;
255 
256  virtual void transpose( EpetraMatrix & dst ) const override;
257 
259  integer const dofsPerNode ) const override;
260 
261  virtual real64 clearRow( globalIndex const row,
262  bool const keepDiag = false,
263  real64 const diagValue = 0.0 ) override;
264 
265  virtual void addEntries( EpetraMatrix const & src,
266  MatrixPatternOp const op,
267  real64 const scale ) override;
268 
269  virtual void addDiagonal( EpetraVector const & src,
270  real64 const scale ) override;
271 
272  virtual void clampEntries( real64 const lo,
273  real64 const hi,
274  bool const excludeDiag ) override;
275 
279  virtual localIndex maxRowLengthLocal() const override;
280 
284  virtual localIndex maxRowLength() const override;
285 
286  virtual localIndex rowLength( globalIndex const globalRowIndex ) const override;
287 
288  virtual void getRowLengths( arrayView1d< localIndex > const & lengths ) const override;
289 
290  virtual void getRowLocalLengths( arrayView1d< localIndex > const & lengths ) const override;
291 
292  virtual void getRowCopy( globalIndex globalRow,
293  arraySlice1d< globalIndex > const & colIndices,
294  arraySlice1d< real64 > const & values ) const override;
295 
296  virtual void extractDiagonal( EpetraVector & dst ) const override;
297 
298  virtual void extract( CRSMatrixView< real64, globalIndex > const & localMat ) const override;
299 
300  virtual void extract( CRSMatrixView< real64, globalIndex const > const & localMat ) const override;
301 
302  virtual void extractLocal( CRSMatrixView< real64, localIndex > const & localMat ) const override;
303 
304  virtual void getRowSums( EpetraVector & dst,
305  RowSumType const rowSumType ) const override;
306 
310  virtual globalIndex numGlobalRows() const override;
311 
315  virtual globalIndex numGlobalCols() const override;
316 
320  virtual localIndex numLocalRows() const override;
321 
325  virtual localIndex numLocalCols() const override;
326 
330  virtual globalIndex ilower() const override;
331 
335  virtual globalIndex iupper() const override;
336 
340  virtual globalIndex jlower() const override;
341 
345  virtual globalIndex jupper() const override;
346 
350  virtual localIndex numLocalNonzeros() const override;
351 
355  virtual globalIndex numGlobalNonzeros() const override;
356 
360  virtual real64 normInf() const override;
361 
365  virtual real64 norm1() const override;
366 
370  virtual real64 normFrobenius() const override;
371 
375  virtual real64 normMax() const override;
376 
377  virtual real64 normMax( arrayView1d< globalIndex const > const & rowIndices ) const override;
378 
379  virtual localIndex getLocalRowID( globalIndex const index ) const override;
380 
381  virtual globalIndex getGlobalRowID( localIndex const index ) const override;
382 
386  virtual MPI_Comm comm() const override;
387 
388  virtual void print( std::ostream & os = std::cout ) const override;
389 
390  virtual void write( string const & filename,
391  LAIOutputFormat const format = LAIOutputFormat::MATRIX_MARKET ) const override;
392 
394 
399  Epetra_FECrsMatrix const & unwrapped() const;
400 
405  Epetra_FECrsMatrix & unwrapped();
406 
407 private:
408 
412  void multiply( bool const transA,
413  EpetraMatrix const & B,
414  bool const transB,
415  EpetraMatrix & C ) const;
416 
421  void create( Epetra_CrsMatrix const & src );
422 
424  std::unique_ptr< Epetra_FECrsMatrix > m_matrix;
425 
427  std::unique_ptr< Epetra_Map > m_src_map;
428 
430  std::unique_ptr< Epetra_Map > m_dst_map;
431 };
432 
433 } // namespace geos
434 
435 #endif /*GEOS_LINEARALGEBRA_INTERFACES_EPETRAMATRIX_HPP_*/
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.
Definition: MatrixBase.hpp:250
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.
Definition: MatrixBase.hpp:76
virtual void multiplyPtAP(Matrix const &P, Matrix &dst) const
Compute the triple product dst = P^T * this * P
Definition: MatrixBase.hpp:648
virtual void residual(Vector const &x, Vector const &b, Vector &r) const override
Compute residual r = b - A * x.
Definition: MatrixBase.hpp:553
void setDofManager(DofManager const *const dofManager)
Associate a DofManager with this matrix.
Definition: MatrixBase.hpp:159
bool assembled() const
Query matrix assembled status.
Definition: MatrixBase.hpp:110
bool closed() const
Query matrix closed status.
Definition: MatrixBase.hpp:104
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:206
bool insertable() const
Query matrix status.
Definition: MatrixBase.hpp:132
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:188
CRSMatrix< real64, globalIndex > extract() const
Extract local part of the matrix.
Definition: MatrixBase.hpp:875
DofManager const * dofManager() const
Definition: MatrixBase.hpp:167
bool ready() const
Query matrix ready status.
Definition: MatrixBase.hpp:117
bool modifiable() const
Query matrix status.
Definition: MatrixBase.hpp:125
virtual void multiplyP1tAP2(Matrix const &P1, Matrix const &P2, Matrix &dst) const
Compute the triple product dst = P1^T * this * P2
Definition: MatrixBase.hpp:628
CRSMatrix< real64, localIndex > extractLocal() const
Extract block of the matrix corresponding to locally owned columns.
Definition: MatrixBase.hpp:898
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.
Definition: MatrixBase.hpp:250
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:179
RowSumType
Type of row sum to compute.
Definition: MatrixBase.hpp:36
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:87
LAIOutputFormat
Definition: common.hpp:142
ArraySlice< T, 2, USD > arraySlice2d
Alias for 2D array slice.
Definition: DataTypes.hpp:199
double real64
64-bit floating point type.
Definition: DataTypes.hpp:98
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:84
ArraySlice< T, 1, USD > arraySlice1d
Alias for 1D array slice.
Definition: DataTypes.hpp:183
LvArray::CRSMatrixView< T, COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:309
int integer
Signed integer type.
Definition: DataTypes.hpp:81
MatrixPatternOp
Describes relationship between and treatment of nonzero patterns of arguments in matrix functions lik...
Definition: MatrixBase.hpp:48