20 #ifndef GEOS_LINEARALGEBRA_UTILITIES_LINEARSOLVERPARAMETERS_HPP_
21 #define GEOS_LINEARALGEBRA_UTILITIES_LINEARSOLVERPARAMETERS_HPP_
124 #if GEOS_USE_HYPRE_DEVICE == GEOS_USE_HYPRE_CUDA || GEOS_USE_HYPRE_DEVICE == GEOS_USE_HYPRE_HIP
259 #if GEOS_USE_HYPRE_DEVICE == GEOS_USE_HYPRE_CUDA || GEOS_USE_HYPRE_DEVICE == GEOS_USE_HYPRE_HIP
390 order.resize( numBlocks );
391 for(
integer i = 0; i < numBlocks; ++i )
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",
594 "lagrangianContactMechanics",
595 "augmentedLagrangianContactMechanics",
596 "lagrangianContactMechanicsBubbleStab",
597 "solidMechanicsEmbeddedFractures" );
660 "modifiedExtendedE" );
672 "modifiedMultipass" );
stdVector< string > string_array
A 1-dimensional array of geos::string types.
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
LvArray::CRSMatrix< T, COL_INDEX, INDEX_TYPE, LvArray::ChaiBuffer > CRSMatrix
Alias for CRS Matrix class.
double real64
64-bit floating point type.
int integer
Signed integer type.
Array< T, 1 > array1d
Alias for 1D array.
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)
@ multipass
Multipass (GPU support)
@ default_
Use LAI's default option.
@ extendedI
Extended+i (GPU support)
@ modifiedExtendedE
Modularized extended+e (GPU support)
@ modifiedExtended
Modularized extended classical (GPU support)
@ direct
Direct (GPU support)
@ directBAMG
Direct with separation of weights (GPU support)
@ modifiedClassical
Modified classical.
@ extended
Extended classical (GPU support)
SmootherType
AMG smoother type.
@ jacobi
Jacobi smoothing.
@ default_
Use LAI's default option.
@ chebyshev
Chebyshev polynomial smoothing.
@ iluk
Incomplete LU with k-level of fill.
@ l1jacobi
l1-Jacobi smoothing
@ 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.
@ default_
Use LAI's default option.
@ 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)
@ HMIS
Hybrid PMIS coarsening.
integer numFunctions
Number of amg functions.
NullSpaceType
Null space type.
@ constantModes
Constant modes.
@ rigidBodyModes
Rigid body modes.
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)
@ multipass
Multipass (GPU support)
@ default_
Use LAI's default option.
@ 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)
@ standardStage2
Standard 2-stage.
CycleType cycleType
AMG cycle type.
PreOrPost
AMG pre/post smoothing option.
@ post
post-smoothing only
@ both
pre- and post-smoothing
integer numCycles
Number of multigrid cycles.
integer aggressiveNumPaths
Number of paths agg. coarsening.
integer aggressiveNumLevels
Number of levels for aggressive coarsening.
CoarseType
AMG coarse solver type.
@ jacobi
Jacobi (GPU support in hypre)
@ default_
Use LAI's default option.
@ 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.
@ l1sgs
l1-Symmetric Gauss-Seidel
@ sgs
Symmetric Gauss-Seidel.
@ gsElimWPivoting
Gaussian Elimination with pivoting direct solver.
@ bgs
Gauss-Seidel smoothing (backward sweep)
@ fgs
Gauss-Seidel (forward sweep)
real64 relaxWeight
Relaxation weight.
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.
@ UserProvided
User-provided scaling.
@ FrobeniusNorm
Equilibrate Frobenius norm of the diagonal blocks.
Shape
Shape of the block preconditioner.
@ LowerTriangular
(LD)^{-1}
@ LowerUpperTriangular
(LDU)^{-1}
@ UpperTriangular
(DU)^{-1}
SchurType schurType
Schur complement type.
Chebyshev iteration/smoothing parameters.
integer order
Chebyshev order.
integer eigNumIter
Number of eigenvalue estimation CG iterations.
Domain decomposition parameters.
integer overlap
Ghost overlap.
Direct solver parameters: used for SuperLU_Dist interface through hypre and PETSc.
integer equilibrate
Whether to scale the rows and columns of the matrix.
RowPerm rowPerm
Rows permutation.
ColPerm colPerm
Columns permutation.
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.
ColPerm
How to permute the columns.
@ 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.
RowPerm
How to permute the rows.
@ mc64
using HSL routine MC64
Incomplete factorization parameters.
real64 threshold
Dropping threshold.
Krylov-method 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)
StrategyType
MGR available strategies.
@ 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
@ hydrofracture
hydrofracture
@ 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)
Method
METIS partitioning algorithm.
Method method
METIS partitioning algorithm.
integer seed
Random number generator seed.
Graph coarsening parameters.
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.
Method method
Graph partitioning method (library)
Method
Graph partitioning method (library) to use.
integer minCommonNodes
Min number of common nodes when constructing a cell connectivity graph.
struct geos::LinearSolverParameters::Multiscale::Coarsening::Graph::Scotch scotch
Scotch parameters.
Structured coarsening parameters.
integer semicoarsening
Use semicoarsenining in Z-direction.
Multiscale coarsening parameters.
PartitionType
Type of cell partitioning.
@ semistructured
For 2.5D (extruded) meshes.
@ cartesian
Cartesian (only with internal mesh)
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)
Multiscale coupled parameters.
integer useBlockSmoother
Whether to use block smoother.
integer updateFrequency
Coarse operator update frequency (num smoothing iterations)
SupportType
Type of support regions.
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
Multiscale smoother parameters.
integer numSweeps
Number of smoother sweeps.
PreconditionerType type
Smoother type.
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)
BasisType
Type of basis functions.
@ 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.
Matrix-scaling parameters.
integer useRowScaling
Apply row scaling.
integer useRowColScaling
Apply row and column scaling (not yet implemented)
Set of parameters for a linear solver or preconditioner.
SolverType
Linear solver type.
@ richardson
Richardson iteration.
@ preconditioner
Preconditioner only.
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.
PreconditionerType
Preconditioner type.
@ amg
Algebraic Multigrid.
@ block
Block preconditioner.
@ jacobi
Jacobi smoothing.
@ chebyshev
Chebyshev polynomial smoothing.
@ mgr
Multigrid reduction (Hypre only)
@ iluk
Incomplete LU with k-level of fill.
@ l1jacobi
l1-Jacobi smoothing
@ 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.
@ multiscale
Multiscale preconditioner.
@ 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.
SolverType solverType
Solver type.
struct geos::LinearSolverParameters::MGR mgr
Multigrid reduction (MGR) parameters.
struct geos::LinearSolverParameters::Krylov krylov
Krylov-method parameter struct.