7.0
general documentation
cs_equation_param.h File Reference

Structure and routines handling the specific settings related to a cs_equation_t structure. More...

#include "cs_advection_field.h"
#include "cs_param_cdo.h"
#include "cs_param_sles.h"
#include "cs_hodge.h"
#include "cs_property.h"
#include "cs_xdef.h"
+ Include dependency graph for cs_equation_param.h:

Go to the source code of this file.

Data Structures

struct  cs_equation_param_t
 Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources. More...
 

Macros

Flags specifying which term is needed for an equation.
#define CS_EQUATION_LOCKED   (1 << 0) /* 1 */
 Parameters for setting an equation are not modifiable anymore. More...
 
#define CS_EQUATION_UNSTEADY   (1 << 1) /* 2 */
 Unsteady term is needed. More...
 
#define CS_EQUATION_CONVECTION   (1 << 2) /* 4 */
 Convection term is needed. More...
 
#define CS_EQUATION_DIFFUSION   (1 << 3) /* 8 */
 Diffusion term is needed. A scalar-/vector-valued Laplacian with div .grad. More...
 
#define CS_EQUATION_CURLCURL   (1 << 4) /* 16 */
 The term related to the curl-curl operator is needed. More...
 
#define CS_EQUATION_GRADDIV   (1 << 5) /* 32 */
 The term related to the grad-div operator is needed. More...
 
#define CS_EQUATION_REACTION   (1 << 6) /* 64 */
 Reaction term is needed. More...
 
#define CS_EQUATION_FORCE_VALUES   (1 << 7) /* 128 */
 Add an algebraic manipulation to set the value of a given set of interior degrees of freedom. More...
 
#define CS_EQUATION_USER_HOOK   (1 << 8) /* 256 */
 Activate a user hook to get a fine control of the discretization process during the cellwise building of the linear system Need to match the cs_equation_user_hook_t prototype. More...
 
Flags specifying which extra operation is needed for an equation.
#define CS_EQUATION_POST_BALANCE   (1 << 0) /* 1 */
 Compute and postprocess the equation balance. More...
 
#define CS_EQUATION_POST_PECLET   (1 << 1) /* 2 */
 Compute and postprocess the Peclet number. More...
 
#define CS_EQUATION_POST_UPWIND_COEF   (1 << 2) /* 4 */
 Postprocess the value of the upwinding coefficient. More...
 
#define CS_EQUATION_POST_NORMAL_FLUX   (1 << 3) /* 8 */
 Postprocess the value of the normal flux across the boundary faces. More...
 
Flags to handle the enforcement of degrees of freedom (DoFs)
#define CS_EQUATION_ENFORCE_BY_CELLS   (1 << 0) /* 1 */
 Definition of a selection of DoFs to enforce using a cell selection. More...
 
#define CS_EQUATION_ENFORCE_BY_DOFS   (1 << 1) /* 2 */
 Definition of a selection of DoFs. More...
 
#define CS_EQUATION_ENFORCE_BY_REFERENCE_VALUE   (1 << 2) /* 4 */
 Assign to all the selected DoFs the same value. This value is stored in enforcement_ref_value. More...
 

Enumerations

enum  cs_equation_type_t {
  CS_EQUATION_TYPE_GROUNDWATER, CS_EQUATION_TYPE_MAXWELL, CS_EQUATION_TYPE_NAVSTO, CS_EQUATION_TYPE_PREDEFINED,
  CS_EQUATION_TYPE_THERMAL, CS_EQUATION_TYPE_SOLIDIFICATION, CS_EQUATION_TYPE_USER, CS_EQUATION_N_TYPES
}
 Type of equations managed by the solver. More...
 
enum  cs_equation_key_t {
  CS_EQKEY_ADV_EXTRAPOL, CS_EQKEY_ADV_FORMULATION, CS_EQKEY_ADV_SCHEME, CS_EQKEY_ADV_STRATEGY,
  CS_EQKEY_ADV_UPWIND_PORTION, CS_EQKEY_AMG_TYPE, CS_EQKEY_BC_ENFORCEMENT, CS_EQKEY_BC_QUADRATURE,
  CS_EQKEY_BC_STRONG_PENA_COEFF, CS_EQKEY_BC_WEAK_PENA_COEFF, CS_EQKEY_DO_LUMPING, CS_EQKEY_DOF_REDUCTION,
  CS_EQKEY_EXTRA_OP, CS_EQKEY_HODGE_DIFF_ALGO, CS_EQKEY_HODGE_DIFF_COEF, CS_EQKEY_HODGE_TIME_ALGO,
  CS_EQKEY_HODGE_REAC_ALGO, CS_EQKEY_ITSOL, CS_EQKEY_ITSOL_EPS, CS_EQKEY_ITSOL_MAX_ITER,
  CS_EQKEY_ITSOL_RESNORM_TYPE, CS_EQKEY_OMP_ASSEMBLY_STRATEGY, CS_EQKEY_PRECOND, CS_EQKEY_SLES_VERBOSITY,
  CS_EQKEY_SOLVER_FAMILY, CS_EQKEY_SPACE_SCHEME, CS_EQKEY_TIME_SCHEME, CS_EQKEY_TIME_THETA,
  CS_EQKEY_VERBOSITY, CS_EQKEY_N_KEYS
}
 List of available keys for setting the parameters of an equation. More...
 

Functions

cs_equation_param_tcs_equation_create_param (const char *name, cs_equation_type_t type, int dim, cs_param_bc_type_t default_bc)
 Create a cs_equation_param_t structure. More...
 
void cs_equation_copy_param_from (const cs_equation_param_t *ref, cs_equation_param_t *dst, bool copy_fid)
 Copy the settings from one cs_equation_param_t structure to another one. The name is not copied. More...
 
void cs_equation_clear_param (cs_equation_param_t *eqp)
 Free the contents of a cs_equation_param_t. More...
 
cs_equation_param_tcs_equation_free_param (cs_equation_param_t *eqp)
 Free a cs_equation_param_t. More...
 
void cs_equation_set_param (cs_equation_param_t *eqp, cs_equation_key_t key, const char *keyval)
 Set a parameter attached to a keyname in a cs_equation_param_t structure. More...
 
void cs_equation_param_set_sles (cs_equation_param_t *eqp)
 Set parameters for initializing SLES structures used for the resolution of the linear system. Settings are related to this equation. More...
 
void cs_equation_param_last_stage (cs_equation_param_t *eqp)
 Last modification of the cs_equation_param_t structure before launching the computation. More...
 
void cs_equation_param_log (const cs_equation_param_t *eqp)
 Print the detail of a cs_equation_param_t structure. More...
 
bool cs_equation_param_has_robin_bc (const cs_equation_param_t *eqp)
 Ask if the parameter settings of the equation has requested the treatment of Robin boundary conditions. More...
 
cs_xdef_tcs_equation_add_ic_by_value (cs_equation_param_t *eqp, const char *z_name, cs_real_t *val)
 Define the initial condition for the unknown related to this equation This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here a constant value is set to all the entities belonging to the given mesh location. More...
 
cs_xdef_tcs_equation_add_ic_by_qov (cs_equation_param_t *eqp, const char *z_name, double quantity)
 Define the initial condition for the unknown related to this equation This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the value related to all the entities belonging to the given mesh location is such that the integral over these cells returns the requested quantity. More...
 
cs_xdef_tcs_equation_add_ic_by_analytic (cs_equation_param_t *eqp, const char *z_name, cs_analytic_func_t *analytic, void *input)
 Define the initial condition for the unknown related to this equation. This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set according to an analytical function. More...
 
void cs_equation_add_xdef_bc (cs_equation_param_t *eqp, cs_xdef_t *xdef)
 Set a boundary condition from an existing cs_xdef_t structure The lifecycle of the cs_xdef_t structure is now managed by the current cs_equation_param_t structure. More...
 
cs_xdef_tcs_equation_add_bc_by_value (cs_equation_param_t *eqp, const cs_param_bc_type_t bc_type, const char *z_name, cs_real_t *values)
 Define and initialize a new structure to set a boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t. More...
 
cs_xdef_tcs_equation_add_bc_by_array (cs_equation_param_t *eqp, const cs_param_bc_type_t bc_type, const char *z_name, cs_flag_t loc, cs_real_t *array, bool is_owner, cs_lnum_t *index)
 Define and initialize a new structure to set a boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t. More...
 
cs_xdef_tcs_equation_add_bc_by_analytic (cs_equation_param_t *eqp, const cs_param_bc_type_t bc_type, const char *z_name, cs_analytic_func_t *analytic, void *input)
 Define and initialize a new structure to set a boundary condition related to the given equation param structure ml_name corresponds to the name of a pre-existing cs_mesh_location_t. More...
 
cs_xdef_tcs_equation_add_bc_by_dof_func (cs_equation_param_t *eqp, const cs_param_bc_type_t bc_type, const char *z_name, cs_flag_t loc_flag, cs_dof_func_t *func, void *input)
 Define and initialize a new structure to set a boundary condition related to the given cs_equation_param_t structure ml_name corresponds to the name of a pre-existing cs_mesh_location_t. More...
 
cs_xdef_tcs_equation_find_bc (cs_equation_param_t *eqp, const char *z_name)
 Return pointer to existing boundary condition definition structure for the given equation param structure and zone. More...
 
void cs_equation_remove_bc (cs_equation_param_t *eqp, const char *z_name)
 Remove boundary condition from the given equation param structure for a given zone. More...
 
void cs_equation_add_sliding_condition (cs_equation_param_t *eqp, const char *z_name)
 Define and initialize a new structure to set a sliding boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t. More...
 
void cs_equation_add_diffusion (cs_equation_param_t *eqp, cs_property_t *property)
 Associate a new term related to the Laplacian operator for the equation associated to the given cs_equation_param_t structure Laplacian means div-grad (either for vector-valued or scalar-valued equations) More...
 
void cs_equation_add_curlcurl (cs_equation_param_t *eqp, cs_property_t *property, int inversion)
 Associate a new term related to the curl-curl operator for the equation associated to the given cs_equation_param_t structure. More...
 
void cs_equation_add_graddiv (cs_equation_param_t *eqp, cs_property_t *property)
 Associate a new term related to the grad-div operator for the equation associated to the given cs_equation_param_t structure. More...
 
void cs_equation_add_time (cs_equation_param_t *eqp, cs_property_t *property)
 Associate a new term related to the time derivative operator for the equation associated to the given cs_equation_param_t structure. More...
 
void cs_equation_add_advection (cs_equation_param_t *eqp, cs_adv_field_t *adv_field)
 Associate a new term related to the advection operator for the equation associated to the given cs_equation_param_t structure. More...
 
void cs_equation_add_advection_scaling_property (cs_equation_param_t *eqp, cs_property_t *property)
 Associate a scaling property to the advection. More...
 
int cs_equation_add_reaction (cs_equation_param_t *eqp, cs_property_t *property)
 Associate a new term related to the reaction operator for the equation associated to the given cs_equation_param_t structure. More...
 
cs_xdef_tcs_equation_add_source_term_by_val (cs_equation_param_t *eqp, const char *z_name, cs_real_t *val)
 Add a new source term by initializing a cs_xdef_t structure. Case of a definition by a constant value. More...
 
cs_xdef_tcs_equation_add_source_term_by_analytic (cs_equation_param_t *eqp, const char *z_name, cs_analytic_func_t *func, void *input)
 Add a new source term by initializing a cs_xdef_t structure. Case of a definition by an analytical function. More...
 
cs_xdef_tcs_equation_add_source_term_by_dof_func (cs_equation_param_t *eqp, const char *z_name, cs_flag_t loc_flag, cs_dof_func_t *func, void *input)
 Add a new source term by initializing a cs_xdef_t structure. Case of a definition by a DoF function. More...
 
cs_xdef_tcs_equation_add_source_term_by_array (cs_equation_param_t *eqp, const char *z_name, cs_flag_t loc, cs_real_t *array, bool is_owner, cs_lnum_t *index)
 Add a new source term by initializing a cs_xdef_t structure. Case of a definition by an array. More...
 
cs_xdef_tcs_equation_add_volume_mass_injection_by_value (cs_equation_param_t *eqp, const char *z_name, double *val)
 Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using a constant value. More...
 
cs_xdef_tcs_equation_add_volume_mass_injection_by_qov (cs_equation_param_t *eqp, const char *z_name, double *quantity)
 Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using a constant quantity distributed over the associated zone's volume. More...
 
cs_xdef_tcs_equation_add_volume_mass_injection_by_analytic (cs_equation_param_t *eqp, const char *z_name, cs_analytic_func_t *func, void *input)
 Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using an analytical function. More...
 
void cs_equation_enforce_vertex_dofs (cs_equation_param_t *eqp, cs_lnum_t n_elts, const cs_lnum_t elt_ids[], const cs_real_t ref_value[], const cs_real_t elt_values[])
 Add an enforcement of the value of degrees of freedom located at mesh vertices. The spatial discretization scheme for the given equation has to be CDO-Vertex based or CDO-Vertex+Cell-based schemes. More...
 
void cs_equation_enforce_value_on_cell_selection (cs_equation_param_t *eqp, cs_lnum_t n_elts, const cs_lnum_t elt_ids[], const cs_real_t ref_value[], const cs_real_t elt_values[])
 Add an enforcement of the value related to the degrees of freedom associated to the list of selected cells. More...
 

Detailed Description

Structure and routines handling the specific settings related to a cs_equation_t structure.

Macro Definition Documentation

◆ CS_EQUATION_CONVECTION

#define CS_EQUATION_CONVECTION   (1 << 2) /* 4 */

Convection term is needed.

◆ CS_EQUATION_CURLCURL

#define CS_EQUATION_CURLCURL   (1 << 4) /* 16 */

The term related to the curl-curl operator is needed.

◆ CS_EQUATION_DIFFUSION

#define CS_EQUATION_DIFFUSION   (1 << 3) /* 8 */

Diffusion term is needed. A scalar-/vector-valued Laplacian with div .grad.

◆ CS_EQUATION_ENFORCE_BY_CELLS

#define CS_EQUATION_ENFORCE_BY_CELLS   (1 << 0) /* 1 */

Definition of a selection of DoFs to enforce using a cell selection.

◆ CS_EQUATION_ENFORCE_BY_DOFS

#define CS_EQUATION_ENFORCE_BY_DOFS   (1 << 1) /* 2 */

Definition of a selection of DoFs.

◆ CS_EQUATION_ENFORCE_BY_REFERENCE_VALUE

#define CS_EQUATION_ENFORCE_BY_REFERENCE_VALUE   (1 << 2) /* 4 */

Assign to all the selected DoFs the same value. This value is stored in enforcement_ref_value.

◆ CS_EQUATION_FORCE_VALUES

#define CS_EQUATION_FORCE_VALUES   (1 << 7) /* 128 */

Add an algebraic manipulation to set the value of a given set of interior degrees of freedom.

◆ CS_EQUATION_GRADDIV

#define CS_EQUATION_GRADDIV   (1 << 5) /* 32 */

The term related to the grad-div operator is needed.

◆ CS_EQUATION_LOCKED

#define CS_EQUATION_LOCKED   (1 << 0) /* 1 */

Parameters for setting an equation are not modifiable anymore.

◆ CS_EQUATION_POST_BALANCE

#define CS_EQUATION_POST_BALANCE   (1 << 0) /* 1 */

Compute and postprocess the equation balance.

◆ CS_EQUATION_POST_NORMAL_FLUX

#define CS_EQUATION_POST_NORMAL_FLUX   (1 << 3) /* 8 */

Postprocess the value of the normal flux across the boundary faces.

◆ CS_EQUATION_POST_PECLET

#define CS_EQUATION_POST_PECLET   (1 << 1) /* 2 */

Compute and postprocess the Peclet number.

◆ CS_EQUATION_POST_UPWIND_COEF

#define CS_EQUATION_POST_UPWIND_COEF   (1 << 2) /* 4 */

Postprocess the value of the upwinding coefficient.

◆ CS_EQUATION_REACTION

#define CS_EQUATION_REACTION   (1 << 6) /* 64 */

Reaction term is needed.

◆ CS_EQUATION_UNSTEADY

#define CS_EQUATION_UNSTEADY   (1 << 1) /* 2 */

Unsteady term is needed.

◆ CS_EQUATION_USER_HOOK

#define CS_EQUATION_USER_HOOK   (1 << 8) /* 256 */

Activate a user hook to get a fine control of the discretization process during the cellwise building of the linear system Need to match the cs_equation_user_hook_t prototype.

Enumeration Type Documentation

◆ cs_equation_key_t

List of available keys for setting the parameters of an equation.

Enumerator
CS_EQKEY_ADV_EXTRAPOL 

Choice in the way to extrapolate the advection field when building the advection operator

CS_EQKEY_ADV_FORMULATION 

Kind of formulation of the advective term. Available choices are:

  • "conservative"
  • "non_conservative"
CS_EQKEY_ADV_SCHEME 

Type of numerical scheme for the advective term. The available choices depend on the space discretization scheme.

CS_EQKEY_ADV_STRATEGY 

Strategy used to handle the advection term

  • "fully_implicit" or "implicit" (default choice)
  • "linearized" or "implicit_linear"
  • "explicit" There is a difference between the two first choices when the advection term induces a non-linearity. In this situation, the first choice implies that a non-linear algorithm has to be used.
CS_EQKEY_ADV_UPWIND_PORTION 

Value between 0 and 1 specifying the portion of upwind added to a centered discretization.

CS_EQKEY_AMG_TYPE 

Specify which type of algebraic multigrid (AMG) to choose. Available choices are:

  • "none" –> (default) No predefined AMG solver
  • "boomer" –> Boomer AMG multigrid from the Hypre library
  • "gamg" –> GAMG multigrid from the PETSc library
  • "pcmg" –> MG multigrid as preconditioner from the PETSc library
  • "v_cycle" –> Code_Saturne's built-in multigrid with a V-cycle strategy
  • "k_cycle" –> Code_Saturne's built-in multigrid with a K-cycle strategy WARNING: For "boomer" and "gamg",one needs to install Code_Saturne with PETSc in this case
CS_EQKEY_BC_ENFORCEMENT 

Set the type of enforcement of the boundary conditions. Available choices are:

  • "algebraic": Modify the linear system so as to add the contribution of the Dirichlet in the right-hand side and replace the line associated to a Dirichlet BC by identity. This is a good choice for pure diffusinon or pure convection problem.
  • "penalization": Add a huge penalization coefficient on the diagonal term of the line related to DoFs associated a Dirichlet BC. The right-hand side is also modified to take into account this penalization. Be aware that it may worsen the matrix conditioning.
  • "weak": weak enforcement using the Nitsche method (there is also penalization term but the scaling is such that a moderate value (1-100) of the penalization coefficient is sufficient). This a good choice for convection/diffusion problem.
  • "weak_sym": Same as the "weak" option but with a modification so as to add a symmetric contribution to the system. If the problem yields a symmetric matrix. This choice is more relevant than "weak". This a good choice for a diffusion problem.

For HHO and CDO-Face based schemes, only the "penalization" and "algebraic" technique is available up to now.

CS_EQKEY_BC_QUADRATURE 

Set the quadrature algorithm used for evaluating integral quantities on faces or volumes. Available choices are:

  • "bary" used the barycenter approximation
  • "higher" used 4 Gauss points for approximating the integral
  • "highest" used 5 Gauss points for approximating the integral

Remark: "higher" and "highest" implies automatically a subdivision into tetrahedra of each cell.

CS_EQKEY_BC_STRONG_PENA_COEFF 

Set the value of the penalization coefficient when "penalization" is activated The default value is 1e12. cf. CS_PARAM_BC_ENFORCE_PENALIZED

CS_EQKEY_BC_WEAK_PENA_COEFF 

Set the value of the penalization coefficient when "weak" or "weak_sym" is activated. The default value is 100. cf. CS_PARAM_BC_ENFORCE_WEAK_NITSCHE or CS_PARAM_BC_ENFORCE_WEAK_SYM

CS_EQKEY_DO_LUMPING 
CS_EQKEY_DOF_REDUCTION 

Set how is defined each degree of freedom (DoF).

  • "de_rham" (default): 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
  • "average": DoF are defined as the average on the element
CS_EQKEY_EXTRA_OP 

Set the additional post-processing to perform. Available choices are:

  • "balance" post-process the balance result in each control volume for each main term of an equation (diffusion, convection, time...)
  • "peclet" post-process an estimation of the Peclet number in each cell
  • "upwind_coef" post-process an estimation of the upwinding coefficient
  • "normal_flux" post-process the normal flux across boundary faces
CS_EQKEY_HODGE_DIFF_ALGO 

Set the algorithm used for building the discrete Hodge operator used in the diffusion term. Available choices are:

  • "voronoi" –> leads to a diagonal discrete Hodge operator but is not consistent for all meshes. Require an "orthogonal" (or admissible) mesh;
  • "cost" –> (default for diffusion) is more robust (i.e. it handles more general meshes but is is less efficient)
  • "wbs" –> is robust and accurate but is limited to the reconstruction of potential-like degrees of freedom and needs a correct computation of the cell barycenter
CS_EQKEY_HODGE_DIFF_COEF 

This key is only useful if CS_EQKEY_HODGE_{TIME, DIFF, REAC}_ALGO is set to "cost". keyval is either a name or a value:

  • "dga" corresponds to the value $ 1./3. $
  • "sushi" corresponds to the value $1./\sqrt(3.)$
  • "gcr" corresponds to the value $1$.
  • or "1.5", "9" for instance
CS_EQKEY_HODGE_TIME_ALGO 

Set the algorithm used for building the discrete Hodge operator used in the unsteady term. Available choices are:

  • "voronoi" –> (default) leads to diagonal discrete Hodge operator. It acts as a mass lumping.
  • "wbs" is robust and accurate but is less efficient. It needs a correct computation of the cell barycenter
CS_EQKEY_HODGE_REAC_ALGO 

Set the algorithm used for building the discrete Hodge operator used in the reaction term. Available choices are:

  • "voronoi" –> leads to diagonal discrete Hodge operator (similar to a lumping).
  • "wbs" –> (default) is robust and accurate but is limited to the reconstruction of potential-like degrees of freedom and needs a correct computation of the cell barycenter
CS_EQKEY_ITSOL 

Specify the iterative solver for solving the linear system related to an equation. Avalaible choices are:

  • "jacobi" –> simpliest iterative solver
  • "gauss_seidel" –> Gauss-Seidel algorithm
  • "sym_gauss_seidel" –> Symmetric version of Gauss-Seidel algorithm; one backward and forward sweep
  • "cg" –> conjuguate gradient algorithm
  • "fcg" –> flexible version of the conjuguate gradient algorithm used when the preconditioner can change iteration by iteration
  • "bicg" –> Bi-CG algorithm (for non-symmetric linear systems)
  • "bicgstab2" –> BiCG-Stab2 algorithm (for non-symmetric linear systems)
  • "cr3" –> a 3-layer conjugate residual solver (when "cs" is chosen as the solver family)
  • "gmres" –> robust iterative solver. Not the best choice if the system is easy to solve
  • "amg" –> algebraic multigrid iterative solver. Good choice for a scalable solver related to symmetric positive definite system.
  • "minres" –> Solver of choice for symmetric indefinite systems
  • "mumps" –> Direct solver (very robust but memory consumming) via PETSc only. LU factorization.
  • "mumps_ldlt" –> Direct solver (very robust but memory consumming) via PETSc only. LDLT factorization.
CS_EQKEY_ITSOL_EPS 

Tolerance factor for stopping the iterative processus for solving the linear system related to an equation

  • Example: "1e-10"
CS_EQKEY_ITSOL_MAX_ITER 

Maximum number of iterations for solving the linear system

  • Example: "2000"
CS_EQKEY_ITSOL_RESNORM_TYPE 

Normalized or not the residual before testing if one continues iterating for solving the linear system. This normalization is performed before applying the boundary conditions to avoid handling the penalization of Dirichlet boundary conditions. If the RHS norm is equal to zero, then the "vol_tot" option is used for rescaling the residual.

Available choices are: "false" or "none" (default) "rhs" "weighted_rhs" or "weighted" "filtered_rhs" or "fieltered_rhs"

CS_EQKEY_OMP_ASSEMBLY_STRATEGY 

Choice of the way to perform the assembly when OpenMP is active Available choices are:

  • "atomic" or "critical"
CS_EQKEY_PRECOND 

Specify the preconditioner associated to an iterative solver. Available choices are:

  • "jacobi" or "diag" or "diagonal": diagonal preconditoner
  • "block_jacobi": Only with PETSc
  • "poly1": Neumann polynomial of order 1 (only with Code_Saturne)
  • "poly2": Neumann polynomial of order 2 (only with Code_Saturne)
  • "ssor": symmetric successive over-relaxation (only with PETSC)
  • "ilu0": incomplete LU factorization (only with PETSc)
  • "icc0": incomplete Cholesky factorization (for symmetric matrices and only with PETSc)
  • "amg": algebraic multigrid
  • "amg_block": algebraic multigrid by block (useful for vector-valued equations)
CS_EQKEY_SLES_VERBOSITY 

Level of details written by the code for the resolution of the linear system

  • Examples: "0", "1", "2" or higher
CS_EQKEY_SOLVER_FAMILY 
CS_EQKEY_SPACE_SCHEME 

Set the space discretization scheme. Available choices are:

  • "cdo_vb" or "cdovb" for CDO vertex-based scheme
  • "cdo_vcb" or "cdovcb" for CDO vertex+cell-based scheme
  • "cdo_fb" or "cdofb" for CDO face-based scheme
  • "cdo_eb" or "cdoeb" for CDO edge-based scheme
  • "hho_p1" for HHO schemes with $\mathbb{P}_1$ polynomial approximation
  • "hho_p2" for HHO schemes with $\mathbb{P}_2$ polynomial approximation
CS_EQKEY_TIME_SCHEME 

Set the scheme for the temporal discretization. Available choices are:

  • "euler_implicit" or "forward_euler" (1st order)
  • "euler_explicit" or "backward_euler" (1st order, conditional stability)
  • "theta_scheme": Time scheme encompassing several other schemes according to the value of theta. One recovers for instance "euler_implicit" with theta equal to "1", "euler_explicit" with "0" or "crank_nicolson" with theta equal to "0.5"
  • "crank_nicolson": Shortcut to set a theta scheme with theta equal to 0.5. This time discretization is second_order in time accurate.
  • "bdf2": 2nd order, implicit time scheme
CS_EQKEY_TIME_THETA 

Set a value betwwen 0 and 1 for the theta parameter when a "theta_scheme" is set for the time discretization. Only useful if CS_EQKEY_TIME_SCHEME is set to "theta_scheme".

  • Example: "0.75"
CS_EQKEY_VERBOSITY 

Set the level of details written by the code for an equation. The higher the more detailed information is displayed.

  • "0" (default)
  • "1" detailed setup resume and coarse grain timer stats
  • "2" fine grain for timer stats
CS_EQKEY_N_KEYS 

◆ cs_equation_type_t

Type of equations managed by the solver.

Enumerator
CS_EQUATION_TYPE_GROUNDWATER 

Equation related to the groundwater flow module

CS_EQUATION_TYPE_MAXWELL 

Equation related to the Maxwell module

CS_EQUATION_TYPE_NAVSTO 

Equation related to the resolution of the Navier-Stokes system

  • Example: momentum, prediction, correction, energy...
CS_EQUATION_TYPE_PREDEFINED 

Predefined equation (most part of the setting is already done)

  • Example: equation for the wall distance or ALE
CS_EQUATION_TYPE_THERMAL 

Equation related to the heat transfer

CS_EQUATION_TYPE_SOLIDIFICATION 

Equation related to the solidification module

CS_EQUATION_TYPE_USER 

User-defined equation

CS_EQUATION_N_TYPES 

Function Documentation

◆ cs_equation_add_advection()

void cs_equation_add_advection ( cs_equation_param_t eqp,
cs_adv_field_t adv_field 
)

Associate a new term related to the advection operator for the equation associated to the given cs_equation_param_t structure.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]adv_fieldpointer to a cs_adv_field_t structure

◆ cs_equation_add_advection_scaling_property()

void cs_equation_add_advection_scaling_property ( cs_equation_param_t eqp,
cs_property_t property 
)

Associate a scaling property to the advection.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_t structure

◆ cs_equation_add_bc_by_analytic()

cs_xdef_t* cs_equation_add_bc_by_analytic ( cs_equation_param_t eqp,
const cs_param_bc_type_t  bc_type,
const char *  z_name,
cs_analytic_func_t analytic,
void *  input 
)

Define and initialize a new structure to set a boundary condition related to the given equation param structure ml_name corresponds to the name of a pre-existing cs_mesh_location_t.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]bc_typetype of boundary condition to add
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]analyticpointer to an analytic function defining the value
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_bc_by_array()

cs_xdef_t* cs_equation_add_bc_by_array ( cs_equation_param_t eqp,
const cs_param_bc_type_t  bc_type,
const char *  z_name,
cs_flag_t  loc,
cs_real_t array,
bool  is_owner,
cs_lnum_t index 
)

Define and initialize a new structure to set a boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]bc_typetype of boundary condition to add
[in]z_namename of the related boundary zone
[in]locinformation to know where are located values
[in]arraypointer to an array
[in]is_ownertransfer the lifecycle to the cs_xdef_t structure (true or false)
[in]indexoptional pointer to the array index
Returns
a pointer to the new allocated cs_xdef_t structure

◆ cs_equation_add_bc_by_dof_func()

cs_xdef_t* cs_equation_add_bc_by_dof_func ( cs_equation_param_t eqp,
const cs_param_bc_type_t  bc_type,
const char *  z_name,
cs_flag_t  loc_flag,
cs_dof_func_t func,
void *  input 
)

Define and initialize a new structure to set a boundary condition related to the given cs_equation_param_t structure ml_name corresponds to the name of a pre-existing cs_mesh_location_t.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]bc_typetype of boundary condition to add
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]loc_flaglocation where values are computed
[in]funcpointer to cs_dof_func_t function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_bc_by_value()

cs_xdef_t* cs_equation_add_bc_by_value ( cs_equation_param_t eqp,
const cs_param_bc_type_t  bc_type,
const char *  z_name,
cs_real_t values 
)

Define and initialize a new structure to set a boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]bc_typetype of boundary condition to add
[in]z_namename of the related boundary zone
[in]valuespointer to a array storing the values
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_curlcurl()

void cs_equation_add_curlcurl ( cs_equation_param_t eqp,
cs_property_t property,
int  inversion 
)

Associate a new term related to the curl-curl operator for the equation associated to the given cs_equation_param_t structure.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_t structure
[in]inversion1: true, false otherwise
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_t structure
[in]inversion> 0: true, false otherwise

◆ cs_equation_add_diffusion()

void cs_equation_add_diffusion ( cs_equation_param_t eqp,
cs_property_t property 
)

Associate a new term related to the Laplacian operator for the equation associated to the given cs_equation_param_t structure Laplacian means div-grad (either for vector-valued or scalar-valued equations)

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_t structure

◆ cs_equation_add_graddiv()

void cs_equation_add_graddiv ( cs_equation_param_t eqp,
cs_property_t property 
)

Associate a new term related to the grad-div operator for the equation associated to the given cs_equation_param_t structure.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_t structure

◆ cs_equation_add_ic_by_analytic()

cs_xdef_t* cs_equation_add_ic_by_analytic ( cs_equation_param_t eqp,
const char *  z_name,
cs_analytic_func_t analytic,
void *  input 
)

Define the initial condition for the unknown related to this equation. This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set according to an analytical function.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]analyticpointer to an analytic function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_ic_by_qov()

cs_xdef_t* cs_equation_add_ic_by_qov ( cs_equation_param_t eqp,
const char *  z_name,
double  quantity 
)

Define the initial condition for the unknown related to this equation This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the value related to all the entities belonging to the given mesh location is such that the integral over these cells returns the requested quantity.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" all cells are considered)
[in]quantityquantity to distribute over the mesh location
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_ic_by_value()

cs_xdef_t* cs_equation_add_ic_by_value ( cs_equation_param_t eqp,
const char *  z_name,
cs_real_t val 
)

Define the initial condition for the unknown related to this equation This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here a constant value is set to all the entities belonging to the given mesh location.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" all cells are considered)
[in]valpointer to the value
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_reaction()

int cs_equation_add_reaction ( cs_equation_param_t eqp,
cs_property_t property 
)

Associate a new term related to the reaction operator for the equation associated to the given cs_equation_param_t structure.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_t structure
Returns
the id related to the reaction term

◆ cs_equation_add_sliding_condition()

void cs_equation_add_sliding_condition ( cs_equation_param_t eqp,
const char *  z_name 
)

Define and initialize a new structure to set a sliding boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the related boundary zone

◆ cs_equation_add_source_term_by_analytic()

cs_xdef_t* cs_equation_add_source_term_by_analytic ( cs_equation_param_t eqp,
const char *  z_name,
cs_analytic_func_t func,
void *  input 
)

Add a new source term by initializing a cs_xdef_t structure. Case of a definition by an analytical function.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]funcpointer to an analytical function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_source_term_by_array()

cs_xdef_t* cs_equation_add_source_term_by_array ( cs_equation_param_t eqp,
const char *  z_name,
cs_flag_t  loc,
cs_real_t array,
bool  is_owner,
cs_lnum_t index 
)

Add a new source term by initializing a cs_xdef_t structure. Case of a definition by an array.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]locinformation to know where are located values
[in]arraypointer to an array
[in]is_ownertransfer the lifecycle to the cs_xdef_t structure (true or false)
[in]indexoptional pointer to the array index
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_source_term_by_dof_func()

cs_xdef_t* cs_equation_add_source_term_by_dof_func ( cs_equation_param_t eqp,
const char *  z_name,
cs_flag_t  loc_flag,
cs_dof_func_t func,
void *  input 
)

Add a new source term by initializing a cs_xdef_t structure. Case of a definition by a DoF function.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]loc_flaglocation of element ids given as parameter
[in]funcpointer to a DoF function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_source_term_by_val()

cs_xdef_t* cs_equation_add_source_term_by_val ( cs_equation_param_t eqp,
const char *  z_name,
cs_real_t val 
)

Add a new source term by initializing a cs_xdef_t structure. Case of a definition by a constant value.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" all cells are considered)
[in]valpointer to the value
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_time()

void cs_equation_add_time ( cs_equation_param_t eqp,
cs_property_t property 
)

Associate a new term related to the time derivative operator for the equation associated to the given cs_equation_param_t structure.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_t structure

◆ cs_equation_add_volume_mass_injection_by_analytic()

cs_xdef_t* cs_equation_add_volume_mass_injection_by_analytic ( cs_equation_param_t eqp,
const char *  z_name,
cs_analytic_func_t func,
void *  input 
)

Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using an analytical function.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]funcpointer to an analytical function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_volume_mass_injection_by_qov()

cs_xdef_t* cs_equation_add_volume_mass_injection_by_qov ( cs_equation_param_t eqp,
const char *  z_name,
double *  quantity 
)

Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using a constant quantity distributed over the associated zone's volume.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]quantitypointer to quantity to distribute over the zone
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_volume_mass_injection_by_value()

cs_xdef_t* cs_equation_add_volume_mass_injection_by_value ( cs_equation_param_t eqp,
const char *  z_name,
double *  val 
)

Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using a constant value.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]valpointer to the value
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_xdef_bc()

void cs_equation_add_xdef_bc ( cs_equation_param_t eqp,
cs_xdef_t xdef 
)

Set a boundary condition from an existing cs_xdef_t structure The lifecycle of the cs_xdef_t structure is now managed by the current cs_equation_param_t structure.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]xdefpointer to the cs_xdef_t structure to transfer

◆ cs_equation_clear_param()

void cs_equation_clear_param ( cs_equation_param_t eqp)

Free the contents of a cs_equation_param_t.

The cs_equation_param_t structure itself is not freed, but all its sub-structures are freed.

This is useful for equation parameters which are accessed through field keywords.

Parameters
[in,out]eqppointer to a cs_equation_param_t

◆ cs_equation_copy_param_from()

void cs_equation_copy_param_from ( const cs_equation_param_t ref,
cs_equation_param_t dst,
bool  copy_fid 
)

Copy the settings from one cs_equation_param_t structure to another one. The name is not copied.

Parameters
[in]refpointer to the reference cs_equation_param_t
[in,out]dstpointer to the cs_equation_param_t to update
[in]copy_fidcopy also the field id or not

◆ cs_equation_create_param()

cs_equation_param_t* cs_equation_create_param ( const char *  name,
cs_equation_type_t  type,
int  dim,
cs_param_bc_type_t  default_bc 
)

Create a cs_equation_param_t structure.

Parameters
[in]namename of the equation
[in]typetype of equation
[in]dimdim of the variable associated to this equation
[in]default_bctype of boundary condition set by default
Returns
a pointer to a new allocated cs_equation_param_t structure

Create a cs_equation_param_t structure.

Parameters
[in]namename of the equation
[in]typetype of equation
[in]dimdim of the variable associated to this equation
[in]default_bctype of boundary condition set by default
Returns
a pointer to a new allocated cs_equation_param_t structure

◆ cs_equation_enforce_value_on_cell_selection()

void cs_equation_enforce_value_on_cell_selection ( cs_equation_param_t eqp,
cs_lnum_t  n_elts,
const cs_lnum_t  elt_ids[],
const cs_real_t  ref_value[],
const cs_real_t  elt_values[] 
)

Add an enforcement of the value related to the degrees of freedom associated to the list of selected cells.

One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]n_eltsnumber of selected cells
[in]elt_idslist of cell ids
[in]ref_valueignored if NULL
[in]elt_valueslist of associated values, ignored if NULL

◆ cs_equation_enforce_vertex_dofs()

void cs_equation_enforce_vertex_dofs ( cs_equation_param_t eqp,
cs_lnum_t  n_elts,
const cs_lnum_t  elt_ids[],
const cs_real_t  ref_value[],
const cs_real_t  elt_values[] 
)

Add an enforcement of the value of degrees of freedom located at mesh vertices. The spatial discretization scheme for the given equation has to be CDO-Vertex based or CDO-Vertex+Cell-based schemes.

One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]n_eltsnumber of vertices to enforce
[in]elt_idslist of vertices
[in]ref_valueignored if NULL
[in]elt_valueslist of associated values, ignored if NULL

◆ cs_equation_find_bc()

cs_xdef_t* cs_equation_find_bc ( cs_equation_param_t eqp,
const char *  z_name 
)

Return pointer to existing boundary condition definition structure for the given equation param structure and zone.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
Returns
a pointer to the cs_xdef_t structure if present, or NULL

◆ cs_equation_free_param()

cs_equation_param_t* cs_equation_free_param ( cs_equation_param_t eqp)

Free a cs_equation_param_t.

Parameters
[in]eqppointer to a cs_equation_param_t
Returns
a NULL pointer
Parameters
[in,out]eqppointer to a cs_equation_param_t
Returns
a NULL pointer

◆ cs_equation_param_has_robin_bc()

bool cs_equation_param_has_robin_bc ( const cs_equation_param_t eqp)

Ask if the parameter settings of the equation has requested the treatment of Robin boundary conditions.

Parameters
[in]eqppointer to a cs_equation_param_t
Returns
true or false

◆ cs_equation_param_last_stage()

void cs_equation_param_last_stage ( cs_equation_param_t eqp)

Last modification of the cs_equation_param_t structure before launching the computation.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure

◆ cs_equation_param_log()

void cs_equation_param_log ( const cs_equation_param_t eqp)

Print the detail of a cs_equation_param_t structure.

Parameters
[in]eqppointer to a cs_equation_param_t structure

◆ cs_equation_param_set_sles()

void cs_equation_param_set_sles ( cs_equation_param_t eqp)

Set parameters for initializing SLES structures used for the resolution of the linear system. Settings are related to this equation.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure

◆ cs_equation_remove_bc()

void cs_equation_remove_bc ( cs_equation_param_t eqp,
const char *  z_name 
)

Remove boundary condition from the given equation param structure for a given zone.

If no matching boundary condition is found, the function returns silently.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)

◆ cs_equation_set_param()

void cs_equation_set_param ( cs_equation_param_t eqp,
cs_equation_key_t  key,
const char *  keyval 
)

Set a parameter attached to a keyname in a cs_equation_param_t structure.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]keykey related to the member of eq to set
[in]keyvalaccessor to the value to set