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 Total, S.A
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 
112  virtual void createWithLocalSize( localIndex const localRows,
113  localIndex const localCols,
114  localIndex const maxEntriesPerRow,
115  MPI_Comm const & comm ) override;
116 
117  virtual void createWithGlobalSize( globalIndex const globalRows,
118  globalIndex const globalCols,
119  localIndex const maxEntriesPerRow,
120  MPI_Comm const & comm ) override;
121 
122  virtual void open() override;
123 
124  virtual void close() override;
125 
129  virtual bool created() const override;
130 
131  virtual void reset() override;
132 
133  virtual void set( real64 const value ) override;
134 
135  virtual void zero() override;
136 
137  virtual void add( globalIndex const rowIndex,
138  globalIndex const colIndex,
139  real64 const value ) override;
140 
141  virtual void set( globalIndex const rowIndex,
142  globalIndex const colIndex,
143  real64 const value ) override;
144 
145  virtual void insert( globalIndex const rowIndex,
146  globalIndex const colIndex,
147  real64 const value ) override;
148 
149  virtual void add( globalIndex const rowIndex,
150  globalIndex const * colIndices,
151  real64 const * values,
152  localIndex const size ) override;
153 
154  virtual void set( globalIndex const rowIndex,
155  globalIndex const * colIndices,
156  real64 const * values,
157  localIndex const size ) override;
158 
159  virtual void insert( globalIndex const rowIndex,
160  globalIndex const * colIndices,
161  real64 const * values,
162  localIndex const size ) override;
163 
164  virtual void add( globalIndex const rowIndex,
165  arraySlice1d< globalIndex const > const & colIndices,
166  arraySlice1d< real64 const > const & values ) override;
167 
168  virtual void set( globalIndex const rowIndex,
169  arraySlice1d< globalIndex const > const & colIndices,
170  arraySlice1d< real64 const > const & values ) override;
171 
172  virtual void insert( globalIndex const rowIndex,
173  arraySlice1d< globalIndex const > const & colIndices,
174  arraySlice1d< real64 const > const & values ) override;
175 
176  virtual void add( arraySlice1d< globalIndex const > const & rowIndices,
177  arraySlice1d< globalIndex const > const & colIndices,
178  arraySlice2d< real64 const > const & values ) override;
179 
180  virtual void set( arraySlice1d< globalIndex const > const & rowIndices,
181  arraySlice1d< globalIndex const > const & colIndices,
182  arraySlice2d< real64 const > const & values ) override;
183 
184  virtual void insert( arraySlice1d< globalIndex const > const & rowIndices,
185  arraySlice1d< globalIndex const > const & colIndices,
186  arraySlice2d< real64 const > const & values ) override;
187 
188  virtual void add( globalIndex const * rowIndices,
189  globalIndex const * colIndices,
190  real64 const * values,
191  localIndex const numRows,
192  localIndex const numCols ) override;
193 
194  virtual void set( globalIndex const * rowIndices,
195  globalIndex const * colIndices,
196  real64 const * values,
197  localIndex const numRows,
198  localIndex const numCols ) override;
199 
200  virtual void insert( globalIndex const * rowIndices,
201  globalIndex const * colIndices,
202  real64 const * values,
203  localIndex const numRows,
204  localIndex const numCols ) override;
205 
206  virtual void insert( arrayView1d< globalIndex const > const & rowIndices,
207  arrayView1d< globalIndex const > const & colIndices,
208  arrayView1d< real64 const > const & values ) override;
209 
210  virtual void apply( EpetraVector const & src,
211  EpetraVector & dst ) const override;
212 
213  virtual void applyTranspose( EpetraVector const & src,
214  EpetraVector & dst ) const override;
215 
216  virtual void multiply( EpetraMatrix const & src,
217  EpetraMatrix & dst ) const override;
218 
219  virtual void leftMultiplyTranspose( EpetraMatrix const & src,
220  EpetraMatrix & dst ) const override;
221 
222  virtual void rightMultiplyTranspose( EpetraMatrix const & src,
223  EpetraMatrix & dst ) const override;
224 
225  virtual void multiplyRAP( EpetraMatrix const & R,
226  EpetraMatrix const & P,
227  EpetraMatrix & dst ) const override;
228 
229  virtual void multiplyPtAP( EpetraMatrix const & P,
230  EpetraMatrix & dst ) const override;
231 
232  virtual void gemv( real64 const alpha,
233  EpetraVector const & x,
234  real64 const beta,
235  EpetraVector & y,
236  bool useTranspose = false ) const override;
237 
238  virtual void scale( real64 const scalingFactor ) override;
239 
240  virtual void leftScale( EpetraVector const & vec ) override;
241 
242  virtual void rightScale( EpetraVector const & vec ) override;
243 
244  virtual void leftRightScale( EpetraVector const & vecLeft,
245  EpetraVector const & vecRight ) override;
246 
247  virtual void rescaleRows( arrayView1d< globalIndex const > const & rowIndices,
248  RowSumType const rowSumType ) override;
249 
250  virtual void transpose( EpetraMatrix & dst ) const override;
251 
253  integer const dofsPerNode ) const override;
254 
255  virtual real64 clearRow( globalIndex const row,
256  bool const keepDiag = false,
257  real64 const diagValue = 0.0 ) override;
258 
259  virtual void addEntries( EpetraMatrix const & src,
260  MatrixPatternOp const op,
261  real64 const scale ) override;
262 
263  virtual void addDiagonal( EpetraVector const & src,
264  real64 const scale ) override;
265 
266  virtual void clampEntries( real64 const lo,
267  real64 const hi,
268  bool const excludeDiag ) override;
269 
273  virtual localIndex maxRowLength() const override;
274 
275  virtual localIndex rowLength( globalIndex const globalRowIndex ) const override;
276 
277  virtual void getRowLengths( arrayView1d< localIndex > const & lengths ) const override;
278 
279  virtual void getRowCopy( globalIndex globalRow,
280  arraySlice1d< globalIndex > const & colIndices,
281  arraySlice1d< real64 > const & values ) const override;
282 
283  virtual void extractDiagonal( EpetraVector & dst ) const override;
284 
285  virtual void getRowSums( EpetraVector & dst,
286  RowSumType const rowSumType ) const override;
287 
291  virtual globalIndex numGlobalRows() const override;
292 
296  virtual globalIndex numGlobalCols() const override;
297 
301  virtual localIndex numLocalRows() const override;
302 
306  virtual localIndex numLocalCols() const override;
307 
311  virtual globalIndex ilower() const override;
312 
316  virtual globalIndex iupper() const override;
317 
321  virtual globalIndex jlower() const override;
322 
326  virtual globalIndex jupper() const override;
327 
331  virtual localIndex numLocalNonzeros() const override;
332 
336  virtual globalIndex numGlobalNonzeros() const override;
337 
341  virtual real64 normInf() const override;
342 
346  virtual real64 norm1() const override;
347 
351  virtual real64 normFrobenius() const override;
352 
356  virtual real64 normMax() const override;
357 
358  virtual real64 normMax( arrayView1d< globalIndex const > const & rowIndices ) const override;
359 
360  virtual localIndex getLocalRowID( globalIndex const index ) const override;
361 
362  virtual globalIndex getGlobalRowID( localIndex const index ) const override;
363 
367  virtual MPI_Comm comm() const override;
368 
369  virtual void print( std::ostream & os = std::cout ) const override;
370 
371  virtual void write( string const & filename,
372  LAIOutputFormat const format = LAIOutputFormat::MATRIX_MARKET ) const override;
373 
375 
380  Epetra_FECrsMatrix const & unwrapped() const;
381 
386  Epetra_FECrsMatrix & unwrapped();
387 
388 private:
389 
393  void multiply( bool const transA,
394  EpetraMatrix const & B,
395  bool const transB,
396  EpetraMatrix & C ) const;
397 
402  void create( Epetra_CrsMatrix const & src );
403 
405  std::unique_ptr< Epetra_FECrsMatrix > m_matrix;
406 
408  std::unique_ptr< Epetra_Map > m_src_map;
409 
411  std::unique_ptr< Epetra_Map > m_dst_map;
412 };
413 
414 } // namespace geos
415 
416 #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
Get the number of global columns.
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 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
Returns the total number of nonzeros in the matrix.
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 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 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 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
Get the number of global rows.
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 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 localIndex numLocalRows() const override
Get the number of local rows.
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
Returns the number of nonzeros in the local portion of the matrix.
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
Get the number of local columns.
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:248
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 residual(Vector const &x, Vector const &b, Vector &r) const override
Compute residual r = b - A * x.
Definition: MatrixBase.hpp:551
void setDofManager(DofManager const *const dofManager)
Associate a DofManager with this matrix.
Definition: MatrixBase.hpp:157
bool assembled() const
Query matrix assembled status.
Definition: MatrixBase.hpp:108
bool closed() const
Query matrix closed status.
Definition: MatrixBase.hpp:102
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:204
bool insertable() const
Query matrix status.
Definition: MatrixBase.hpp:130
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:186
DofManager const * dofManager() const
Definition: MatrixBase.hpp:165
bool ready() const
Query matrix ready status.
Definition: MatrixBase.hpp:115
bool modifiable() const
Query matrix status.
Definition: MatrixBase.hpp:123
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:248
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
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:88
LAIOutputFormat
Definition: common.hpp:142
ArraySlice< T, 2, USD > arraySlice2d
Alias for 2D array slice.
Definition: DataTypes.hpp:200
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:85
ArraySlice< T, 1, USD > arraySlice1d
Alias for 1D array slice.
Definition: DataTypes.hpp:184
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
MatrixPatternOp
Describes relationship between and treatment of nonzero patterns of arguments in matrix functions lik...
Definition: MatrixBase.hpp:48