GEOSX
EpetraMatrix.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_EPETRAMATRIX_HPP_
20 #define GEOSX_LINEARALGEBRA_INTERFACES_EPETRAMATRIX_HPP_
21 
22 #include "common/DataTypes.hpp"
26 
27 class Epetra_Map;
28 
29 class Epetra_CrsMatrix;
30 
31 class Epetra_FECrsMatrix;
32 
33 namespace geosx
34 {
35 
39 class EpetraMatrix final : public virtual LinearOperator< EpetraVector >,
40  private MatrixBase< EpetraMatrix, EpetraVector >
41 {
42 public:
43 
46 
50 
57  EpetraMatrix();
58 
63  EpetraMatrix( EpetraMatrix const & src );
64 
68  virtual ~EpetraMatrix() override;
69 
71 
75 
79  using MatrixBase::create;
80  using MatrixBase::closed;
84  using MatrixBase::ready;
86 
87  virtual void createWithLocalSize( localIndex const localRows,
88  localIndex const localCols,
89  localIndex const maxEntriesPerRow,
90  MPI_Comm const & comm ) override;
91 
92  virtual void createWithGlobalSize( globalIndex const globalRows,
93  globalIndex const globalCols,
94  localIndex const maxEntriesPerRow,
95  MPI_Comm const & comm ) override;
96 
97  virtual void open() override;
98 
99  virtual void close() override;
100 
104  virtual bool created() const override;
105 
106  virtual void reset() override;
107 
108  virtual void set( real64 const value ) override;
109 
110  virtual void zero() override;
111 
112  virtual void add( globalIndex const rowIndex,
113  globalIndex const colIndex,
114  real64 const value ) override;
115 
116  virtual void set( globalIndex const rowIndex,
117  globalIndex const colIndex,
118  real64 const value ) override;
119 
120  virtual void insert( globalIndex const rowIndex,
121  globalIndex const colIndex,
122  real64 const value ) override;
123 
124  virtual void add( globalIndex const rowIndex,
125  globalIndex const * colIndices,
126  real64 const * values,
127  localIndex const size ) override;
128 
129  virtual void set( globalIndex const rowIndex,
130  globalIndex const * colIndices,
131  real64 const * values,
132  localIndex const size ) override;
133 
134  virtual void insert( globalIndex const rowIndex,
135  globalIndex const * colIndices,
136  real64 const * values,
137  localIndex const size ) override;
138 
139  virtual void add( globalIndex const rowIndex,
140  arraySlice1d< globalIndex const > const & colIndices,
141  arraySlice1d< real64 const > const & values ) override;
142 
143  virtual void set( globalIndex const rowIndex,
144  arraySlice1d< globalIndex const > const & colIndices,
145  arraySlice1d< real64 const > const & values ) override;
146 
147  virtual void insert( globalIndex const rowIndex,
148  arraySlice1d< globalIndex const > const & colIndices,
149  arraySlice1d< real64 const > const & values ) override;
150 
151  virtual void add( arraySlice1d< globalIndex const > const & rowIndices,
152  arraySlice1d< globalIndex const > const & colIndices,
154 
155  virtual void set( arraySlice1d< globalIndex const > const & rowIndices,
156  arraySlice1d< globalIndex const > const & colIndices,
158 
159  virtual void insert( arraySlice1d< globalIndex const > const & rowIndices,
160  arraySlice1d< globalIndex const > const & colIndices,
162 
163  virtual void add( arraySlice1d< globalIndex const > const & rowIndices,
164  arraySlice1d< globalIndex const > const & colIndices,
166 
167  virtual void set( arraySlice1d< globalIndex const > const & rowIndices,
168  arraySlice1d< globalIndex const > const & colIndices,
170 
171  virtual void insert( arraySlice1d< globalIndex const > const & rowIndices,
172  arraySlice1d< globalIndex const > const & colIndices,
174 
175  virtual void add( globalIndex const * rowIndices,
176  globalIndex const * colIndices,
177  real64 const * values,
178  localIndex const numRows,
179  localIndex const numCols ) override;
180 
181  virtual void set( globalIndex const * rowIndices,
182  globalIndex const * colIndices,
183  real64 const * values,
184  localIndex const numRows,
185  localIndex const numCols ) override;
186 
187  virtual void insert( globalIndex const * rowIndices,
188  globalIndex const * colIndices,
189  real64 const * values,
190  localIndex const numRows,
191  localIndex const numCols ) override;
192 
193  virtual void apply( EpetraVector const & src,
194  EpetraVector & dst ) const override;
195 
196  virtual void applyTranspose( EpetraVector const & src,
197  EpetraVector & dst ) const override;
198 
199  virtual void multiply( EpetraMatrix const & src,
200  EpetraMatrix & dst ) const override;
201 
202  virtual void leftMultiplyTranspose( EpetraMatrix const & src,
203  EpetraMatrix & dst ) const override;
204 
205  virtual void rightMultiplyTranspose( EpetraMatrix const & src,
206  EpetraMatrix & dst ) const override;
207 
208  virtual void multiplyRAP( EpetraMatrix const & R,
209  EpetraMatrix const & P,
210  EpetraMatrix & dst ) const override;
211 
212  virtual void multiplyPtAP( EpetraMatrix const & P,
213  EpetraMatrix & dst ) const override;
214 
215  virtual void gemv( real64 const alpha,
216  EpetraVector const & x,
217  real64 const beta,
218  EpetraVector & y,
219  bool useTranspose = false ) const override;
220 
221  virtual void scale( real64 const scalingFactor ) override;
222 
223  virtual void leftScale( EpetraVector const & vec ) override;
224 
225  virtual void rightScale( EpetraVector const & vec ) override;
226 
227  virtual void leftRightScale( EpetraVector const & vecLeft,
228  EpetraVector const & vecRight ) override;
229 
230  virtual void transpose( EpetraMatrix & dst ) const override;
231 
232  virtual real64 clearRow( globalIndex const row,
233  bool const keepDiag = false,
234  real64 const diagValue = 0.0 ) override;
235 
236  virtual void addEntries( EpetraMatrix const & src,
237  real64 const scale = 1.0,
238  bool const samePattern = true ) override;
239 
240  virtual void addDiagonal( EpetraVector const & src ) override;
241 
245  virtual localIndex maxRowLength() const override;
246 
247  virtual localIndex localRowLength( localIndex localRowIndex ) const override;
248 
249  virtual localIndex globalRowLength( globalIndex globalRowIndex ) const override;
250 
251  virtual void getRowCopy( globalIndex globalRow,
252  arraySlice1d< globalIndex > const & colIndices,
253  arraySlice1d< real64 > const & values ) const override;
254 
255  virtual real64 getDiagValue( globalIndex globalRow ) const override;
256 
257  virtual void extractDiagonal( EpetraVector & dst ) const override;
258 
262  virtual globalIndex numGlobalRows() const override;
263 
267  virtual globalIndex numGlobalCols() const override;
268 
272  virtual localIndex numLocalRows() const override;
273 
277  virtual localIndex numLocalCols() const override;
278 
282  virtual globalIndex ilower() const override;
283 
287  virtual globalIndex iupper() const override;
288 
292  virtual globalIndex jlower() const override;
293 
297  virtual globalIndex jupper() const override;
298 
302  virtual localIndex numLocalNonzeros() const override;
303 
307  virtual globalIndex numGlobalNonzeros() const override;
308 
312  virtual real64 normInf() const override;
313 
317  virtual real64 norm1() const override;
318 
322  virtual real64 normFrobenius() const override;
323 
324  virtual localIndex getLocalRowID( globalIndex const index ) const override;
325 
326  virtual globalIndex getGlobalRowID( localIndex const index ) const override;
327 
331  virtual MPI_Comm getComm() const override;
332 
333  virtual void print( std::ostream & os = std::cout ) const override;
334 
335  virtual void write( string const & filename,
336  LAIOutputFormat const format = LAIOutputFormat::MATRIX_MARKET ) const override;
337 
339 
344  Epetra_FECrsMatrix const & unwrapped() const;
345 
350  Epetra_FECrsMatrix & unwrapped();
351 
352 private:
353 
357  void multiply( bool const transA,
358  EpetraMatrix const & B,
359  bool const transB,
360  EpetraMatrix & C ) const;
361 
366  void create( Epetra_CrsMatrix const & src );
367 
369  std::unique_ptr< Epetra_FECrsMatrix > m_matrix;
370 
372  std::unique_ptr< Epetra_Map > m_src_map;
373 
375  std::unique_ptr< Epetra_Map > m_dst_map;
376 };
377 
378 } // namespace geosx
379 
380 #endif /*GEOSX_LINEARALGEBRA_INTERFACES_EPETRAMATRIX_HPP_*/
Wrapper class for Epetra&#39;s CrsMatrix.
virtual void applyTranspose(EpetraVector const &src, EpetraVector &dst) const override
Apply transpose of the matrix to a vector.
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 bool created() const override
Query matrix creation status.
virtual localIndex maxRowLength() const override
Returns the number of nonzero entries in the longest row of the matrix.
virtual localIndex localRowLength(localIndex localRowIndex) const override
Get row length via local row index.
long long int globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:128
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 void zero() override
Set all elements to zero.
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 reset() override
Reset the matrix to default state.
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 real64 norm1() const override
Returns the one norm of the matrix.
virtual globalIndex ilower() const override
Returns the index of the first global row owned by that processor.
bool ready() const
Query matrix ready status.
Definition: MatrixBase.hpp:125
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 globalIndex jlower() const override
Returns the index of the first global col owned by that processor.
virtual void leftMultiplyTranspose(EpetraMatrix const &src, EpetraMatrix &dst) const override
Matrix/Matrix transpose multiplication.
virtual void extractDiagonal(EpetraVector &dst) const override
Extract diagonal values into a vector.
virtual void addEntries(EpetraMatrix 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.
double real64
64-bit floating point type.
Definition: DataTypes.hpp:136
virtual real64 getDiagValue(globalIndex globalRow) const override
get diagonal element value on a given row
virtual void add(globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
Add to one element.
virtual void transpose(EpetraMatrix &dst) const override
Matrix transposition.
virtual void rightMultiplyTranspose(EpetraMatrix const &src, EpetraMatrix &dst) const override
Matrix/Matrix transpose multiplication.
virtual localIndex numLocalNonzeros() const override
Returns the number of nonzeros in the local portion of the matrix.
virtual real64 normFrobenius() const override
Returns the Frobenius norm of the matrix.
virtual ~EpetraMatrix() override
Virtual destructor.
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 ...
EpetraMatrix()
Empty matrix constructor.
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 void leftScale(EpetraVector const &vec) override
Pre-multiplies (left) with diagonal matrix consisting of the values in vec.
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 localIndex numLocalRows() const override
Return the local number of columns on each processor.
virtual void multiply(EpetraMatrix const &src, EpetraMatrix &dst) const override
Matrix/Matrix multiplication.
virtual globalIndex numGlobalRows() const override
Returns the number of global rows.
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 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 open() override
Open matrix for adding new entries.
LAIOutputFormat
Definition: common.hpp:139
virtual void rightScale(EpetraVector const &vec) override
Post-multiplies (right) with diagonal matrix consisting of the values in vec.
virtual globalIndex numGlobalCols() const override
Returns the number of global columns.
bool closed() const
Query matrix closed status.
Definition: MatrixBase.hpp:112
virtual MPI_Comm getComm() const override
Get the MPI communicator the matrix was created with.
virtual localIndex globalRowLength(globalIndex globalRowIndex) const override
Get row length via global row index.
virtual void close() override
Assemble and compress the matrix.
virtual void print(std::ostream &os=std::cout) const override
Print the matrix in Trilinos format to a stream.
bool modifiable() const
Query matrix status.
Definition: MatrixBase.hpp:133
virtual void multiplyPtAP(EpetraMatrix const &P, EpetraMatrix &dst) const override
Compute the triple product dst = P^T * this * P
virtual void multiplyRAP(EpetraMatrix const &R, EpetraMatrix const &P, EpetraMatrix &dst) const override
Compute the triple product dst = R * this * P
virtual globalIndex getGlobalRowID(localIndex const index) const override
Map a local row index to global row index.
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 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 numGlobalNonzeros() const override
Returns the total number of nonzeros in the matrix.
virtual void residual(Vector const &x, Vector const &b, Vector &r) const override
Compute residual r = Ax - b.
Definition: MatrixBase.hpp:546
virtual localIndex getLocalRowID(globalIndex const index) const override
Map a global row index to local row index.
virtual globalIndex iupper() const override
Returns index one past the last global row owned by that processor.
virtual void scale(real64 const scalingFactor) override
Multiply all elements by scalingFactor.
virtual void insert(globalIndex const rowIndex, globalIndex const colIndex, real64 const value) override
Insert one element.
This class creates and provides basic support for the Epetra_FEVector vector object type used in Tril...
Epetra_FECrsMatrix const & unwrapped() const
Returns a const pointer to the underlying matrix.
Abstract base class for linear operators.
virtual globalIndex jupper() const override
Returns index one past the last global col owned by that processor.
virtual void apply(EpetraVector const &src, EpetraVector &dst) const override
Apply operator to a vector.
virtual localIndex numLocalCols() const override
Return the local number of columns on each processor.
virtual void addDiagonal(EpetraVector const &src) override
Add entries of a vector to the diagonal of this matrix.
Common base template for all matrix wrapper types.
Definition: MatrixBase.hpp:49