GEOS
LinearSolverParameters.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_UTILITIES_LINEARSOLVERPARAMETERS_HPP_
21 #define GEOS_LINEARALGEBRA_UTILITIES_LINEARSOLVERPARAMETERS_HPP_
22 
24 #include "common/DataTypes.hpp"
25 
26 namespace geos
27 {
28 
36 {
40  enum class SolverType : integer
41  {
42  direct,
43  cg,
44  gmres,
45  fgmres,
46  bicgstab,
47  richardson,
49  };
50 
55  {
56  none,
57  jacobi,
58  l1jacobi,
59  fgs,
60  sgs,
61  l1sgs,
62  chebyshev,
63  iluk,
64  ilut,
65  ick,
66  ict,
67  amg,
68  mgr,
69  block,
70  direct,
71  bgs,
72  multiscale
73  };
74 
77  bool isSymmetric = false;
79 
82 
84  struct Direct
85  {
89  enum class ColPerm : integer
90  {
91  none,
92  MMD_AtplusA,
93  MMD_AtA,
94  colAMD,
95  metis,
96  parmetis
97  };
98 
102  enum class RowPerm : integer
103  {
104  none,
105  mc64
106  };
107 
116  }
118 
120  struct Krylov
121  {
124 #if GEOS_USE_HYPRE_DEVICE == GEOS_USE_HYPRE_CUDA || GEOS_USE_HYPRE_DEVICE == GEOS_USE_HYPRE_HIP
126 #else
127  integer maxRestart = 200;
128 #endif
134  }
136 
138  struct Relaxation
139  {
140  real64 weight = 2.0 / 3.0;
141  }
143 
145  struct Chebyshev
146  {
149  }
151 
153  struct Scaling
154  {
157  }
159 
161  struct AMG
162  {
164  enum class CycleType : integer
165  {
166  V,
167  W
168  };
169 
171  enum class PreOrPost : integer
172  {
173  pre,
174  post,
175  both
176  };
177 
179  enum class SmootherType : integer
180  {
181  default_,
182  jacobi,
183  l1jacobi,
184  fgs,
185  bgs,
186  sgs,
187  l1sgs,
188  chebyshev,
189  iluk,
190  ilut,
191  ick,
192  ict
193  };
194 
196  enum class CoarseType : integer
197  {
198  default_,
199  jacobi,
200  l1jacobi,
201  fgs,
202  sgs,
203  l1sgs,
204  chebyshev,
205  direct,
206  bgs,
209  };
210 
212  enum class CoarseningType : integer
213  {
214  default_,
215  CLJP,
216  RugeStueben,
217  Falgout,
218  PMIS,
219  HMIS
220  };
221 
223  enum class InterpType : integer
224  {
225  default_,
227  direct,
228  multipass,
229  extendedI,
230  standard,
231  extended,
232  directBAMG,
236  };
237 
239  enum class AggInterpType : integer
240  {
241  default_,
245  multipass,
250  };
251 
253  enum class NullSpaceType : integer
254  {
255  constantModes,
257  };
258 
259 #if GEOS_USE_HYPRE_DEVICE == GEOS_USE_HYPRE_CUDA || GEOS_USE_HYPRE_DEVICE == GEOS_USE_HYPRE_HIP
262 #else
265 #endif
266 
286  }
287  amg;
288 
290  struct MGR
291  {
295  enum class StrategyType : integer
296  {
297  invalid,
318  hydrofracture,
323  };
324 
329  }
330  mgr;
331 
333  struct IFact
334  {
335  integer fill = 0;
337  }
339 
341  struct DD
342  {
344  }
345  dd;
346 
348  struct Block
349  {
351  enum class Shape
352  {
353  Diagonal,
357  };
358 
360  enum class SchurType
361  {
362  None,
366  };
367 
369  enum class Scaling
370  {
371  None,
372  FrobeniusNorm,
373  UserProvided
374  };
375 
379 
382 
387  void resize( integer const numBlocks )
388  {
389  subParams.resize( numBlocks );
390  order.resize( numBlocks );
391  for( integer i = 0; i < numBlocks; ++i )
392  {
393  order[i] = i;
394  }
395  }
396  }
398 
400  struct Multiscale
401  {
403  enum class BasisType
404  {
405  msrsb,
406  };
407 
408  // Core parameters
410  string fieldName;
413  real64 droptol = 0.0;
417 
418  // Debugging/user-display settings
419  string label;
421 
423  struct Smoother
424  {
428  }
430 
432  struct Coupled
433  {
435  }
437 
439  struct Coarsening
440  {
442  enum class PartitionType
443  {
444  graph,
445  cartesian,
447  };
448 
452 
454  struct Structured
455  {
457  }
459 
461  struct Graph
462  {
464  enum class Method
465  {
466  metis,
467  scotch,
468  };
469 
470  Method method = Method::metis;
475 
477  struct Metis
478  {
480  enum class Method
481  {
482  kway,
483  recursive
484  };
485 
486  Method method = Method::kway;
488  integer seed = 2020;
489  }
491 
493  struct Scotch
494  {
495  // TODO
496  } scotch;
497 
498  } graph;
499 
501 
503  struct MsRSB
504  {
506  enum class SupportType : integer
507  {
508  layers,
509  matching
510  };
511 
512  SupportType support = SupportType::matching;
513 
515  integer maxIter = 100;
516  real64 tolerance = 1e-3;
517  real64 relaxation = 2.0 / 3.0;
520  }
522  }
524 };
525 
528  "direct",
529  "cg",
530  "gmres",
531  "fgmres",
532  "bicgstab",
533  "richardson",
534  "preconditioner" );
535 
538  "none",
539  "jacobi",
540  "l1jacobi",
541  "fgs",
542  "sgs",
543  "l1sgs",
544  "chebyshev",
545  "iluk",
546  "ilut",
547  "ick",
548  "ict",
549  "amg",
550  "mgr",
551  "block",
552  "direct",
553  "bgs",
554  "multiscale" );
555 
558  "none",
559  "MMD_AtplusA",
560  "MMD_AtA",
561  "colAMD",
562  "metis",
563  "parmetis" );
564 
567  "none",
568  "mc64" );
569 
572  "invalid",
573  "singlePhaseReservoirFVM",
574  "singlePhaseHybridFVM",
575  "singlePhaseReservoirHybridFVM",
576  "singlePhasePoromechanics",
577  "thermalSinglePhasePoromechanics",
578  "hybridSinglePhasePoromechanics",
579  "singlePhasePoromechanicsEmbeddedFractures",
580  "singlePhasePoromechanicsConformingFractures",
581  "singlePhasePoromechanicsReservoirFVM",
582  "compositionalMultiphaseFVM",
583  "compositionalMultiphaseHybridFVM",
584  "compositionalMultiphaseReservoirFVM",
585  "compositionalMultiphaseReservoirHybridFVM",
586  "immiscibleMultiphaseFVM",
587  "reactiveCompositionalMultiphaseOBL",
588  "thermalCompositionalMultiphaseFVM",
589  "thermalCompositionalMultiphaseReservoirFVM",
590  "multiphasePoromechanics",
591  "multiphasePoromechanicsReservoirFVM",
592  "thermalMultiphasePoromechanics",
593  "hydrofracture",
594  "lagrangianContactMechanics",
595  "augmentedLagrangianContactMechanics",
596  "lagrangianContactMechanicsBubbleStab",
597  "solidMechanicsEmbeddedFractures" );
598 
601  "V",
602  "W" );
603 
606  "pre",
607  "post",
608  "both" );
609 
612  "default",
613  "jacobi",
614  "l1jacobi",
615  "fgs",
616  "bgs",
617  "sgs",
618  "l1sgs",
619  "chebyshev",
620  "iluk",
621  "ilut",
622  "ick",
623  "ict" );
624 
627  "default",
628  "jacobi",
629  "l1jacobi",
630  "fgs",
631  "sgs",
632  "l1sgs",
633  "chebyshev",
634  "direct",
635  "bgs",
636  "gsElimWPivoting",
637  "gsElimWInverse" );
638 
641  "default",
642  "CLJP",
643  "RugeStueben",
644  "Falgout",
645  "PMIS",
646  "HMIS" );
647 
650  "default",
651  "modifiedClassical",
652  "direct",
653  "multipass",
654  "extendedI",
655  "standard",
656  "extended",
657  "directBAMG",
658  "modifiedExtended",
659  "modifiedExtendedI",
660  "modifiedExtendedE" );
661 
664  "default",
665  "extendedIStage2",
666  "standardStage2",
667  "extendedStage2",
668  "multipass",
669  "modifiedExtended",
670  "modifiedExtendedI",
671  "modifiedExtendedE",
672  "modifiedMultipass" );
673 
676  "constantModes",
677  "rigidBodyModes" );
678 
681  "D",
682  "DU",
683  "LD",
684  "LDU" );
685 
688  "none",
689  "diagonal",
690  "probing",
691  "user" );
692 
695  "none",
696  "frobenius",
697  "user" );
698 
701  "msrsb" );
702 
705  "graph",
706  "cartesian",
707  "semistructured" );
708 
711  "metis",
712  "scotch" );
713 
716  "kway",
717  "recursive" );
718 
721  "layers",
722  "matching" );
723 
724 } /* namespace geos */
725 
726 #endif /*GEOS_LINEARALGEBRA_UTILITIES_LINEARSOLVERPARAMETERS_HPP_ */
stdVector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:361
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:87
LvArray::CRSMatrix< T, COL_INDEX, INDEX_TYPE, LvArray::ChaiBuffer > CRSMatrix
Alias for CRS Matrix class.
Definition: DataTypes.hpp:305
double real64
64-bit floating point type.
Definition: DataTypes.hpp:98
int integer
Signed integer type.
Definition: DataTypes.hpp:81
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:175
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "richardson", "preconditioner")
Declare strings associated with enumeration values.
Algebraic multigrid parameters.
integer numSweeps
Number of smoother sweeps.
InterpType
AMG interpolation type (HYPRE only)
@ modifiedExtendedI
Modularized extended+i (GPU support)
@ modifiedExtendedE
Modularized extended+e (GPU support)
@ modifiedExtended
Modularized extended classical (GPU support)
@ directBAMG
Direct with separation of weights (GPU support)
@ extended
Extended classical (GPU support)
@ chebyshev
Chebyshev polynomial smoothing.
@ iluk
Incomplete LU with k-level of fill.
@ ick
Incomplete Cholesky with k-level of fill.
@ l1sgs
l1-Symmetric Gauss-Seidel smoothing
@ ilut
Incomplete LU with thresholding.
@ sgs
Symmetric Gauss-Seidel smoothing.
@ bgs
Gauss-Seidel smoothing (backward sweep)
@ ict
Incomplete Cholesky with thresholding.
@ fgs
Gauss-Seidel smoothing (forward sweep)
integer aggressiveInterpMaxNonZeros
Aggressive Interpolation - Max. nonzeros/row.
CoarseningType coarseningType
Coarsening algorithm (GPUs)
integer separateComponents
Apply a separate component filter before AMG construction.
integer maxLevels
Maximum number of coarsening levels.
integer interpolationMaxNonZeros
Interpolation - Max. nonzeros/row.
InterpType interpolationType
Interpolation algorithm.
CoarseType coarseType
Coarse-level solver/smoother.
CoarseningType
AMG coarsening types (HYPRE only)
@ CLJP
A parallel coarsening algorithm using independent sets.
@ Falgout
Ruge-Stueben followed by CLJP.
@ RugeStueben
Classical Ruge-Stueben on each processor, followed by a third pass.
@ PMIS
Parallel coarsening as CLJP but with lower complexities (GPU support)
integer numFunctions
Number of amg functions.
NullSpaceType nullSpaceType
Null space type [constantModes,rigidBodyModes].
SmootherType smootherType
Smoother type (GPUs)
integer maxCoarseSize
Threshold for coarse grid size.
PreOrPost preOrPostSmoothing
Pre and/or post smoothing.
AggInterpType
AMG interpolation type for aggressive coarsening levels (HYPRE only)
@ modifiedExtendedI
Modularized Extended+i (GPU support)
@ extendedIStage2
Extended+i 2-stage (GPU support)
@ modifiedExtendedE
Modularized Extended+e (GPU support)
@ modifiedExtended
Modularized Extended (GPU support)
@ modifiedMultipass
Modularized Multipass (GPU support)
@ extendedStage2
Extended 2-stage (GPU support)
PreOrPost
AMG pre/post smoothing option.
integer numCycles
Number of multigrid cycles.
integer aggressiveNumPaths
Number of paths agg. coarsening.
integer aggressiveNumLevels
Number of levels for aggressive coarsening.
@ jacobi
Jacobi (GPU support in hypre)
@ chebyshev
Chebyshev polynomial (GPU support in hypre)
@ l1jacobi
l1-Jacobi (GPU support in hypre)
@ gsElimWInverse
Direct inverse with Gaussian Elimination.
@ direct
Direct solver as preconditioner.
@ gsElimWPivoting
Gaussian Elimination with pivoting direct solver.
@ bgs
Gauss-Seidel smoothing (backward sweep)
AggInterpType aggressiveInterpType
Interp. type for agg. coarsening.
Block preconditioner parameters.
array1d< integer > order
Order of application of sub-problem solvers.
array1d< LinearSolverParameters const * > subParams
Pointers to parameters for sub-problems.
Shape shape
Block preconditioner shape.
Scaling scaling
Type of system scaling to use.
SchurType
Type of Schur complement approximation used.
@ FirstBlockDiagonal
Approximate first block with its diagonal.
@ None
No Schur complement - just block-GS/block-Jacobi preconditioner.
@ FirstBlockUserDefined
User defined preconditioner for the first block.
@ RowsumDiagonalProbing
Rowsum-preserving diagonal approximation constructed with probing.
void resize(integer const numBlocks)
Set the number of blocks and resize arrays accordingly.
Scaling
Type of block row scaling to apply.
@ FrobeniusNorm
Equilibrate Frobenius norm of the diagonal blocks.
Shape
Shape of the block preconditioner.
SchurType schurType
Schur complement type.
Chebyshev iteration/smoothing parameters.
integer eigNumIter
Number of eigenvalue estimation CG iterations.
Domain decomposition parameters.
Direct solver parameters: used for SuperLU_Dist interface through hypre and PETSc.
integer equilibrate
Whether to scale the rows and columns of the matrix.
integer iterativeRefine
Whether to perform iterative refinement.
integer replaceTinyPivot
Whether to replace tiny pivots by sqrt(epsilon)*norm(A)
integer parallel
Whether to use a parallel solver (instead of a serial one)
integer reuseFactorization
Whether to reuse the LU factorization or not.
@ MMD_AtplusA
multiple minimum degree on At+A
@ colAMD
approximate minimum degree on columns
@ MMD_AtA
multiple minimum degree on At*A (heavy)
integer checkResidual
Whether to check the linear system solution residual.
Incomplete factorization parameters.
real64 adaptiveGamma
Gamma parameter for adaptive method.
integer maxRestart
Max number of vectors in Krylov basis before restarting (GPUs)
real64 weakestTol
Weakest allowed tolerance when using adaptive method.
real64 relTolerance
Relative convergence tolerance for iterative solvers.
integer maxIterations
Max iterations before declaring convergence failure.
real64 adaptiveExponent
Exponent parameter for adaptive method.
real64 strongestTol
Strongest allowed tolerance when using adaptive method.
integer useAdaptiveTol
Use Eisenstat-Walker adaptive tolerance.
Multigrid reduction parameters.
integer separateComponents
Apply a separate displacement component (SDC) filter before AMG construction.
StrategyType strategy
Predefined MGR solution strategy (solver specific)
@ hybridSinglePhasePoromechanics
single phase poromechanics with hybrid finite volume single phase flow
@ compositionalMultiphaseFVM
finite volume compositional multiphase flow
@ singlePhasePoromechanics
single phase poromechanics with finite volume single phase flow
@ compositionalMultiphaseHybridFVM
hybrid finite volume compositional multiphase flow
@ solidMechanicsEmbeddedFractures
Embedded fractures mechanics.
@ singlePhaseReservoirHybridFVM
hybrid finite volume single-phase flow with wells
@ multiphasePoromechanicsReservoirFVM
multiphase poromechanics with finite volume compositional multiphase flow with wells
@ thermalCompositionalMultiphaseReservoirFVM
finite volume thermal compositional multiphase flow
@ thermalCompositionalMultiphaseFVM
finite volume thermal compositional multiphase flow
@ singlePhaseReservoirFVM
finite volume single-phase flow with wells
@ augmentedLagrangianContactMechanics
Augmented Lagrangian contact mechanics.
@ multiphasePoromechanics
multiphase poromechanics with finite volume compositional multiphase flow
@ singlePhasePoromechanicsReservoirFVM
single phase poromechanics with finite volume single phase flow with wells
@ singlePhaseHybridFVM
hybrid finite volume single-phase flow
@ thermalMultiphasePoromechanics
thermal multiphase poromechanics with finite volume compositional multiphase flow
@ singlePhasePoromechanicsConformingFractures
single phase poromechanics with conforming fractures
@ immiscibleMultiphaseFVM
finite volume immiscible multiphase flow
@ singlePhasePoromechanicsEmbeddedFractures
single phase poromechanics with FV embedded fractures
@ reactiveCompositionalMultiphaseOBL
finite volume reactive compositional flow with OBL
@ thermalSinglePhasePoromechanics
thermal single phase poromechanics with finite volume single phase flow
@ lagrangianContactMechanicsBubbleStab
Lagrangian contact mechanics with bubble stabilization.
@ compositionalMultiphaseReservoirHybridFVM
hybrid finite volume compositional multiphase flow with wells
@ lagrangianContactMechanics
Lagrangian contact mechanics.
@ compositionalMultiphaseReservoirFVM
finite volume compositional multiphase flow with wells
@ invalid
default value, to ensure solver sets something
integer ufactor
METIS UFACTOR option (allowed partition imbalance)
integer preserveRegions
Attempt to keep cells from the same region in one aggregate.
CRSMatrix< real64, localIndex > const * localMatrix
Local matrix to use for weights.
struct geos::LinearSolverParameters::Multiscale::Coarsening::Graph::Metis metis
METIS parameters.
integer matrixWeights
If >0, specifies matrix weight multiplier.
integer minCommonNodes
Min number of common nodes when constructing a cell connectivity graph.
struct geos::LinearSolverParameters::Multiscale::Coarsening::Graph::Scotch scotch
Scotch parameters.
struct geos::LinearSolverParameters::Multiscale::Coarsening::Graph graph
Graph coarsening parameters.
struct geos::LinearSolverParameters::Multiscale::Coarsening::Structured structured
Structured coarsening parameters.
PartitionType partitionType
Partitioning/coarsening method.
array1d< real64 > ratio
Coarsening ratio (cell-based), total or per-dimension.
globalIndex maxCoarseDof
Limit of coarsening globally (trims the grid hierarchy)
integer useBlockSmoother
Whether to use block smoother.
integer updateFrequency
Coarse operator update frequency (num smoothing iterations)
real64 tolerance
Smoothing iteration convergence tolerance.
integer maxIter
Max number of smoothing iterations.
integer numLayers
Number of layers to grow (for support = layers)
real64 relaxation
Relaxation parameter for Jacobi smoothing.
integer checkFrequency
Convergence check frequency (num smoothing iterations)
SupportType support
algorithm used to construct support regions
AMG::PreOrPost preOrPost
Pre and/or post smoothing [pre,post,both].
Multiscale preconditioner parameters.
struct geos::LinearSolverParameters::Multiscale::Coupled coupled
Multiscale coupled parameters.
integer galerkin
Whether to use Galerkin definition R = P^T (otherwise R = P_0^T)
integer maxLevels
Limit on total number of grid levels.
struct geos::LinearSolverParameters::Multiscale::MsRSB msrsb
MsRSB parameters.
struct geos::LinearSolverParameters::Multiscale::Coarsening coarsening
Multiscale coarsening parameters.
integer debugLevel
Flag for enabling addition debugging output.
string fieldName
DofManager field name, populated by the physics solver.
string label
User-displayed label of the scheme.
PreconditionerType coarseType
Coarse solver type.
real64 droptol
Dropping tolerance for coarse matrix values (relative to row max)
@ msrsb
Restricted Smoothing Basis Multiscale.
integer separateComponents
Apply a separate component filter before multiscale.
struct geos::LinearSolverParameters::Multiscale::Smoother smoother
Multiscale smoother parameters.
string_array boundarySets
List of boundary node set names (needed for coarse node detection)
BasisType basisType
Type of basis functions.
Relaxation/stationary iteration parameters (Richardson, damped Jacobi, etc.)
real64 weight
Relaxation weight (omega) for stationary iterations.
integer useRowColScaling
Apply row and column scaling (not yet implemented)
Set of parameters for a linear solver or preconditioner.
struct geos::LinearSolverParameters::Direct direct
direct solver parameter struct
integer dofsPerNode
Dofs per node (or support location) for non-scalar problems.
struct geos::LinearSolverParameters::Chebyshev chebyshev
Chebyshev smoother parameters.
struct geos::LinearSolverParameters::IFact ifact
Incomplete factorization parameter struct.
struct geos::LinearSolverParameters::Relaxation relaxation
Relaxation method parameters.
bool isSymmetric
Whether input matrix is symmetric (may affect choice of scheme)
struct geos::LinearSolverParameters::Scaling scaling
Matrix-scaling parameter struct.
struct geos::LinearSolverParameters::DD dd
Domain decomposition parameter struct.
@ chebyshev
Chebyshev polynomial smoothing.
@ mgr
Multigrid reduction (Hypre only)
@ iluk
Incomplete LU with k-level of fill.
@ ick
Incomplete Cholesky with k-level of fill.
@ direct
Direct solver as preconditioner.
@ l1sgs
l1-Symmetric Gauss-Seidel smoothing
@ ilut
Incomplete LU with thresholding.
@ sgs
Symmetric Gauss-Seidel smoothing.
@ bgs
Gauss-Seidel smoothing (backward sweep)
@ ict
Incomplete Cholesky with thresholding.
@ fgs
Gauss-Seidel smoothing (forward sweep)
integer stopIfError
Whether to stop the simulation if the linear solver reports an error.
struct geos::LinearSolverParameters::Multiscale multiscale
Multiscale preconditioner parameters.
PreconditionerType preconditionerType
Preconditioner type.
struct geos::LinearSolverParameters::AMG amg
Algebraic Multigrid (AMG) parameters.
integer logLevel
Output level [0=none, 1=basic, 2=everything].
struct geos::LinearSolverParameters::Block block
Block preconditioner parameters.
struct geos::LinearSolverParameters::MGR mgr
Multigrid reduction (MGR) parameters.
struct geos::LinearSolverParameters::Krylov krylov
Krylov-method parameter struct.