GEOS
MatrixBase.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_MATRIXBASE_HPP_
21 #define GEOS_LINEARALGEBRA_INTERFACES_MATRIXBASE_HPP_
22 
25 #include "LvArray/src/output.hpp"
26 
27 namespace geos
28 {
29 
30 class DofManager;
31 
35 enum class RowSumType
36 {
37  SumValues,
38  SumAbsValues,
39  SumSqrValues,
40  MaxAbsValues
41 };
42 
47 enum class MatrixPatternOp
48 {
49  Same, // Caller guarantees patterns of arguments are exactly the same
50  Subset, // Caller guarantees pattern of second argument is a subset of the first
51  Restrict, // Restrict pattern of second argument, ignoring any entries that don't exist in first
52  Extend // Extend pattern of the first argument with entries from the second
53 };
54 
74 template< typename MATRIX, typename VECTOR >
75 class MatrixBase : public virtual LinearOperator< VECTOR >
76 {
77 protected:
78 
81 
83  using Matrix = MATRIX;
84 
86  using Vector = VECTOR;
87 
88  using Base::numLocalRows;
89  using Base::numLocalCols;
90  using Base::numGlobalRows;
91  using Base::numGlobalCols;
92 
97 
102  inline bool closed() const { return m_closed; }
103 
108  inline bool assembled() const { return m_assembled; }
109 
115  inline bool ready() const { return closed() && assembled(); }
116 
123  inline bool modifiable() const { return !closed() && assembled(); }
124 
130  inline bool insertable() const { return !closed() && !assembled(); }
131 
136  virtual bool created() const = 0;
137 
139 
141 
157  void setDofManager( DofManager const * const dofManager )
158  {
160  }
161 
165  DofManager const * dofManager() const
166  {
167  return m_dofManager;
168  }
169 
171 
176 
185  virtual void
186  createWithLocalSize( localIndex const localSize,
187  localIndex const maxEntriesPerRow,
188  MPI_Comm const & comm )
189  {
190  createWithLocalSize( localSize, localSize, maxEntriesPerRow, comm );
191  }
192 
203  virtual void
205  localIndex const maxEntriesPerRow,
206  MPI_Comm const & comm )
207  {
208  createWithGlobalSize( globalSize, globalSize, maxEntriesPerRow, comm );
209  }
210 
219  virtual void
220  createWithLocalSize( localIndex const localRows,
221  localIndex const localCols,
222  localIndex const maxEntriesPerRow,
223  MPI_Comm const & comm ) = 0;
224 
233  virtual void
235  globalIndex const globalCols,
236  localIndex const maxEntriesPerRow,
237  MPI_Comm const & comm ) = 0;
238 
248  virtual void create( CRSMatrixView< real64 const, globalIndex const > const & localMatrix,
249  localIndex const numLocalColumns,
250  MPI_Comm const & comm )
251  {
252  localMatrix.move( hostMemorySpace, false );
253 
254  localIndex maxEntriesPerRow = 0;
255  for( localIndex i = 0; i < localMatrix.numRows(); ++i )
256  {
257  maxEntriesPerRow = std::max( maxEntriesPerRow, localMatrix.numNonZeros( i ) );
258  }
259 
260  createWithLocalSize( localMatrix.numRows(),
261  numLocalColumns,
262  maxEntriesPerRow,
263  comm );
264 
265  globalIndex const rankOffset = ilower();
266 
267  open();
268  for( localIndex localRow = 0; localRow < localMatrix.numRows(); ++localRow )
269  {
270  insert( localRow + rankOffset, localMatrix.getColumns( localRow ), localMatrix.getEntries( localRow ) );
271  }
272  close();
273  }
274 
276 
281 
288  virtual void open() = 0;
289 
296  virtual void close() = 0;
297 
301  virtual void reset()
302  {
303  m_assembled = false;
304  m_closed = true;
305  }
306 
308 
313 
318  virtual void set( real64 const value ) = 0;
319 
323  virtual void zero() = 0;
324 
326 
338 
345  virtual void add( globalIndex const rowIndex,
346  globalIndex const colIndex,
347  real64 const value ) = 0;
348 
355  virtual void set( globalIndex const rowIndex,
356  globalIndex const colIndex,
357  real64 const value ) = 0;
358 
365  virtual void insert( globalIndex const rowIndex,
366  globalIndex const colIndex,
367  real64 const value ) = 0;
368 
376  virtual void add( globalIndex const rowIndex,
377  globalIndex const * colIndices,
378  real64 const * values,
379  localIndex const size ) = 0;
380 
388  virtual void set( globalIndex const rowIndex,
389  globalIndex const * colIndices,
390  real64 const * values,
391  localIndex const size ) = 0;
392 
400  virtual void insert( globalIndex const rowIndex,
401  globalIndex const * colIndices,
402  real64 const * values,
403  localIndex const size ) = 0;
404 
411  virtual void add( globalIndex const rowIndex,
412  arraySlice1d< globalIndex const > const & colIndices,
413  arraySlice1d< real64 const > const & values ) = 0;
414 
421  virtual void set( globalIndex const rowIndex,
422  arraySlice1d< globalIndex const > const & colIndices,
423  arraySlice1d< real64 const > const & values ) = 0;
424 
431  virtual void insert( globalIndex const rowIndex,
432  arraySlice1d< globalIndex const > const & colIndices,
433  arraySlice1d< real64 const > const & values ) = 0;
434 
441  virtual void add( arraySlice1d< globalIndex const > const & rowIndices,
442  arraySlice1d< globalIndex const > const & colIndices,
443  arraySlice2d< real64 const > const & values ) = 0;
444 
451  virtual void set( arraySlice1d< globalIndex const > const & rowIndices,
452  arraySlice1d< globalIndex const > const & colIndices,
453  arraySlice2d< real64 const > const & values ) = 0;
454 
461  virtual void insert( arraySlice1d< globalIndex const > const & rowIndices,
462  arraySlice1d< globalIndex const > const & colIndices,
463  arraySlice2d< real64 const > const & values ) = 0;
464 
475  virtual void add( globalIndex const * rowIndices,
476  globalIndex const * colIndices,
477  real64 const * values,
478  localIndex const numRows,
479  localIndex const numCols ) = 0;
480 
491  virtual void set( globalIndex const * rowIndices,
492  globalIndex const * colIndices,
493  real64 const * values,
494  localIndex const numRows,
495  localIndex const numCols ) = 0;
496 
507  virtual void insert( globalIndex const * rowIndices,
508  globalIndex const * colIndices,
509  real64 const * values,
510  localIndex const numRows,
511  localIndex const numCols ) = 0;
512 
520  virtual void insert( arrayView1d< globalIndex const > const & rowIndices,
521  arrayView1d< globalIndex const > const & colIndices,
522  arrayView1d< real64 const > const & values ) = 0;
523 
525 
530 
531  virtual void apply( Vector const & src, Vector & dst ) const override = 0;
532 
538  virtual void applyTranspose( Vector const & src, Vector & dst ) const = 0;
539 
551  virtual void residual( Vector const & x, Vector const & b, Vector & r ) const override
552  {
553  if( &b != &r )
554  {
555  r.copy( b );
556  }
557  gemv( -1.0, x, 1.0, r );
558  }
559 
570  virtual void multiply( Matrix const & src,
571  Matrix & dst ) const = 0;
572 
583  virtual void leftMultiplyTranspose( Matrix const & src,
584  Matrix & dst ) const = 0;
585 
596  virtual void rightMultiplyTranspose( Matrix const & src,
597  Matrix & dst ) const = 0;
598 
605  virtual void multiplyRAP( Matrix const & R,
606  Matrix const & P,
607  Matrix & dst ) const
608  {
609  GEOS_LAI_ASSERT( ready() );
610  GEOS_LAI_ASSERT( R.ready() );
611  GEOS_LAI_ASSERT( P.ready() );
612  GEOS_LAI_ASSERT_EQ( numGlobalRows(), R.numGlobalCols() );
613  GEOS_LAI_ASSERT_EQ( numGlobalCols(), P.numGlobalRows() );
614 
615  Matrix AP;
616  multiply( P, AP );
617  R.multiply( AP, dst );
618  }
619 
625  virtual void multiplyPtAP( Matrix const & P,
626  Matrix & dst ) const
627  {
628  GEOS_LAI_ASSERT( ready() );
629  GEOS_LAI_ASSERT( P.ready() );
630  GEOS_LAI_ASSERT_EQ( numGlobalRows(), P.numGlobalRows() );
631  GEOS_LAI_ASSERT_EQ( numGlobalCols(), P.numGlobalRows() );
632 
633  Matrix AP;
634  multiply( P, AP );
635  P.leftMultiplyTranspose( AP, dst );
636  }
637 
643  virtual void multiplyRARt( Matrix const & R,
644  Matrix & dst ) const
645  {
646  Matrix P;
647  R.transpose( P );
648  multiplyPtAP( P, dst );
649  }
650 
664  virtual void gemv( real64 const alpha,
665  Vector const & x,
666  real64 const beta,
667  Vector & y,
668  bool useTranspose = false ) const = 0;
669 
675  virtual void separateComponentFilter( Matrix & dst,
676  integer const dofPerPoint ) const = 0;
677 
678 
679 
684  virtual void scale( real64 const scalingFactor ) = 0;
685 
690  virtual void leftScale( Vector const & vec ) = 0;
691 
696  virtual void rightScale( Vector const & vec ) = 0;
697 
704  virtual void leftRightScale( Vector const & vecLeft,
705  Vector const & vecRight ) = 0;
706 
712  virtual void rescaleRows( arrayView1d< globalIndex const > const & rowIndices,
713  RowSumType const rowSumType ) = 0;
714 
723  virtual void transpose( Matrix & dst ) const = 0;
724 
734  virtual real64 clearRow( globalIndex const row,
735  bool const keepDiag = false,
736  real64 const diagValue = 0.0 ) = 0;
737 
747  virtual void addEntries( Matrix const & src,
748  MatrixPatternOp const op,
749  real64 const scale ) = 0;
750 
759  virtual void addDiagonal( Vector const & src,
760  real64 const scale ) = 0;
761 
770  virtual void clampEntries( real64 const lo,
771  real64 const hi,
772  bool const excludeDiag ) = 0;
773 
775 
780 
787  virtual localIndex maxRowLength() const = 0;
788 
794  virtual localIndex rowLength( globalIndex const globalRowIndex ) const = 0;
795 
801  virtual void getRowLengths( arrayView1d< localIndex > const & lengths ) const = 0;
802 
809  virtual void getRowCopy( globalIndex const globalRow,
810  arraySlice1d< globalIndex > const & colIndices,
811  arraySlice1d< real64 > const & values ) const = 0;
812 
817  virtual void extractDiagonal( Vector & dst ) const = 0;
818 
824  virtual void getRowSums( Vector & dst, RowSumType const rowSumType ) const = 0;
825 
830  virtual globalIndex ilower() const = 0;
831 
838  virtual globalIndex iupper() const = 0;
839 
849  virtual globalIndex jlower() const = 0;
850 
858  virtual globalIndex jupper() const = 0;
859 
864  virtual localIndex numLocalNonzeros() const = 0;
865 
870  virtual globalIndex numGlobalNonzeros() const = 0;
871 
876  virtual real64 normInf() const = 0;
877 
882  virtual real64 norm1() const = 0;
883 
888  virtual real64 normFrobenius() const = 0;
889 
894  virtual real64 normMax() const = 0;
895 
901  virtual real64 normMax( arrayView1d< globalIndex const > const & rowIndices ) const = 0;
902 
908  virtual localIndex getLocalRowID( globalIndex const index ) const = 0;
909 
915  virtual globalIndex getGlobalRowID( localIndex const index ) const = 0;
916 
918 
923 
928  virtual void print( std::ostream & os = std::cout ) const = 0;
929 
939  virtual void write( string const & filename,
940  LAIOutputFormat const format = LAIOutputFormat::MATRIX_MARKET ) const = 0;
941 
943 
950  friend std::ostream & operator<<( std::ostream & os, Matrix const & matrix )
951  {
952  matrix.print( os );
953  return os;
954  }
955 
957  bool m_closed = true;
958 
960  bool m_assembled = false;
961 
964 
965 };
966 
967 } // namespace geos
968 
969 #endif //GEOS_LINEARALGEBRA_INTERFACES_MATRIXBASE_HPP_
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
Definition: DofManager.hpp:44
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:76
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:625
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:551
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:605
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:950
void setDofManager(DofManager const *const dofManager)
Associate a DofManager with this matrix.
Definition: MatrixBase.hpp:157
MATRIX Matrix
Type alias for actual derived matrix class.
Definition: MatrixBase.hpp:83
bool assembled() const
Query matrix assembled status.
Definition: MatrixBase.hpp:108
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: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
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:957
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:130
bool m_assembled
Flag indicating whether the matrix (sparsity pattern) has been assembled.
Definition: MatrixBase.hpp:960
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:186
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:963
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:643
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:165
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:115
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:123
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:301
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:248
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:49
#define GEOS_LAI_ASSERT(expr)
Definition: common.hpp:35
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
RowSumType
Type of row sum to compute.
Definition: MatrixBase.hpp:36
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:310
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