My Project
programmer's documentation
Data Structures | Enumerations | Functions
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.h"
#include "cs_param_cdo.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. More...
 
#define CS_EQUATION_REACTION   (1 << 4) /* 16 */
 Reaction term is needed. More...
 
#define CS_EQUATION_FORCE_VALUES   (1 << 5) /* 32 */
 Add an algebraic manipulation to set the value of a given set of interior degrees of freedom. 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...
 

Enumerations

enum  cs_equation_type_t {
  CS_EQUATION_TYPE_USER, CS_EQUATION_TYPE_GROUNDWATER, CS_EQUATION_TYPE_NAVSTO, CS_EQUATION_TYPE_PREDEFINED,
  CS_EQUATION_N_TYPES
}
 Type of equations managed by the solver. More...
 
enum  cs_equation_key_t {
  CS_EQKEY_ADV_FORMULATION, CS_EQKEY_ADV_SCHEME, 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

static void cs_equation_param_set_flag (cs_equation_param_t *eqp, cs_flag_t flag)
 Update the flag related to a cs_equation_param_t structure. More...
 
static bool cs_equation_param_has_diffusion (const cs_equation_param_t *eqp)
 Ask if the parameters of the equation needs a diffusion term. More...
 
static bool cs_equation_param_has_convection (const cs_equation_param_t *eqp)
 Ask if the parameters of the equation needs a convection term. More...
 
static bool cs_equation_param_has_reaction (const cs_equation_param_t *eqp)
 Ask if the parameters of the equation needs a reaction term. More...
 
static bool cs_equation_param_has_time (const cs_equation_param_t *eqp)
 Ask if the parameters of the equation needs an unsteady term. More...
 
static bool cs_equation_param_has_sourceterm (const cs_equation_param_t *eqp)
 Ask if the parameters of the equation needs a source term. More...
 
static bool cs_equation_param_has_internal_enforcement (const cs_equation_param_t *eqp)
 Ask if the parameters of the equation has an internal enforcement of the degrees of freedom. More...
 
static bool cs_equation_param_has_name (cs_equation_param_t *eqp, const char *name)
 Check if a cs_equation_param_t structure has its name member equal to the given name. More...
 
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_param_update_from (const cs_equation_param_t *ref, cs_equation_param_t *dst)
 Copy the settings from one cs_equation_param_t structure to another one. 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, int field_id)
 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_summary_param (const cs_equation_param_t *eqp)
 Summary of a cs_equation_param_t structure. 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...
 
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)
 Define and initialize a new structure to store parameters related to a diffusion term. More...
 
void cs_equation_add_time (cs_equation_param_t *eqp, cs_property_t *property)
 Define and initialize a new structure to store parameters related to an unsteady term. More...
 
void cs_equation_add_advection (cs_equation_param_t *eqp, cs_adv_field_t *adv_field)
 Define and initialize a new structure to store parameters related to an advection term. More...
 
int cs_equation_add_reaction (cs_equation_param_t *eqp, cs_property_t *property)
 Define and initialize a new structure to store parameters related to a reaction term. More...
 
cs_xdef_tcs_equation_add_source_term_by_val (cs_equation_param_t *eqp, const char *z_name, cs_real_t *val)
 Define a new source term structure and initialize it by value. More...
 
cs_xdef_tcs_equation_add_source_term_by_analytic (cs_equation_param_t *eqp, const char *z_name, cs_analytic_func_t *ana, void *input)
 Define a new source term structure and initialize it by an analytical 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)
 Define a new source term defined by an array. 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 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. We assume that values are interlaced (if eqp->dim > 1) 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_DIFFUSION

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

Diffusion term is needed.

◆ CS_EQUATION_FORCE_VALUES

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

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

◆ 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 << 4) /* 16 */

Reaction term is needed.

◆ CS_EQUATION_UNSTEADY

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

Unsteady term is needed.

Enumeration Type Documentation

◆ cs_equation_key_t

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

Enumerator
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_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 preconditionner from the PETSc library
  • "v_cycle" --> Code_Saturne's in house multigrid with a V-cycle strategy
  • "k_cycle" --> Code_Saturne's in house 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 algorithm
  • "gauss_seidel" --> Gauss-Seidel algorithm
  • "sym_gauss_seidel" --> Symmetric version of Gauss-Seidel algorithm; one backward and forward sweep
  • "cg" --> (default) the standard 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" --> a robust iterative solver but slower as previous one if the system is not difficult to solve
  • "amg" --> an algebraic multigrid iterative solver. Good choice for a symmetric positive definite system.
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" "vol_tot" "weighted_rhs" (default) "matrix_diag"

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 preconditionner associated to an iterative solver. Available choices are:

  • "jacobi" --> 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" for CDO vertex-based scheme
  • "cdo_vcb" for CDO vertex+cell-based scheme
  • "cdo_fb" for CDO face-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": first-order in time (inconditionnally stable)
  • "euler_explicit":
  • "crank_nicolson": second_order in time
  • "theta_scheme": generic time scheme. One recovers "euler_implicit" with theta equal to "1", "explicit" with "0", "crank_nicolson" with "0.5"
CS_EQKEY_TIME_THETA 

Set the value of theta. Only useful if CS_EQKEY_TIME_SCHEME is set to "theta_scheme"

  • Example: "0.75" (keyval must be between 0 and 1)
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_USER 

User-defined equation

CS_EQUATION_TYPE_GROUNDWATER 

Equation related to the groundwater flow 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_N_TYPES 

Function Documentation

◆ cs_equation_add_advection()

void cs_equation_add_advection ( cs_equation_param_t eqp,
cs_adv_field_t adv_field 
)

Define and initialize a new structure to store parameters related to an advection term.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]adv_fieldpointer to a cs_adv_field_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_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_diffusion()

void cs_equation_add_diffusion ( cs_equation_param_t eqp,
cs_property_t property 
)

Define and initialize a new structure to store parameters related to a diffusion term.

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 
)

Define and initialize a new structure to store parameters related to a reaction term.

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 ana,
void *  input 
)

Define a new source term structure and initialize it 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]anapointer to an analytical function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_source_term_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 
)

Define a new source term defined 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_source_term_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 
)

Define a new source term structure and initialize it by 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 
)

Define and initialize a new structure to store parameters related to an unsteady term.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_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_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_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  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. We assume that values are interlaced (if eqp->dim > 1)

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

◆ 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_convection()

static bool cs_equation_param_has_convection ( const cs_equation_param_t eqp)
inlinestatic

Ask if the parameters of the equation needs a convection term.

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

◆ cs_equation_param_has_diffusion()

static bool cs_equation_param_has_diffusion ( const cs_equation_param_t eqp)
inlinestatic

Ask if the parameters of the equation needs a diffusion term.

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

◆ cs_equation_param_has_internal_enforcement()

static bool cs_equation_param_has_internal_enforcement ( const cs_equation_param_t eqp)
inlinestatic

Ask if the parameters of the equation has an internal enforcement of the degrees of freedom.

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

◆ cs_equation_param_has_name()

static bool cs_equation_param_has_name ( cs_equation_param_t eqp,
const char *  name 
)
inlinestatic

Check if a cs_equation_param_t structure has its name member equal to the given name.

Parameters
[in]eqppointer to a cs_equation_param_t structure
[in]namename of the equation
Returns
true if the given eqp has the same name the one given as parameter otherwise false

◆ cs_equation_param_has_reaction()

static bool cs_equation_param_has_reaction ( const cs_equation_param_t eqp)
inlinestatic

Ask if the parameters of the equation needs a reaction term.

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

◆ cs_equation_param_has_sourceterm()

static bool cs_equation_param_has_sourceterm ( const cs_equation_param_t eqp)
inlinestatic

Ask if the parameters of the equation needs a source term.

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

◆ cs_equation_param_has_time()

static bool cs_equation_param_has_time ( const cs_equation_param_t eqp)
inlinestatic

Ask if the parameters of the equation needs an unsteady term.

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_set_flag()

static void cs_equation_param_set_flag ( cs_equation_param_t eqp,
cs_flag_t  flag 
)
inlinestatic

Update the flag related to a cs_equation_param_t structure.

Parameters
[in,out]eqppointer to a cs_equation_param_t
[in]flagflag to add

◆ cs_equation_param_set_sles()

void cs_equation_param_set_sles ( cs_equation_param_t eqp,
int  field_id 
)

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

Parameters
[in]eqppointer to a cs_equation_param_t struct.
[in]field_idid of the cs_field_t struct. for this equation

◆ cs_equation_param_update_from()

void cs_equation_param_update_from ( const cs_equation_param_t ref,
cs_equation_param_t dst 
)

Copy the settings from one cs_equation_param_t structure to another one.

Parameters
[in]refpointer to the reference cs_equation_param_t
[in,out]dstpointer to the cs_equation_param_t to update

◆ 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

◆ cs_equation_summary_param()

void cs_equation_summary_param ( const cs_equation_param_t eqp)

Summary of a cs_equation_param_t structure.

Parameters
[in]eqppointer to a cs_equation_param_t structure