8.3
general documentation
cs_param_types.h File Reference
#include "cs_defs.h"
+ Include dependency graph for cs_param_types.h:

Go to the source code of this file.

Data Structures

struct  cs_param_convergence_t
 Set of parameters to check the convergence (or the divergence) of an iterative process (tolerances or max. number of iterations) More...
 

Macros

#define CS_DOF_VTX_SCAL   0
 
#define CS_DOF_VTX_VECT   1
 
#define CS_DOF_FACE_SCAL   2
 
#define CS_DOF_FACE_VECT   3
 
#define CS_DOF_FACE_SCAP1   3
 
#define CS_DOF_FACE_SCAP2   4
 
#define CS_DOF_FACE_VECP0   3
 
#define CS_DOF_FACE_VECP1   5
 
#define CS_DOF_FACE_VECP2   6
 
#define CS_DOF_EDGE_SCAL   7
 
#define CS_N_DOF_CASES   8
 
#define CS_ISOTROPIC_DIFFUSION   (1 << 0)
 
#define CS_ORTHOTROPIC_DIFFUSION   (1 << 1)
 
#define CS_ANISOTROPIC_LEFT_DIFFUSION   (1 << 2)
 
#define CS_ANISOTROPIC_RIGHT_DIFFUSION   (1 << 3)
 
#define CS_ANISOTROPIC_DIFFUSION   ((1 << 2) + (1 << 3))
 

Typedefs

typedef void() cs_analytic_func_t(cs_real_t time, cs_lnum_t n_elts, const cs_lnum_t *elt_ids, const cs_real_t *coords, bool dense_output, void *input, cs_real_t *retval)
 Generic function pointer for an evaluation relying on an analytic function. More...
 
typedef void() cs_dof_func_t(cs_lnum_t n_elts, const cs_lnum_t *elt_ids, bool dense_output, void *input, cs_real_t *retval)
 Generic function pointer for computing a quantity at predefined locations such as degrees of freedom (DoF): cells, faces, edges or or vertices. More...
 
typedef void() cs_time_func_t(double time, void *input, cs_real_t *retval)
 Function which defines the evolution of a quantity according to the current time and any structure given as a parameter. More...
 

Enumerations

enum  cs_param_space_scheme_t {
  CS_SPACE_SCHEME_LEGACY , CS_SPACE_SCHEME_CDOVB , CS_SPACE_SCHEME_CDOVCB , CS_SPACE_SCHEME_CDOEB ,
  CS_SPACE_SCHEME_CDOFB , CS_SPACE_SCHEME_CDOCB , CS_SPACE_SCHEME_HHO_P0 , CS_SPACE_SCHEME_HHO_P1 ,
  CS_SPACE_SCHEME_HHO_P2 , CS_SPACE_SCHEME_MACFB , CS_SPACE_N_SCHEMES
}
 Type of numerical scheme for the discretization in space. More...
 
enum  cs_param_dof_reduction_t { CS_PARAM_REDUCTION_DERHAM , CS_PARAM_REDUCTION_AVERAGE , CS_PARAM_N_REDUCTIONS }
 

Functions

bool cs_param_space_scheme_is_face_based (cs_param_space_scheme_t scheme)
 Return true if the space scheme has degrees of freedom on faces, otherwise false. More...
 
const char * cs_param_get_space_scheme_name (cs_param_space_scheme_t scheme)
 Get the name of the space discretization scheme. More...
 
const char * cs_param_get_time_scheme_name (cs_param_time_scheme_t scheme)
 Get the name of the time discretization scheme. More...
 
const char * cs_param_get_advection_form_name (cs_param_advection_form_t adv_form)
 Get the label associated to the advection formulation. More...
 
const char * cs_param_get_advection_scheme_name (cs_param_advection_scheme_t scheme)
 Get the label of the advection scheme. More...
 
const char * cs_param_get_advection_strategy_name (cs_param_advection_strategy_t adv_stra)
 Get the label associated to the advection strategy. More...
 
const char * cs_param_get_advection_extrapol_name (cs_param_advection_extrapol_t extrapol)
 Get the label associated to the extrapolation used for the advection field. More...
 
const char * cs_param_get_bc_name (cs_param_bc_type_t bc)
 Get the name of the type of boundary condition. More...
 
const char * cs_param_get_bc_enforcement_name (cs_param_bc_enforce_t type)
 Get the name of the type of enforcement of the boundary condition. More...
 
const char * cs_param_get_nl_algo_name (cs_param_nl_algo_t algo)
 Get the name of the non-linear algorithm. More...
 
const char * cs_param_get_nl_algo_label (cs_param_nl_algo_t algo)
 Get the label (short name) of the non-linear algorithm. More...
 
const char * cs_param_get_dotprod_type_name (cs_param_dotprod_type_t dp_type)
 Get the name of the type of dot product to apply. More...
 
const char * cs_param_get_solver_name (cs_param_solver_type_t solver)
 Get the name of the solver. More...
 
const char * cs_param_get_precond_name (cs_param_precond_type_t precond)
 Get the name of the preconditioner. More...
 
const char * cs_param_get_precond_block_name (cs_param_precond_block_t type)
 Get the name of the type of block preconditioning. More...
 

Variables

const char cs_sep_h1 [80]
 
const char cs_sep_h2 [80]
 
const char cs_sepline [80]
 
const char cs_med_sepline [50]
 

Settings for the boundary conditions

#define CS_PARAM_BC_HMG_DIRICHLET   CS_BC_HMG_DIRICHLET
 
#define CS_PARAM_BC_DIRICHLET   CS_BC_DIRICHLET
 
#define CS_PARAM_BC_HMG_NEUMANN   CS_BC_SYMMETRY
 
#define CS_PARAM_BC_NEUMANN   CS_BC_NEUMANN
 
#define CS_PARAM_BC_NEUMANN_FULL   CS_BC_NEUMANN_FULL
 
#define CS_PARAM_BC_ROBIN   CS_BC_ROBIN
 
#define CS_PARAM_BC_SLIDING   CS_BC_SYMMETRY
 
#define CS_PARAM_BC_CIRCULATION   CS_BC_CIRCULATION
 
#define CS_PARAM_BC_WALL_PRECRIBED   CS_BC_WALL_MODELLED
 
enum  cs_param_bc_type_t {
  CS_BC_HMG_DIRICHLET = 0 , CS_BC_DIRICHLET = 1 , CS_BC_RADIATIVE_OUTLET = 2 , CS_BC_NEUMANN = 3 ,
  CS_BC_SYMMETRY = 4 , CS_BC_WALL_MODELLED = 5 , CS_BC_ROUGH_WALL_MODELLED = 6 , CS_BC_NEUMANN_FULL = 7 ,
  CS_BC_ROBIN = 9 , CS_BC_CIRCULATION = 11 , CS_BC_IMPOSED_TOT_FLUX = 13 , CS_BC_GENERALIZED_SYM = 14 ,
  CS_BC_IMPOSED_EXCHANGE_COEF = 15 , CS_PARAM_N_BC_TYPES
}
 
enum  cs_param_bc_enforce_t {
  CS_PARAM_BC_ENFORCE_ALGEBRAIC , CS_PARAM_BC_ENFORCE_PENALIZED , CS_PARAM_BC_ENFORCE_WEAK_NITSCHE , CS_PARAM_BC_ENFORCE_WEAK_SYM ,
  CS_PARAM_N_BC_ENFORCEMENTS
}
 

Numerical settings for time discretization

enum  cs_param_time_scheme_t {
  CS_TIME_SCHEME_STEADY , CS_TIME_SCHEME_EULER_IMPLICIT , CS_TIME_SCHEME_EULER_EXPLICIT , CS_TIME_SCHEME_CRANKNICO ,
  CS_TIME_SCHEME_THETA , CS_TIME_SCHEME_BDF2 , CS_TIME_N_SCHEMES
}
 

Settings for the advection term

enum  cs_param_advection_form_t { CS_PARAM_ADVECTION_FORM_CONSERV , CS_PARAM_ADVECTION_FORM_NONCONS , CS_PARAM_ADVECTION_FORM_SKEWSYM , CS_PARAM_N_ADVECTION_FORMULATIONS }
 
enum  cs_param_advection_scheme_t {
  CS_PARAM_ADVECTION_SCHEME_CENTERED , CS_PARAM_ADVECTION_SCHEME_CIP , CS_PARAM_ADVECTION_SCHEME_CIP_CW , CS_PARAM_ADVECTION_SCHEME_HYBRID_CENTERED_UPWIND ,
  CS_PARAM_ADVECTION_SCHEME_SAMARSKII , CS_PARAM_ADVECTION_SCHEME_SG , CS_PARAM_ADVECTION_SCHEME_UPWIND , CS_PARAM_N_ADVECTION_SCHEMES
}
 
enum  cs_param_advection_strategy_t { CS_PARAM_ADVECTION_IMPLICIT_FULL , CS_PARAM_ADVECTION_IMPLICIT_LINEARIZED , CS_PARAM_ADVECTION_EXPLICIT , CS_PARAM_N_ADVECTION_STRATEGIES }
 Choice of how to handle the advection term in an equation. More...
 
enum  cs_param_advection_extrapol_t { CS_PARAM_ADVECTION_EXTRAPOL_NONE , CS_PARAM_ADVECTION_EXTRAPOL_TAYLOR_2 , CS_PARAM_ADVECTION_EXTRAPOL_ADAMS_BASHFORTH_2 , CS_PARAM_N_ADVECTION_EXTRAPOLATIONS }
 Choice of how to extrapolate the advection field in the advection term. More...
 

Settings for non-linear algorithms

enum  cs_param_nl_algo_t {
  CS_PARAM_NL_ALGO_NONE , CS_PARAM_NL_ALGO_PICARD , CS_PARAM_NL_ALGO_MODIFIED_PICARD , CS_PARAM_NL_ALGO_ANDERSON ,
  CS_PARAM_N_NL_ALGOS
}
 Class of non-linear iterative algorithm. More...
 

Settings for the linear solvers or SLES (Sparse Linear Equation Solver)

enum  cs_param_solver_class_t {
  CS_PARAM_SOLVER_CLASS_CS , CS_PARAM_SOLVER_CLASS_HYPRE , CS_PARAM_SOLVER_CLASS_MUMPS , CS_PARAM_SOLVER_CLASS_PETSC ,
  CS_PARAM_N_SOLVER_CLASSES
}
 Class of iterative solvers to consider for solver the linear system. More...
 
enum  cs_param_precond_block_t {
  CS_PARAM_PRECOND_BLOCK_NONE , CS_PARAM_PRECOND_BLOCK_DIAG , CS_PARAM_PRECOND_BLOCK_LOWER_TRIANGULAR , CS_PARAM_PRECOND_BLOCK_SYM_GAUSS_SEIDEL ,
  CS_PARAM_PRECOND_BLOCK_UPPER_TRIANGULAR , CS_PARAM_N_PCD_BLOCK_TYPES
}
 
enum  cs_param_precond_type_t {
  CS_PARAM_PRECOND_NONE , CS_PARAM_PRECOND_AMG , CS_PARAM_PRECOND_BJACOB_ILU0 , CS_PARAM_PRECOND_BJACOB_SGS ,
  CS_PARAM_PRECOND_DIAG , CS_PARAM_PRECOND_GKB_CG , CS_PARAM_PRECOND_GKB_GMRES , CS_PARAM_PRECOND_LU ,
  CS_PARAM_PRECOND_ILU0 , CS_PARAM_PRECOND_ICC0 , CS_PARAM_PRECOND_MUMPS , CS_PARAM_PRECOND_HPDDM ,
  CS_PARAM_PRECOND_POLY1 , CS_PARAM_PRECOND_POLY2 , CS_PARAM_PRECOND_SSOR , CS_PARAM_N_PRECOND_TYPES
}
 
enum  cs_param_solver_type_t {
  CS_PARAM_SOLVER_NONE , CS_PARAM_SOLVER_AMG , CS_PARAM_SOLVER_BICGS , CS_PARAM_SOLVER_BICGS2 ,
  CS_PARAM_SOLVER_CG , CS_PARAM_SOLVER_CR3 , CS_PARAM_SOLVER_FCG , CS_PARAM_SOLVER_FGMRES ,
  CS_PARAM_SOLVER_GAUSS_SEIDEL , CS_PARAM_SOLVER_GCR , CS_PARAM_SOLVER_GMRES , CS_PARAM_SOLVER_JACOBI ,
  CS_PARAM_SOLVER_MUMPS , CS_PARAM_SOLVER_SYM_GAUSS_SEIDEL , CS_PARAM_SOLVER_USER_DEFINED , CS_PARAM_N_SOLVER_TYPES
}
 
enum  cs_param_resnorm_type_t {
  CS_PARAM_RESNORM_NONE , CS_PARAM_RESNORM_NORM2_RHS , CS_PARAM_RESNORM_WEIGHTED_RHS , CS_PARAM_RESNORM_FILTERED_RHS ,
  CS_PARAM_N_RESNORM_TYPES
}
 
enum  cs_param_dotprod_type_t { CS_PARAM_DOTPROD_EUCLIDEAN , CS_PARAM_DOTPROD_CDO , CS_PARAM_N_DOTPROD_TYPES }
 

Macro Definition Documentation

◆ CS_PARAM_BC_CIRCULATION

#define CS_PARAM_BC_CIRCULATION   CS_BC_CIRCULATION

Homogeneous Dirichlet boundary conditions. The value of a variable is set to zero.

◆ CS_PARAM_BC_DIRICHLET

#define CS_PARAM_BC_DIRICHLET   CS_BC_DIRICHLET

Homogeneous Dirichlet boundary conditions. The value of a variable is set to zero.

◆ CS_PARAM_BC_HMG_DIRICHLET

#define CS_PARAM_BC_HMG_DIRICHLET   CS_BC_HMG_DIRICHLET

Homogeneous Dirichlet boundary conditions. The value of a variable is set to zero.

◆ CS_PARAM_BC_HMG_NEUMANN

#define CS_PARAM_BC_HMG_NEUMANN   CS_BC_SYMMETRY

Homogeneous Dirichlet boundary conditions. The value of a variable is set to zero.

◆ CS_PARAM_BC_NEUMANN

#define CS_PARAM_BC_NEUMANN   CS_BC_NEUMANN

Homogeneous Dirichlet boundary conditions. The value of a variable is set to zero.

◆ CS_PARAM_BC_NEUMANN_FULL

#define CS_PARAM_BC_NEUMANN_FULL   CS_BC_NEUMANN_FULL

Homogeneous Dirichlet boundary conditions. The value of a variable is set to zero.

◆ CS_PARAM_BC_ROBIN

#define CS_PARAM_BC_ROBIN   CS_BC_ROBIN

Homogeneous Dirichlet boundary conditions. The value of a variable is set to zero.

◆ CS_PARAM_BC_SLIDING

#define CS_PARAM_BC_SLIDING   CS_BC_SYMMETRY

Homogeneous Dirichlet boundary conditions. The value of a variable is set to zero.

◆ CS_PARAM_BC_WALL_PRECRIBED

#define CS_PARAM_BC_WALL_PRECRIBED   CS_BC_WALL_MODELLED

Homogeneous Dirichlet boundary conditions. The value of a variable is set to zero.

Typedef Documentation

◆ cs_analytic_func_t

typedef void() cs_analytic_func_t(cs_real_t time, cs_lnum_t n_elts, const cs_lnum_t *elt_ids, const cs_real_t *coords, bool dense_output, void *input, cs_real_t *retval)

Generic function pointer for an evaluation relying on an analytic function.

For the calling function, elt_ids is optional. If non-null, the coords array should be accessed with an indirection. The same indirection can be applied to fill retval if dense_output is set to false.

Parameters
[in]timewhen ?
[in]n_eltsnumber of elements to consider
[in]elt_idslist of elements ids (in coords and retval)
[in]coordswhere ?
[in]dense_outputperform an indirection in retval or not
[in]inputnull or pointer to a structure cast on-the-fly
[in,out]retvalresulting value(s). Must be allocated.

◆ cs_dof_func_t

typedef void() cs_dof_func_t(cs_lnum_t n_elts, const cs_lnum_t *elt_ids, bool dense_output, void *input, cs_real_t *retval)

Generic function pointer for computing a quantity at predefined locations such as degrees of freedom (DoF): cells, faces, edges or or vertices.

For the calling function, elt_ids is optional. If non-null, array(s) should be accessed with an indirection. The same indirection can be applied to fill retval if dense_output is set to false.

Parameters
[in]n_eltsnumber of elements to consider
[in]elt_idslist of elements ids
[in]dense_outputperform an indirection in retval or not
[in]inputnull or pointer to a structure cast on-the-fly
[in,out]retvalresulting value(s). Must be allocated.

◆ cs_time_func_t

typedef void() cs_time_func_t(double time, void *input, cs_real_t *retval)

Function which defines the evolution of a quantity according to the current time and any structure given as a parameter.

Parameters
[in]timevalue of the time at the end of the last iteration
[in]inputnull or pointer to a structure cast on-the-fly
[in]retvalresult of the evaluation

Enumeration Type Documentation

◆ cs_param_advection_extrapol_t

Choice of how to extrapolate the advection field in the advection term.

Enumerator
CS_PARAM_ADVECTION_EXTRAPOL_NONE 

The advection field is not extrapolated. The last known advection field is considered $ \phi^n $ if one computes $ u^(n+1) $ knowing $ u^n. $ In case of a a non-linearity, this is $ \phi^(n+1,k) $ if one computes $ u^(n+1,k+1) $ knowing $  u^(n+1,k) $ and $ u^n.$ The initial step is such that $ u^(n+1,0) = u^n $. This is a good choice with a first-order forward or backward time scheme.

CS_PARAM_ADVECTION_EXTRAPOL_TAYLOR_2 

The advection field is extrapolated with a 2nd order Taylor expansion yielding $ \phi^extrap = 2 \phi^n - \phi^(n-1) $. This corresponds to an estimation of $ \phi $ at n+1. Thus, this is a good choice when associated with a BDF2 time scheme.

CS_PARAM_ADVECTION_EXTRAPOL_ADAMS_BASHFORTH_2 

The advection field is extrapolated with a 2nd order Adams-Bashforth technique yielding $ \phi^extrap = 3/2 \phi^n - 1/2 \phi^(n-1) $. This corresponds to an estimation of $ \phi $ at n+1/2. Thus, this ) is a good choice when associated with a Crank-Nilcolson time scheme.

CS_PARAM_N_ADVECTION_EXTRAPOLATIONS 

◆ cs_param_advection_form_t

Type of formulation for the advection term

Enumerator
CS_PARAM_ADVECTION_FORM_CONSERV 

conservative formulation (i.e. divergence of a flux)

CS_PARAM_ADVECTION_FORM_NONCONS 

non-conservative formation (i.e advection field times the gradient)

CS_PARAM_ADVECTION_FORM_SKEWSYM 

skew-symmetric form

CS_PARAM_N_ADVECTION_FORMULATIONS 

◆ cs_param_advection_scheme_t

Type of numerical scheme for the advection term

Enumerator
CS_PARAM_ADVECTION_SCHEME_CENTERED 

centered discretization

CS_PARAM_ADVECTION_SCHEME_CIP 

Continuous Interior Penalty discretization. Only available for CS_SPACE_SCHEME_CDOVCB space discretization.

CS_PARAM_ADVECTION_SCHEME_CIP_CW 

Continuous Interior Penalty discretization. Only available for CS_SPACE_SCHEME_CDOVCB space discretization. Variant with a cellwise constant approximation of the advection field

CS_PARAM_ADVECTION_SCHEME_HYBRID_CENTERED_UPWIND 

Centered discretization with a portion between [0,1] of upwinding. The portion is specified thanks to CS_EQKEY_ADV_UPWIND_PORTION If the portion is equal to 0, then one recovers CS_PARAM_ADVECTION_SCHEME_CENTERED. If the portion is equal to 1, then one recovers CS_PARAM_ADVECTION_SCHEME_UPWIND

CS_PARAM_ADVECTION_SCHEME_SAMARSKII 

Weighting between an upwind and a centered discretization relying on the Peclet number. Weighting function = Samarskii

CS_PARAM_ADVECTION_SCHEME_SG 

Weighting between an upwind and a centered discretization relying on the Peclet number. Weighting function = Scharfetter-Gummel

CS_PARAM_ADVECTION_SCHEME_UPWIND 

Low order upwind discretization

CS_PARAM_N_ADVECTION_SCHEMES 

◆ cs_param_advection_strategy_t

Choice of how to handle the advection term in an equation.

Enumerator
CS_PARAM_ADVECTION_IMPLICIT_FULL 

The advection term is implicitly treated. The non-linearity stemming from this term has to be solved using a specific algorithm such as the Picard (fixed-point) technique or more elaborated techniques.

CS_PARAM_ADVECTION_IMPLICIT_LINEARIZED 

The advection term is implicitly treated. The non-linearity stemming from this term is simplified. Namely, one assumes a linearized advection. This is equivalent to a one-step Picard technique.

CS_PARAM_ADVECTION_EXPLICIT 

The advection term is treated explicitly. One keeps the non-linearity stemming from this term at the right hand-side. An extrapolation can be used for the advection field (cf. cs_param_advection_extrapol_t)

CS_PARAM_N_ADVECTION_STRATEGIES 

◆ cs_param_bc_enforce_t

Type of method for enforcing the boundary conditions. According to the type of numerical scheme, some enforcements are not available.

Enumerator
CS_PARAM_BC_ENFORCE_ALGEBRAIC 

Weak enforcement of the boundary conditions (i.e. one keeps the degrees of freedom in the algebraic system) with an algebraic manipulation of the linear system

CS_PARAM_BC_ENFORCE_PENALIZED 

Weak enforcement of the boundary conditions (i.e. one keeps the degrees of freedom in the algebraic system) with a penalization technique using a huge value.

CS_PARAM_BC_ENFORCE_WEAK_NITSCHE 

Weak enforcement of the boundary conditions (i.e. one keeps the degrees of freedom in the algebraic system) with a Nitsche-like penalization technique. This technique does not increase the conditioning number as much as CS_PARAM_BC_ENFORCE_PENALIZED but is more computationally intensive. The computation of the diffusion term should be activated with this choice.

CS_PARAM_BC_ENFORCE_WEAK_SYM 

Weak enforcement of the boundary conditions (i.e. one keeps the degrees of freedom in the algebraic system) with a Nitsche-like penalization technique. This variant enables to keep the symmetry of the algebraic contrary to CS_PARAM_BC_ENFORCE_WEAK_NITSCHE. The computation of the diffusion term should be activated with this choice.

CS_PARAM_N_BC_ENFORCEMENTS 

◆ cs_param_bc_type_t

Type of boundary condition to enforce.

Enumerator
CS_BC_HMG_DIRICHLET 

Homogeneous Dirichlet boundary conditions. The value of a variable is set to zero.

CS_BC_DIRICHLET 

Dirichlet boundary conditions. The value of the variable is set to the user requirements.

CS_BC_RADIATIVE_OUTLET 
CS_BC_NEUMANN 

Neumann conditions (along the face normal). The value of the flux of variable is set to the user requirements.

CS_BC_SYMMETRY 

Homogeneous Neumann conditions. The value of the flux of variable is set to zero.

Sliding conditions. Homogeneous Dirichlet for the normal component and homogeneous Neumann for the tangential components. Only available for vector-valued equations.

CS_BC_WALL_MODELLED 

Hybrid way to enforce the value at a wall. Use a wall function to prescribe the wall value knowing the value at the associated cell center for instance

CS_BC_ROUGH_WALL_MODELLED 
CS_BC_NEUMANN_FULL 

Neumann conditions with full definition relative to all directions, i.e. not only the normal direction. The value of the flux of variable is set to the user requirements.

CS_BC_ROBIN 

Robin conditions.

CS_BC_CIRCULATION 

Set the tangential part of a vector-valued field. This is the part lying on the boundary of a part of the computatinal domain. Nothing is prescribed for the normal part of the vector-valued field.

CS_BC_IMPOSED_TOT_FLUX 
CS_BC_GENERALIZED_SYM 
CS_BC_IMPOSED_EXCHANGE_COEF 
CS_PARAM_N_BC_TYPES 

◆ cs_param_dof_reduction_t

Enumerator
CS_PARAM_REDUCTION_DERHAM 

Evaluation at vertices for potentials, integral along a line for circulations, integral across the normal component of a face for fluxes and integral inside a cell for densities

CS_PARAM_REDUCTION_AVERAGE 

Degrees of freedom are defined as the average on the element

CS_PARAM_N_REDUCTIONS 

◆ cs_param_dotprod_type_t

Way to compute a dot product or its related norm

Enumerator
CS_PARAM_DOTPROD_EUCLIDEAN 

Usual Eucliean 2-norm (no weight)

CS_PARAM_DOTPROD_CDO 

Weighted 2-norm from CDO quantities

CS_PARAM_N_DOTPROD_TYPES 

◆ cs_param_nl_algo_t

Class of non-linear iterative algorithm.

Enumerator
CS_PARAM_NL_ALGO_NONE 

No specific algorithm. This mean a linear technique.

CS_PARAM_NL_ALGO_PICARD 

Picard (also called fixed point) algorithm

CS_PARAM_NL_ALGO_MODIFIED_PICARD 

modified Picard algorithm where the unsteady term is treated specifically in a manner close to what is done in the Newton algorithm.

CS_PARAM_NL_ALGO_ANDERSON 

Anderson acceleration

CS_PARAM_N_NL_ALGOS 

No specific algorithm. This mean a linear technique.

◆ cs_param_precond_block_t

Type of preconditioning by block

Enumerator
CS_PARAM_PRECOND_BLOCK_NONE 

No block preconditioner is requested (default)

CS_PARAM_PRECOND_BLOCK_DIAG 

Only the diagonal blocks are considered in the preconditioner

CS_PARAM_PRECOND_BLOCK_LOWER_TRIANGULAR 

The diagonal blocks and the lower blocks are considered in the preconditioner

CS_PARAM_PRECOND_BLOCK_SYM_GAUSS_SEIDEL 

A symmetric Gauss-Seidel block preconditioning is considered (cf. Y. Notay, "A new analysis of block preconditioners for saddle-point problems" (2014), SIAM J. Matrix. Anal. Appl.)

CS_PARAM_PRECOND_BLOCK_UPPER_TRIANGULAR 

The diagonal blocks and the upper blocks are considered in the preconditioner

CS_PARAM_N_PCD_BLOCK_TYPES 

◆ cs_param_precond_type_t

Type of preconditioner to use with the iterative solver. Some of the mentionned preconditioners are available only if the PETSc library is linked with code_saturne

Enumerator
CS_PARAM_PRECOND_NONE 

No preconditioner

CS_PARAM_PRECOND_AMG 

Algebraic multigrid preconditioner (additional options may be set using cs_param_amg_type_t)

CS_PARAM_PRECOND_BJACOB_ILU0 

Only with PETSc

Block Jacobi with an ILU zero fill-in in each block

CS_PARAM_PRECOND_BJACOB_SGS 

Only with PETSc

Block Jacobi with a symmetric Gauss-Seidel in each block (rely on Eisenstat's trick with PETSc)

CS_PARAM_PRECOND_DIAG 

Diagonal (also Jacobi) preconditioner. The cheapest one but not the most efficient one.

CS_PARAM_PRECOND_GKB_CG 

Only with PETSc

Golub-Kahan Bidiagonalization solver used as a preconditioner. Only useful if one has to solve a saddle-point system (such systems arise when solving Stokes or Navier-Stokes in a fully couple manner). Variant with CG as inner solver.

CS_PARAM_PRECOND_GKB_GMRES 

Only with PETSc

Golub-Kahan Bidiagonalization solver used as a preconditioner. Only useful if one has to solve a saddle-point system (such systems arise when solving Stokes or Navier-Stokes in a fully couple manner). Variant with GMRES as inner solver.

CS_PARAM_PRECOND_LU 

Only with PETSc

LU factorization (direct solver) with PETSc

CS_PARAM_PRECOND_ILU0 

Only with PETSc

Incomplete LU factorization (fill-in coefficient set to 0)

CS_PARAM_PRECOND_ICC0 

Only with PETSc

Incomplete Cholesky factorization (fill-in coefficient set to 0). This is variant of the ILU0 preconditioner dedicated to symmetric positive definite system

CS_PARAM_PRECOND_MUMPS 

Only with MUMPS

Sparse direct solver available with the MUMPS library and used as a preconditioner

CS_PARAM_PRECOND_HPDDM 

Only with PETSc

Use HPDDM and Geneo approach for coarse space with PETSc

CS_PARAM_PRECOND_POLY1 

Neumann polynomial preconditioning. Polynoms of order 1.

CS_PARAM_PRECOND_POLY2 

Neumann polynomial preconditioning. Polynoms of order 2.

CS_PARAM_PRECOND_SSOR 

Symmetric Successive OverRelaxations (can be seen as a symmetric Gauss-Seidel preconditioner)

CS_PARAM_N_PRECOND_TYPES 

◆ cs_param_resnorm_type_t

Way to renormalize (or not) the residual arising during the resolution of linear systems

Enumerator
CS_PARAM_RESNORM_NONE 

No renormalization

CS_PARAM_RESNORM_NORM2_RHS 

Renormalization based on the Euclidean norm of the right-hand side

CS_PARAM_RESNORM_WEIGHTED_RHS 

Renormalization based on a weighted Euclidean norm of the right-hand side

CS_PARAM_RESNORM_FILTERED_RHS 

Renormalization based on an Euclidean norm of a selection of the right-hand side (penalized terms are filtered)

CS_PARAM_N_RESNORM_TYPES 

◆ cs_param_solver_class_t

Class of iterative solvers to consider for solver the linear system.

Enumerator
CS_PARAM_SOLVER_CLASS_CS 

Iterative solvers available in code_saturne

CS_PARAM_SOLVER_CLASS_HYPRE 

Solvers available in HYPRE through the PETSc library

CS_PARAM_SOLVER_CLASS_MUMPS 

Solvers available with MUMPS (without the PETSc interface)

CS_PARAM_SOLVER_CLASS_PETSC 

Solvers available in PETSc. Please notice that the MUMPS solver can be handled within PETSc if the installation of PETSc includes the MUMPS library

CS_PARAM_N_SOLVER_CLASSES 

Iterative solvers available in code_saturne

◆ cs_param_solver_type_t

Type of solver used to solve a linear system. Some of the mentionned solver are available only if the PETSc library is linked with code_saturne.

Enumerator
CS_PARAM_SOLVER_NONE 

No iterative solver (equivalent to a "preonly" choice in PETSc)

CS_PARAM_SOLVER_AMG 

Algebraic multigrid solver (additional options may be set using cs_param_amg_type_t)

CS_PARAM_SOLVER_BICGS 

Bi-Conjugate gradient (useful for non-symmetric systems)

CS_PARAM_SOLVER_BICGS2 

Stabilized Bi-Conjugate gradient (useful for non-symmetric systems)

CS_PARAM_SOLVER_CG 

Conjugate Gradient (solver of choice for symmetric positive definite systems)

CS_PARAM_SOLVER_CR3 

3-layer conjugate residual (can handle non-symmetric systems)

CS_PARAM_SOLVER_FCG 

Flexible Conjugate Gradient (variant of the CG when the preconditioner may change from one iteration to another. For instance when using an AMG as preconditioner)

CS_PARAM_SOLVER_FGMRES 

Only with PETSc

Flexible Generalized Minimal RESidual

CS_PARAM_SOLVER_GAUSS_SEIDEL 

Gauss-Seidel

CS_PARAM_SOLVER_GCR 

Generalized conjugate residual (flexible iterative solver for symmetric or non-symmetric system)

CS_PARAM_SOLVER_GMRES 

Only with PETSc

Generalized Minimal RESidual

CS_PARAM_SOLVER_JACOBI 

Jacobi (diagonal relaxation)

CS_PARAM_SOLVER_MUMPS 

Only with PETSc or MUMPS

MUMPS sparse direct solver

CS_PARAM_SOLVER_SYM_GAUSS_SEIDEL 

Symmetric Gauss-Seidel

CS_PARAM_SOLVER_USER_DEFINED 

User-defined iterative solver. It relies on the implementation of the the function cs_user_sles_it_solver()

CS_PARAM_N_SOLVER_TYPES 

◆ cs_param_space_scheme_t

Type of numerical scheme for the discretization in space.

How is defined the degree of freedom.

Enumerator
CS_SPACE_SCHEME_LEGACY 

Cell-centered Finite Volume Two-Point Flux

CS_SPACE_SCHEME_CDOVB 

CDO scheme with vertex-based positioning

CS_SPACE_SCHEME_CDOVCB 

CDO scheme with vertex+cell-based positioning

CS_SPACE_SCHEME_CDOEB 

CDO scheme with edge-based positioning

CS_SPACE_SCHEME_CDOFB 

CDO scheme with face-based positioning

CS_SPACE_SCHEME_CDOCB 

CDO scheme with cell-based positioning (potential at cells and fluxes at faces)

CS_SPACE_SCHEME_HHO_P0 

Hybrid High Order (HHO) schemes HHO scheme with face-based positioning (lowest order)

CS_SPACE_SCHEME_HHO_P1 

Hybrid High Order (HHO) schemes HHO scheme with face-based positioning (k=1 up to order 3)

CS_SPACE_SCHEME_HHO_P2 

Hybrid High Order (HHO) schemes HHO scheme with face-based positioning (k=2 up to order 4)

CS_SPACE_SCHEME_MACFB 

MAC scheme with face-based positioning faces)

CS_SPACE_N_SCHEMES 

◆ cs_param_time_scheme_t

Type of numerical scheme for the discretization in time

Enumerator
CS_TIME_SCHEME_STEADY 

No time scheme. Steady-state computation.

CS_TIME_SCHEME_EULER_IMPLICIT 

fully implicit (forward Euler/theta-scheme = 1)

CS_TIME_SCHEME_EULER_EXPLICIT 

fully explicit (backward Euler/theta-scheme = 0)

CS_TIME_SCHEME_CRANKNICO 

Crank-Nicolson (theta-scheme = 0.5)

CS_TIME_SCHEME_THETA 

theta-scheme

CS_TIME_SCHEME_BDF2 

theta-scheme

CS_TIME_N_SCHEMES 

Function Documentation

◆ cs_param_get_advection_extrapol_name()

const char * cs_param_get_advection_extrapol_name ( cs_param_advection_extrapol_t  extrapol)

Get the label associated to the extrapolation used for the advection field.

Parameters
[in]extrapoltype of extrapolation for the advection field
Returns
the associated label

◆ cs_param_get_advection_form_name()

const char * cs_param_get_advection_form_name ( cs_param_advection_form_t  adv_form)

Get the label associated to the advection formulation.

Parameters
[in]adv_formtype of advection formulation
Returns
the associated label

◆ cs_param_get_advection_scheme_name()

const char * cs_param_get_advection_scheme_name ( cs_param_advection_scheme_t  scheme)

Get the label of the advection scheme.

Parameters
[in]schemetype of advection scheme
Returns
the associated advection scheme label

◆ cs_param_get_advection_strategy_name()

const char * cs_param_get_advection_strategy_name ( cs_param_advection_strategy_t  adv_stra)

Get the label associated to the advection strategy.

Parameters
[in]adv_stratype of advection strategy
Returns
the associated label

◆ cs_param_get_bc_enforcement_name()

const char * cs_param_get_bc_enforcement_name ( cs_param_bc_enforce_t  type)

Get the name of the type of enforcement of the boundary condition.

Parameters
[in]typetype of enforcement of boundary conditions
Returns
the associated name

◆ cs_param_get_bc_name()

const char * cs_param_get_bc_name ( cs_param_bc_type_t  type)

Get the name of the type of boundary condition.

Parameters
[in]typetype of boundary condition
Returns
the associated bc name

◆ cs_param_get_dotprod_type_name()

const char * cs_param_get_dotprod_type_name ( cs_param_dotprod_type_t  dp_type)

Get the name of the type of dot product to apply.

Parameters
[in]dp_typetype of dot product
Returns
the associated type name

◆ cs_param_get_nl_algo_label()

const char * cs_param_get_nl_algo_label ( cs_param_nl_algo_t  algo)

Get the label (short name) of the non-linear algorithm.

Parameters
[in]algotype of algorithm
Returns
the associated algorithm label

◆ cs_param_get_nl_algo_name()

const char * cs_param_get_nl_algo_name ( cs_param_nl_algo_t  algo)

Get the name of the non-linear algorithm.

Parameters
[in]algotype of algorithm
Returns
the associated algorithm name

◆ cs_param_get_precond_block_name()

const char * cs_param_get_precond_block_name ( cs_param_precond_block_t  type)

Get the name of the type of block preconditioning.

Parameters
[in]typetype of block preconditioning
Returns
the associated type name

◆ cs_param_get_precond_name()

const char * cs_param_get_precond_name ( cs_param_precond_type_t  precond)

Get the name of the preconditioner.

Parameters
[in]precondtype of preconditioner
Returns
the associated preconditioner name

◆ cs_param_get_solver_name()

const char * cs_param_get_solver_name ( cs_param_solver_type_t  solver)

Get the name of the solver.

Parameters
[in]solvertype of iterative solver
Returns
the associated solver name

◆ cs_param_get_space_scheme_name()

const char * cs_param_get_space_scheme_name ( cs_param_space_scheme_t  scheme)

Get the name of the space discretization scheme.

Parameters
[in]schemetype of space scheme
Returns
the associated space scheme name

◆ cs_param_get_time_scheme_name()

const char * cs_param_get_time_scheme_name ( cs_param_time_scheme_t  scheme)

Get the name of the time discretization scheme.

Parameters
[in]schemetype of time scheme
Returns
the associated time scheme name

◆ cs_param_space_scheme_is_face_based()

bool cs_param_space_scheme_is_face_based ( cs_param_space_scheme_t  scheme)

Return true if the space scheme has degrees of freedom on faces, otherwise false.

Parameters
[in]schemetype of space scheme
Returns
true or false

Variable Documentation

◆ cs_med_sepline

const char cs_med_sepline[50]
extern

◆ cs_sep_h1

const char cs_sep_h1[80]
extern

◆ cs_sep_h2

const char cs_sep_h2[80]
extern

◆ cs_sepline

const char cs_sepline[80]
extern