My Project
programmer's documentation
Data Structures | Typedefs | Enumerations
cs_param.h File Reference
#include "cs_defs.h"
Include dependency graph for cs_param.h:

Go to the source code of this file.

Data Structures

struct  cs_param_sles_t
 Structure storing all metadata related to the resolution of a linear system with an iterative solver. More...
 

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 compact, void *input, cs_real_t *retval)
 Generic function pointer for an analytic function elt_ids is optional. If not NULL, it enables to access in coords at the right location and the same thing to fill retval if compact is set to false. More...
 
typedef void() cs_time_func_t(int time_iter, double time, void *input, cs_real_t *retval)
 Function which defines the evolution of a quantity according to the number of iteration already done, 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_CDOFB,
  CS_SPACE_SCHEME_HHO_P0, CS_SPACE_SCHEME_HHO_P1, CS_SPACE_SCHEME_HHO_P2, 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 }
 

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_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_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_MIX_CENTERED_UPWIND,
  CS_PARAM_ADVECTION_SCHEME_SAMARSKII, CS_PARAM_ADVECTION_SCHEME_SG, CS_PARAM_ADVECTION_SCHEME_UPWIND, CS_PARAM_N_ADVECTION_SCHEMES
}
 

Settings for the boundary conditions

enum  cs_param_bc_type_t {
  CS_PARAM_BC_HMG_DIRICHLET, CS_PARAM_BC_DIRICHLET, CS_PARAM_BC_HMG_NEUMANN, CS_PARAM_BC_NEUMANN,
  CS_PARAM_BC_ROBIN, CS_PARAM_BC_SLIDING, 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
}
 

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

enum  cs_param_sles_class_t { CS_PARAM_SLES_CLASS_CS, CS_PARAM_SLES_CLASS_PETSC, CS_PARAM_SLES_N_CLASSES }
 Class of iterative solvers to consider for solver the linear system. More...
 
enum  cs_param_amg_type_t {
  CS_PARAM_AMG_NONE, CS_PARAM_AMG_HYPRE_BOOMER, CS_PARAM_AMG_PETSC_GAMG, CS_PARAM_AMG_PETSC_PCMG,
  CS_PARAM_AMG_HOUSE_V, CS_PARAM_AMG_HOUSE_K, CS_PARAM_N_AMG_TYPES
}
 
enum  cs_param_precond_type_t {
  CS_PARAM_PRECOND_NONE, CS_PARAM_PRECOND_DIAG, CS_PARAM_PRECOND_BJACOB, CS_PARAM_PRECOND_POLY1,
  CS_PARAM_PRECOND_POLY2, CS_PARAM_PRECOND_SSOR, CS_PARAM_PRECOND_ILU0, CS_PARAM_PRECOND_ICC0,
  CS_PARAM_PRECOND_AMG, CS_PARAM_PRECOND_AMG_BLOCK, CS_PARAM_PRECOND_AS, CS_PARAM_N_PRECOND_TYPES
}
 
enum  cs_param_itsol_type_t {
  CS_PARAM_ITSOL_AMG, CS_PARAM_ITSOL_BICG, CS_PARAM_ITSOL_BICGSTAB2, CS_PARAM_ITSOL_CG,
  CS_PARAM_ITSOL_CR3, CS_PARAM_ITSOL_FCG, CS_PARAM_ITSOL_GAUSS_SEIDEL, CS_PARAM_ITSOL_GMRES,
  CS_PARAM_ITSOL_JACOBI, CS_PARAM_ITSOL_MINRES, CS_PARAM_ITSOL_SYM_GAUSS_SEIDEL, CS_PARAM_N_ITSOL_TYPES
}
 
enum  cs_param_resnorm_type_t {
  CS_PARAM_RESNORM_NONE, CS_PARAM_RESNORM_VOLTOT, CS_PARAM_RESNORM_WEIGHTED_RHS, CS_PARAM_RESNORM_MAT_DIAG,
  CS_PARAM_N_RESNORM_TYPES
}
 
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_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_solver_name (cs_param_itsol_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_amg_type_name (cs_param_amg_type_t type)
 Get the name of the type of algebraic multigrid (AMG) More...
 

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 compact, void *input, cs_real_t *retval)

Generic function pointer for an analytic function elt_ids is optional. If not NULL, it enables to access in coords at the right location and the same thing to fill retval if compact is set to false.

Parameters
[in]timewhen ?
[in]n_eltsnumber of elements to consider
[in]elt_idslist of elements ids (to access coords and fill)
[in]coordswhere ?
[in]compacttrue:no indirection, false:indirection for filling
[in]inputpointer to a structure cast on-the-fly (may be NULL)
[in,out]retvalresult of the function

◆ cs_time_func_t

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

Function which defines the evolution of a quantity according to the number of iteration already done, the current time and any structure given as a parameter.

Parameters
[in]time_itercurrent number of iterations
[in]timevalue of the time at the end of the last iteration
[in]inputpointer to a structure cast on-the-fly
[in]retvalresult of the evaluation

Enumeration Type Documentation

◆ 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_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

CS_PARAM_ADVECTION_SCHEME_CIP_CW 

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

CS_PARAM_ADVECTION_SCHEME_MIX_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_amg_type_t

Type of AMG (Algebraic MultiGrid) algorithm to use (either as a preconditionnerwith or a solver).

Enumerator
CS_PARAM_AMG_NONE 

No specified algorithm

CS_PARAM_AMG_HYPRE_BOOMER 

Boomer algorithm from Hypre library

CS_PARAM_AMG_PETSC_GAMG 

GAMG algorithm from PETSc

CS_PARAM_AMG_PETSC_PCMG 

preconditionned MG algorithm from PETSc

CS_PARAM_AMG_HOUSE_V 

In-house algorithm with V-cycle

CS_PARAM_AMG_HOUSE_K 

In-house algorithm with K-cycle

CS_PARAM_N_AMG_TYPES 

◆ 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_PARAM_BC_HMG_DIRICHLET 

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

CS_PARAM_BC_DIRICHLET 

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

CS_PARAM_BC_HMG_NEUMANN 

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

CS_PARAM_BC_NEUMANN 

Neumann conditions. The value of the flux of variable is set to the user requirements.

CS_PARAM_BC_ROBIN 

Robin conditions.

CS_PARAM_BC_SLIDING 

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

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_itsol_type_t

Type of iterative solver to use to inverse the linear system. CS_PARAM_ITSOL_CR3 is available only inside the Code_Saturne framework.

Enumerator
CS_PARAM_ITSOL_AMG 

Algebraic MultiGrid

CS_PARAM_ITSOL_BICG 

Bi-Conjuguate gradient

CS_PARAM_ITSOL_BICGSTAB2 

Stabilized Bi-Conjuguate gradient

CS_PARAM_ITSOL_CG 

Conjuguate Gradient

CS_PARAM_ITSOL_CR3 

3-layer conjugate residual

CS_PARAM_ITSOL_FCG 

Flexible Conjuguate Gradient

CS_PARAM_ITSOL_GAUSS_SEIDEL 

Gauss-Seidel

CS_PARAM_ITSOL_GMRES 

Generalized Minimal RESidual

CS_PARAM_ITSOL_JACOBI 

Jacobi

CS_PARAM_ITSOL_MINRES 

Mininal Residual

CS_PARAM_ITSOL_SYM_GAUSS_SEIDEL 

Symetric Gauss-Seidel

CS_PARAM_N_ITSOL_TYPES 

◆ cs_param_precond_type_t

Type of preconditionner to use with the iterative solver. Some preconditionners as CS_PARAM_PRECOND_ILU0, CS_PARAM_PRECOND_ICC0, CS_PARAM_PRECOND_AS and CS_PARAM_PRECOND_AMG_BLOCK are available only with the PETSc interface.

Enumerator
CS_PARAM_PRECOND_NONE 

No preconditioning

CS_PARAM_PRECOND_DIAG 

Diagonal (or Jacobi) preconditioning

CS_PARAM_PRECOND_BJACOB 

Block Jacobi

CS_PARAM_PRECOND_POLY1 

Neumann polynomial preconditioning (Order 1)

CS_PARAM_PRECOND_POLY2 

Neumann polynomial preconditioning (Order 2)

CS_PARAM_PRECOND_SSOR 

Symmetric Successive OverRelaxations

CS_PARAM_PRECOND_ILU0 

Incomplete LU factorization

CS_PARAM_PRECOND_ICC0 

Incomplete Cholesky factorization

CS_PARAM_PRECOND_AMG 

Algebraic MultiGrid

CS_PARAM_PRECOND_AMG_BLOCK 

Algebraic MultiGrid by block

CS_PARAM_PRECOND_AS 

Additive Schwarz method

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_VOLTOT 

Renormalization based on the volume of the computational domain

CS_PARAM_RESNORM_WEIGHTED_RHS 

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

CS_PARAM_RESNORM_MAT_DIAG 
CS_PARAM_N_RESNORM_TYPES 

◆ cs_param_sles_class_t

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

Enumerator
CS_PARAM_SLES_CLASS_CS 

Iterative solvers available in Code_Saturne

CS_PARAM_SLES_CLASS_PETSC 

Solvers available in Code_Saturne

CS_PARAM_SLES_N_CLASSES 

◆ 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 positionning

CS_SPACE_SCHEME_CDOVCB 

CDO scheme with vertex+cell-based positionning

CS_SPACE_SCHEME_CDOFB 

CDO scheme with face-based positionning

CS_SPACE_SCHEME_HHO_P0 

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

CS_SPACE_SCHEME_HHO_P1 

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

CS_SPACE_SCHEME_HHO_P2 

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

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_N_SCHEMES 

Function Documentation

◆ cs_param_get_amg_type_name()

const char* cs_param_get_amg_type_name ( cs_param_amg_type_t  type)

Get the name of the type of algebraic multigrid (AMG)

Parameters
[in]typetype of AMG
Returns
the associated type name

◆ 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_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_itsol_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