GEOSX
MatrixBase.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 TotalEnergies
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 GEOS_LINEARALGEBRA_INTERFACES_MATRIXBASE_HPP_
20 #define GEOS_LINEARALGEBRA_INTERFACES_MATRIXBASE_HPP_
21 
24 #include "LvArray/src/output.hpp"
25 
26 namespace geos
27 {
28 
29 class DofManager;
30 
34 enum class RowSumType
35 {
36  SumValues,
37  SumAbsValues,
38  SumSqrValues,
39  MaxAbsValues
40 };
41 
46 enum class MatrixPatternOp
47 {
48  Same, // Caller guarantees patterns of arguments are exactly the same
49  Subset, // Caller guarantees pattern of second argument is a subset of the first
50  Restrict, // Restrict pattern of second argument, ignoring any entries that don't exist in first
51  Extend // Extend pattern of the first argument with entries from the second
52 };
53 
73 template< typename MATRIX, typename VECTOR >
74 class MatrixBase : public virtual LinearOperator< VECTOR >
75 {
76 protected:
77 
80 
82  using Matrix = MATRIX;
83 
85  using Vector = VECTOR;
86 
87  using Base::numLocalRows;
88  using Base::numLocalCols;
89  using Base::numGlobalRows;
90  using Base::numGlobalCols;
91 
96 
101  inline bool closed() const { return m_closed; }
102 
107  inline bool assembled() const { return m_assembled; }
108 
114  inline bool ready() const { return closed() && assembled(); }
115 
122  inline bool modifiable() const { return !closed() && assembled(); }
123 
129  inline bool insertable() const { return !closed() && !assembled(); }
130 
135  virtual bool created() const = 0;
136 
138 
140 
156  void setDofManager( DofManager const * const dofManager )
157  {
159  }
160 
164  DofManager const * dofManager() const
165  {
166  return m_dofManager;
167  }
168 
170 
175 
184  virtual void
185  createWithLocalSize( localIndex const localSize,
186  localIndex const maxEntriesPerRow,
187  MPI_Comm const & comm )
188  {
189  createWithLocalSize( localSize, localSize, maxEntriesPerRow, comm );
190  }
191 
202  virtual void
204  localIndex const maxEntriesPerRow,
205  MPI_Comm const & comm )
206  {
207  createWithGlobalSize( globalSize, globalSize, maxEntriesPerRow, comm );
208  }
209 
218  virtual void
219  createWithLocalSize( localIndex const localRows,
220  localIndex const localCols,
221  localIndex const maxEntriesPerRow,
222  MPI_Comm const & comm ) = 0;
223 
232  virtual void
234  globalIndex const globalCols,
235  localIndex const maxEntriesPerRow,
236  MPI_Comm const & comm ) = 0;
237 
247  virtual void create( CRSMatrixView< real64 const, globalIndex const > const & localMatrix,
248  localIndex const numLocalColumns,
249  MPI_Comm const & comm )
250  {
251  localMatrix.move( hostMemorySpace, false );
252 
253  localIndex maxEntriesPerRow = 0;
254  for( localIndex i = 0; i < localMatrix.numRows(); ++i )
255  {
256  maxEntriesPerRow = std::max( maxEntriesPerRow, localMatrix.numNonZeros( i ) );
257  }
258 
259  createWithLocalSize( localMatrix.numRows(),
260  numLocalColumns,
261  maxEntriesPerRow,
262  comm );
263 
264  globalIndex const rankOffset = ilower();
265 
266  open();
267  for( localIndex localRow = 0; localRow < localMatrix.numRows(); ++localRow )
268  {
269  insert( localRow + rankOffset, localMatrix.getColumns( localRow ), localMatrix.getEntries( localRow ) );
270  }
271  close();
272  }
273 
275 
280 
287  virtual void open() = 0;
288 
295  virtual void close() = 0;
296 
300  virtual void reset()
301  {
302  m_assembled = false;
303  m_closed = true;
304  }
305 
307 
312 
317  virtual void set( real64 const value ) = 0;
318 
322  virtual void zero() = 0;
323 
325 
337 
344  virtual void add( globalIndex const rowIndex,
345  globalIndex const colIndex,
346  real64 const value ) = 0;
347 
354  virtual void set( globalIndex const rowIndex,
355  globalIndex const colIndex,
356  real64 const value ) = 0;
357 
364  virtual void insert( globalIndex const rowIndex,
365  globalIndex const colIndex,
366  real64 const value ) = 0;
367 
375  virtual void add( globalIndex const rowIndex,
376  globalIndex const * colIndices,
377  real64 const * values,
378  localIndex const size ) = 0;
379 
387  virtual void set( globalIndex const rowIndex,
388  globalIndex const * colIndices,
389  real64 const * values,
390  localIndex const size ) = 0;
391 
399  virtual void insert( globalIndex const rowIndex,
400  globalIndex const * colIndices,
401  real64 const * values,
402  localIndex const size ) = 0;
403 
410  virtual void add( globalIndex const rowIndex,
411  arraySlice1d< globalIndex const > const & colIndices,
412  arraySlice1d< real64 const > const & values ) = 0;
413 
420  virtual void set( globalIndex const rowIndex,
421  arraySlice1d< globalIndex const > const & colIndices,
422  arraySlice1d< real64 const > const & values ) = 0;
423 
430  virtual void insert( globalIndex const rowIndex,
431  arraySlice1d< globalIndex const > const & colIndices,
432  arraySlice1d< real64 const > const & values ) = 0;
433 
440  virtual void add( arraySlice1d< globalIndex const > const & rowIndices,
441  arraySlice1d< globalIndex const > const & colIndices,
442  arraySlice2d< real64 const > const & values ) = 0;
443 
450  virtual void set( arraySlice1d< globalIndex const > const & rowIndices,
451  arraySlice1d< globalIndex const > const & colIndices,
452  arraySlice2d< real64 const > const & values ) = 0;
453 
460  virtual void insert( arraySlice1d< globalIndex const > const & rowIndices,
461  arraySlice1d< globalIndex const > const & colIndices,
462  arraySlice2d< real64 const > const & values ) = 0;
463 
474  virtual void add( globalIndex const * rowIndices,
475  globalIndex const * colIndices,
476  real64 const * values,
477  localIndex const numRows,
478  localIndex const numCols ) = 0;
479 
490  virtual void set( globalIndex const * rowIndices,
491  globalIndex const * colIndices,
492  real64 const * values,
493  localIndex const numRows,
494  localIndex const numCols ) = 0;
495 
506  virtual void insert( globalIndex const * rowIndices,
507  globalIndex const * colIndices,
508  real64 const * values,
509  localIndex const numRows,
510  localIndex const numCols ) = 0;
511 
519  virtual void insert( arrayView1d< globalIndex const > const & rowIndices,
520  arrayView1d< globalIndex const > const & colIndices,
521  arrayView1d< real64 const > const & values ) = 0;
522 
524 
529 
530  virtual void apply( Vector const & src, Vector & dst ) const override = 0;
531 
537  virtual void applyTranspose( Vector const & src, Vector & dst ) const = 0;
538 
550  virtual void residual( Vector const & x, Vector const & b, Vector & r ) const override
551  {
552  if( &b != &r )
553  {
554  r.copy( b );
555  }
556  gemv( -1.0, x, 1.0, r );
557  }
558 
569  virtual void multiply( Matrix const & src,
570  Matrix & dst ) const = 0;
571 
582  virtual void leftMultiplyTranspose( Matrix const & src,
583  Matrix & dst ) const = 0;
584 
595  virtual void rightMultiplyTranspose( Matrix const & src,
596  Matrix & dst ) const = 0;
597 
604  virtual void multiplyRAP( Matrix const & R,
605  Matrix const & P,
606  Matrix & dst ) const
607  {
608  GEOS_LAI_ASSERT( ready() );
609  GEOS_LAI_ASSERT( R.ready() );
610  GEOS_LAI_ASSERT( P.ready() );
611  GEOS_LAI_ASSERT_EQ( numGlobalRows(), R.numGlobalCols() );
612  GEOS_LAI_ASSERT_EQ( numGlobalCols(), P.numGlobalRows() );
613 
614  Matrix AP;
615  multiply( P, AP );
616  R.multiply( AP, dst );
617  }
618 
624  virtual void multiplyPtAP( Matrix const & P,
625  Matrix & dst ) const
626  {
627  GEOS_LAI_ASSERT( ready() );
628  GEOS_LAI_ASSERT( P.ready() );
629  GEOS_LAI_ASSERT_EQ( numGlobalRows(), P.numGlobalRows() );
630  GEOS_LAI_ASSERT_EQ( numGlobalCols(), P.numGlobalRows() );
631 
632  Matrix AP;
633  multiply( P, AP );
634  P.leftMultiplyTranspose( AP, dst );
635  }
636 
642  virtual void multiplyRARt( Matrix const & R,
643  Matrix & dst ) const
644  {
645  Matrix P;
646  R.transpose( P );
647  multiplyPtAP( P, dst );
648  }
649 
663  virtual void gemv( real64 const alpha,
664  Vector const & x,
665  real64 const beta,
666  Vector & y,
667  bool useTranspose = false ) const = 0;
668 
674  virtual void separateComponentFilter( Matrix & dst,
675  integer const dofPerPoint ) const = 0;
676 
677 
678 
683  virtual void scale( real64 const scalingFactor ) = 0;
684 
689  virtual void leftScale( Vector const & vec ) = 0;
690 
695  virtual void rightScale( Vector const & vec ) = 0;
696 
703  virtual void leftRightScale( Vector const & vecLeft,
704  Vector const & vecRight ) = 0;
705 
711  virtual void rescaleRows( arrayView1d< globalIndex const > const & rowIndices,
712  RowSumType const rowSumType ) = 0;
713 
722  virtual void transpose( Matrix & dst ) const = 0;
723 
733  virtual real64 clearRow( globalIndex const row,
734  bool const keepDiag = false,
735  real64 const diagValue = 0.0 ) = 0;
736 
746  virtual void addEntries( Matrix const & src,
747  MatrixPatternOp const op,
748  real64 const scale ) = 0;
749 
758  virtual void addDiagonal( Vector const & src,
759  real64 const scale ) = 0;
760 
769  virtual void clampEntries( real64 const lo,
770  real64 const hi,
771  bool const excludeDiag ) = 0;
772 
774 
779 
786  virtual localIndex maxRowLength() const = 0;
787 
793  virtual localIndex rowLength( globalIndex const globalRowIndex ) const = 0;
794 
800  virtual void getRowLengths( arrayView1d< localIndex > const & lengths ) const = 0;
801 
808  virtual void getRowCopy( globalIndex const globalRow,
809  arraySlice1d< globalIndex > const & colIndices,
810  arraySlice1d< real64 > const & values ) const = 0;
811 
816  virtual void extractDiagonal( Vector & dst ) const = 0;
817 
823  virtual void getRowSums( Vector & dst, RowSumType const rowSumType ) const = 0;
824 
829  virtual globalIndex ilower() const = 0;
830 
837  virtual globalIndex iupper() const = 0;
838 
848  virtual globalIndex jlower() const = 0;
849 
857  virtual globalIndex jupper() const = 0;
858 
863  virtual localIndex numLocalNonzeros() const = 0;
864 
869  virtual globalIndex numGlobalNonzeros() const = 0;
870 
875  virtual real64 normInf() const = 0;
876 
881  virtual real64 norm1() const = 0;
882 
887  virtual real64 normFrobenius() const = 0;
888 
893  virtual real64 normMax() const = 0;
894 
900  virtual real64 normMax( arrayView1d< globalIndex const > const & rowIndices ) const = 0;
901 
907  virtual localIndex getLocalRowID( globalIndex const index ) const = 0;
908 
914  virtual globalIndex getGlobalRowID( localIndex const index ) const = 0;
915 
917 
922 
927  virtual void print( std::ostream & os = std::cout ) const = 0;
928 
938  virtual void write( string const & filename,
939  LAIOutputFormat const format = LAIOutputFormat::MATRIX_MARKET ) const = 0;
940 
942 
949  friend std::ostream & operator<<( std::ostream & os, Matrix const & matrix )
950  {
951  matrix.print( os );
952  return os;
953  }
954 
956  bool m_closed = true;
957 
959  bool m_assembled = false;
960 
963 
964 };
965 
966 } // namespace geos
967 
968 #endif //GEOS_LINEARALGEBRA_INTERFACES_MATRIXBASE_HPP_
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
Definition: DofManager.hpp:43
Abstract base class for linear operators.
virtual globalIndex numGlobalCols() const =0
Get the number of global columns.
VECTOR Vector
Alias for template parameter.
virtual globalIndex numGlobalRows() const =0
Get the number of global rows.
virtual localIndex numLocalRows() const =0
Get the number of local rows.
virtual localIndex numLocalCols() const =0
Get the number of local columns.
virtual MPI_Comm comm() const =0
Get the MPI communicator the matrix was created with.
Common base template for all matrix wrapper types.
Definition: MatrixBase.hpp:75
virtual void multiply(Matrix const &src, Matrix &dst) const =0
Matrix/Matrix multiplication.
virtual void rightMultiplyTranspose(Matrix const &src, Matrix &dst) const =0
Matrix/Matrix transpose multiplication.
virtual real64 norm1() const =0
Returns the one norm of the matrix.
virtual bool created() const =0
Query matrix creation status.
virtual void multiplyPtAP(Matrix const &P, Matrix &dst) const
Compute the triple product dst = P^T * this * P
Definition: MatrixBase.hpp:624
virtual globalIndex ilower() const =0
Returns the index of the first global row owned by that processor.
virtual void add(globalIndex const *rowIndices, globalIndex const *colIndices, real64 const *values, localIndex const numRows, localIndex const numCols)=0
Add a dense block of values.
virtual real64 normMax(arrayView1d< globalIndex const > const &rowIndices) const =0
Returns the max norm of the matrix on a subset of rows.
virtual void insert(arraySlice1d< globalIndex const > const &rowIndices, arraySlice1d< globalIndex const > const &colIndices, arraySlice2d< real64 const > const &values)=0
Insert a dense block of values.
virtual void residual(Vector const &x, Vector const &b, Vector &r) const override
Compute residual r = b - A * x.
Definition: MatrixBase.hpp:550
virtual void createWithGlobalSize(globalIndex const globalRows, globalIndex const globalCols, localIndex const maxEntriesPerRow, MPI_Comm const &comm)=0
Create a rectangular matrix from number of rows/columns.
virtual void multiplyRAP(Matrix const &R, Matrix const &P, Matrix &dst) const
Compute the triple product dst = R * this * P
Definition: MatrixBase.hpp:604
virtual void set(globalIndex const rowIndex, globalIndex const *colIndices, real64 const *values, localIndex const size)=0
Set elements to one row using c-style arrays.
friend std::ostream & operator<<(std::ostream &os, Matrix const &matrix)
Stream insertion operator for all matrix types.
Definition: MatrixBase.hpp:949
void setDofManager(DofManager const *const dofManager)
Associate a DofManager with this matrix.
Definition: MatrixBase.hpp:156
MATRIX Matrix
Type alias for actual derived matrix class.
Definition: MatrixBase.hpp:82
bool assembled() const
Query matrix assembled status.
Definition: MatrixBase.hpp:107
virtual void write(string const &filename, LAIOutputFormat const format=LAIOutputFormat::MATRIX_MARKET) const =0
Write the matrix to filename in a matlab-compatible format.
virtual void getRowCopy(globalIndex const globalRow, arraySlice1d< globalIndex > const &colIndices, arraySlice1d< real64 > const &values) const =0
Returns a copy of the data in row globalRow.
virtual void transpose(Matrix &dst) const =0
Matrix transposition.
virtual void set(real64 const value)=0
Set all non-zero elements to a value.
virtual void leftMultiplyTranspose(Matrix const &src, Matrix &dst) const =0
Matrix/Matrix transpose multiplication.
virtual localIndex maxRowLength() const =0
Returns the number of nonzero entries in the longest row of the matrix.
virtual void add(globalIndex const rowIndex, globalIndex const colIndex, real64 const value)=0
Add to one element.
virtual localIndex numLocalNonzeros() const =0
Returns the number of nonzeros in the local portion of the matrix.
virtual void leftScale(Vector const &vec)=0
Pre-multiplies (left) with diagonal matrix consisting of the values in vec.
bool closed() const
Query matrix closed status.
Definition: MatrixBase.hpp:101
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:203
virtual globalIndex jlower() const =0
Returns the index of the first global col owned by that processor.
virtual void apply(Vector const &src, Vector &dst) const override=0
Apply operator to a vector, dst = this(src).
virtual void clampEntries(real64 const lo, real64 const hi, bool const excludeDiag)=0
Clamp each matrix value between values of lo and hi.
bool m_closed
Flag indicating whether the matrix is currently open for adding new entries.
Definition: MatrixBase.hpp:956
virtual void rescaleRows(arrayView1d< globalIndex const > const &rowIndices, RowSumType const rowSumType)=0
Rescales selected rows of matrix using row sum reciprocal as a factor.
virtual void separateComponentFilter(Matrix &dst, integer const dofPerPoint) const =0
Apply a separate component approximation (filter) to this matrix.
virtual void leftRightScale(Vector const &vecLeft, Vector const &vecRight)=0
Post-multiplies (right) with diagonal matrix consisting of the values in vecRight and pre-multiplies ...
bool insertable() const
Query matrix status.
Definition: MatrixBase.hpp:129
bool m_assembled
Flag indicating whether the matrix (sparsity pattern) has been assembled.
Definition: MatrixBase.hpp:959
virtual real64 normFrobenius() const =0
Returns the Frobenius norm of the matrix.
virtual void applyTranspose(Vector const &src, Vector &dst) const =0
Apply transpose of the matrix to a vector.
virtual void scale(real64 const scalingFactor)=0
Multiply all elements by scalingFactor.
virtual void getRowSums(Vector &dst, RowSumType const rowSumType) const =0
Populate a vector with row sums of this.
virtual globalIndex getGlobalRowID(localIndex const index) const =0
Map a local row index to global row index.
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:185
virtual globalIndex numGlobalRows() const=0
Get the number of global rows.
virtual real64 normInf() const =0
Returns the infinity norm of the matrix.
virtual void insert(globalIndex const rowIndex, globalIndex const colIndex, real64 const value)=0
Insert one element.
virtual void zero()=0
Set all elements to zero.
virtual void open()=0
Open matrix for adding new entries.
virtual void addDiagonal(Vector const &src, real64 const scale)=0
Add (scaled) entries of a vector to the diagonal of this matrix.
virtual real64 normMax() const =0
Returns the max norm of the matrix (the largest absolute element value).
virtual void insert(globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values)=0
Insert elements of one row using array1d.
virtual void add(globalIndex const rowIndex, globalIndex const *colIndices, real64 const *values, localIndex const size)=0
Add elements to one row using c-style arrays.
virtual void insert(globalIndex const *rowIndices, globalIndex const *colIndices, real64 const *values, localIndex const numRows, localIndex const numCols)=0
Insert dense matrix.
virtual void extractDiagonal(Vector &dst) const =0
Extract diagonal values into a vector.
virtual void insert(globalIndex const rowIndex, globalIndex const *colIndices, real64 const *values, localIndex const size)=0
Insert elements to one row using c-style arrays.
virtual void add(globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values)=0
Add elements to one row using array1d.
DofManager const * m_dofManager
(optional) DofManager associated with this matrix
Definition: MatrixBase.hpp:962
virtual void set(globalIndex const *rowIndices, globalIndex const *colIndices, real64 const *values, localIndex const numRows, localIndex const numCols)=0
Set a dense block of values.
virtual globalIndex numGlobalNonzeros() const =0
Returns the total number of nonzeros in the matrix.
virtual void set(arraySlice1d< globalIndex const > const &rowIndices, arraySlice1d< globalIndex const > const &colIndices, arraySlice2d< real64 const > const &values)=0
Set a dense block of values.
virtual void print(std::ostream &os=std::cout) const =0
Print the matrix in Trilinos format to a stream.
virtual void multiplyRARt(Matrix const &R, Matrix &dst) const
Compute the triple product dst = R * this * R^T
Definition: MatrixBase.hpp:642
virtual void getRowLengths(arrayView1d< localIndex > const &lengths) const =0
Get the row lengths of every local row.
virtual void add(arraySlice1d< globalIndex const > const &rowIndices, arraySlice1d< globalIndex const > const &colIndices, arraySlice2d< real64 const > const &values)=0
Add a dense block of values.
virtual globalIndex iupper() const =0
Returns index one past the last global row owned by that processor.
DofManager const * dofManager() const
Definition: MatrixBase.hpp:164
virtual void rightScale(Vector const &vec)=0
Post-multiplies (right) with diagonal matrix consisting of the values in vec.
virtual globalIndex numGlobalCols() const=0
Get the number of global columns.
bool ready() const
Query matrix ready status.
Definition: MatrixBase.hpp:114
virtual void set(globalIndex const rowIndex, arraySlice1d< globalIndex const > const &colIndices, arraySlice1d< real64 const > const &values)=0
Set elements of one row using array1d.
virtual real64 clearRow(globalIndex const row, bool const keepDiag=false, real64 const diagValue=0.0)=0
Clear a row, and optionally set diagonal element to diagValue.
bool modifiable() const
Query matrix status.
Definition: MatrixBase.hpp:122
virtual localIndex getLocalRowID(globalIndex const index) const =0
Map a global row index to local row index.
virtual void reset()
Reset the matrix to default state.
Definition: MatrixBase.hpp:300
virtual void gemv(real64 const alpha, Vector const &x, real64 const beta, Vector &y, bool useTranspose=false) const =0
Compute gemv y = alpha*A*x + beta*y.
virtual localIndex rowLength(globalIndex const globalRowIndex) const =0
Get row length via global row index.
virtual void close()=0
Assemble and compress the matrix.
virtual void createWithLocalSize(localIndex const localRows, localIndex const localCols, localIndex const maxEntriesPerRow, MPI_Comm const &comm)=0
Create a rectangular matrix from number of rows/columns.
virtual void addEntries(Matrix const &src, MatrixPatternOp const op, real64 const scale)=0
Add entries of another matrix to this.
virtual globalIndex jupper() const =0
Returns index one past the last global col owned by that processor.
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:247
virtual void set(globalIndex const rowIndex, globalIndex const colIndex, real64 const value)=0
Set one element.
virtual void insert(arrayView1d< globalIndex const > const &rowIndices, arrayView1d< globalIndex const > const &colIndices, arrayView1d< real64 const > const &values)=0
Insert values stored in 3 linear vectors.
#define GEOS_LAI_ASSERT_EQ(lhs, rhs)
Definition: common.hpp:47
#define GEOS_LAI_ASSERT(expr)
Definition: common.hpp:33
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:220
RowSumType
Type of row sum to compute.
Definition: MatrixBase.hpp:35
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:350
LAIOutputFormat
Definition: common.hpp:140
ArraySlice< T, 2, USD > arraySlice2d
Alias for 2D array slice.
Definition: DataTypes.hpp:240
double real64
64-bit floating point type.
Definition: DataTypes.hpp:139
ArraySlice< T, 1, USD > arraySlice1d
Alias for 1D array slice.
Definition: DataTypes.hpp:224
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:122
GEOSX_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:128
GEOSX_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:125
MatrixPatternOp
Describes relationship between and treatment of nonzero patterns of arguments in matrix functions lik...
Definition: MatrixBase.hpp:47