My Project
programmer's documentation
Input of calculation parameters (C functions in cs_user_parameters.c)

Introduction

C user functions for definition of model options and calculation parameters. These subroutines are called in all cases.

If the Code_Saturne GUI is used, this file is not required (but may be used to override parameters entered through the GUI, and to set parameters not accessible through the GUI).

Several functions are present in the file, each destined to defined specific parameters.

Base model related options

Definition of user variables or properties should be defined here, if not already done throught the GUI.

General options (\ref cs_user_parameters)

Most definitions should be done in the cs_user_parameters function (Fortran equivalent is usipsu subroutine).

Choose a turbulent model among the available models

/* Example: Chose a turbulence model
* CS_TURB_NONE: no turbulence model (laminar flow)
* CS_TURB_MIXING_LENGTH: mixing length model
* CS_TURB_K_EPSILON: standard k-epsilon model
* CS_TURB_K_EPSILON_LIN_PROD: k-epsilon model with
* Linear Production (LP) correction
* CS_TURB_K_EPSILON_LS: Launder-Sharma low Re k-epsilon model
* CS_TURB_RIJ_EPSILON_LRR: Rij-epsilon (LRR)
* CS_TURB_RIJ_EPSILON_SSG: Rij-epsilon (SSG)
* CS_TURB_RIJ_EPSILON_EBRSM: Rij-epsilon (EBRSM)
* CS_TURB_LES_SMAGO_CONST: LES (constant Smagorinsky model)
* CS_TURB_LES_SMAGO_DYN: LES ("classical" dynamic Smagorisky model)
* CS_TURB_LES_WALE: LES (WALE)
* CS_TURB_V2F_PHI: v2f phi-model
* CS_TURB_V2F_BL_V2K: v2f BL-v2-k
* CS_TURB_K_OMEGA: k-omega SST
* CS_TURB_SPALART_ALLMARAS: Spalart-Allmaras model */

Coupled solver for Rij components (when iturb=30, 31 or 32): 0 to switch off, 1 to switch on

/* Example: Coupled solver for Rij components (when iturb=30, 31 or 32)
* 0: switch off
* 1: switch on (default)
*/
rans_model->irijco = 1;

To set a thermal model (0: none, 1: temperature, 2: entyhalpy, 3: total energy)

/* Example: Choose a thermal model
*
* 0: none
* 1: temperature
* 2: enthalpy
* 3: total energy (only for compressible module)
*
* For temperature, the temperature scale may be set later using itpscl
* (1 for Kelvin, 2 for Celsius).
*
* Warning: When using specific physics, this value is
* set automatically by the physics model.
*
*/
thermal_model->itherm = 1;

Volume of Fluid model with mass transfer Merkle model to take into account vaporization / condensation

To activate ALE (Arbitrary Lagrangian Eulerian) method (CS_ALE_NONE: switch off, CS_ALE_LEGACY: legacy solver, CS_ALE_CDO: CDO solver)

/* Example: activate ALE (Arbitrary Lagrangian Eulerian) method
* CS_ALE_NONE: switch off
* CS_ALE_LEGACY: legacy solver
* CS_ALE_CDO: CDO solver
*/

The user can add a scalar to be solved

/* Example: add 2 scalar variables ("species" in the GUI nomenclature).
*
* Note that at this (very early) stage of the data setup, fields are
* not defined yet. Associated fields will be defined later (after
* model-defined fields) in the same order as that used here, and
* after user-defined variables defined throught the GUI, if used.
*
* Currently, only 1d (scalar) fields are handled.
*
* parameters for cs_parameters_add_variable():
* name <-- name of variable and associated field
* dim <-- variable dimension
*/
cs_parameters_add_variable("species_1", 1);

After adding a scalar, the user can add the variance of this scalar

/* Example: add the variance of a user variable.
*
* parameters for cs_parameters_add_variable_variance():
* name <-- name of variance and associated field
* variable_name <-- name of associated variable
*/
"species_1");

Add a user property defined on a mesh location (cells, interior faces, boundary faces or vertices).

/* Example: add a user property defined on boundary faces.
*
* parameters for cs_parameters_add_property():
* name <-- name of property and associated field
* dim <-- property dimension
* location_id <-- id of associated mesh location, which must be one of:
* CS_MESH_LOCATION_CELLS
* CS_MESH_LOCATION_INTERIOR_FACES
* CS_MESH_LOCATION_BOUNDARY_FACES
* CS_MESH_LOCATION_VERTICES
*/
cs_parameters_add_property("user_b_property_1",
1,
}

Choose a time step option

/* Time stepping (0 : uniform and constant
1 : variable in time, uniform in space
2 : variable in time and space
-1 : steady algorithm) */
time_opt->idtvar = 0;

Choose a reference time step

/* Reference time step dt_ref
The example given below is probably not adapted to your case. */
cs_real_t dt_ref = 0.005; // was 0.05
domain->time_step->dt_ref = dt_ref;

To set a duration

/* Duration
nt_max is absolute number of the last time step required
if we have already run 10 time steps and want to
run 10 more, nt_max must be set to 10 + 10 = 20 */
domain->time_step->nt_max = (int) (40./dt_ref);

For example, to change the log (run_solver.log) verbosity of all the variables:

{
int n_fields = cs_field_n_fields();
for (int f_id = 0; f_id < n_fields; f_id++) {
if (f->type & CS_FIELD_VARIABLE) {
int key_cal_opt_id = cs_field_key_id("var_cal_opt");
cs_field_get_key_struct(f, key_cal_opt_id, &vcopt);
vcopt.iwarni = 2;
cs_field_set_key_struct(f, key_cal_opt_id, &vcopt);
}
}
}

For example, to change limiters for the convective scheme for a given scalar:

{
/* retrieve scalar field by its name */
cs_field_t *sca1 = cs_field_by_name("scalar1");
/* ischcv is the type of convective scheme:
0: second order linear upwind
1: centered
2: pure upwind gradient in SOLU */
/* isstpc:
0: swich on the slope test
1: swich off the slope test (default)
2: continuous limiter ensuring boundedness (beta limiter)
3: NVD/TVD Scheme */
int key_cal_opt_id = cs_field_key_id("var_cal_opt");
cs_field_get_key_struct(sca1, key_cal_opt_id, &vcopt);
vcopt.ischcv = 0;
vcopt.isstpc = 3;
cs_field_set_key_struct(sca1, key_cal_opt_id, &vcopt);
/* Min/Max limiter or NVD/TVD limiters
* then "limiter_choice" keyword must be set:
* 0: Gamma
* 1: SMART
* 2: CUBISTA
* 3: SUPERBEE
* 4: MUSCL
* 5: MINMOD
* 6: CLAM
* 7: STOIC
* 8: OSHER
* 9: WASEB
* --- VOF scheme ---
* 10: M-HRIC
* 11: M-CICSAM */
int key_lim_id = cs_field_key_id("limiter_choice");
/* Get the Key for the Sup and Inf for the convective scheme */
int kccmin = cs_field_key_id("min_scalar");
int kccmax = cs_field_key_id("max_scalar");
/* Set the Value for the Sup and Inf of the studied scalar
* for the Gamma beta limiter for the temperature */
}

One can also choose the percentage of upwind blending when using the slope test

{
/* retrieve scalar field by its name */
cs_field_t *sca1 = cs_field_by_name("scalar1");
/* blend_st (can be between 0 and 1):
0: full upwind (default)
1: scheme without upwind */
int key_cal_opt_id = cs_field_key_id("var_cal_opt");
cs_field_get_key_struct(sca1, key_cal_opt_id, &vcopt);
vcopt.blend_st = 0.1;
cs_field_set_key_struct(sca1, key_cal_opt_id, &vcopt);
}

If one wants to declare a scalar as buoyant (i.e. the density depends on this scalar through a given equation of state) and add it in the velocity-pressure optional sub-iterations (PISO), one can set the dedicated keyword:

{
/* retrieve scalar field by its name */
cs_field_t *sca1 = cs_field_by_name("scalar1");
int key_is_buoyant = cs_field_key_id("is_buoyant");
cs_field_set_key_int(sca1, key_is_buoyant, 1);
}

If one wants to activate drift terms on a transported scalar:

{
/* retrieve scalar field by its name */
cs_field_t *sca1 = cs_field_by_name("scalar1");
/* Key id for drift scalar */
int key_drift = cs_field_key_id("drift_scalar_model");
int drift = CS_DRIFT_SCALAR_ON;
if (false)
if (false)
if (false)
if (false)
/* Set the key word "drift_scalar_model" into the field structure */
cs_field_set_key_int(sca1, key_drift, drift);
}

To set options for the solving of the Stokes problem (velocity pressure coupling see cs_stokes_model_t structure):

{
stokes->arak = 0.;
}

To set options for the PISO-like sub iterations (velocity pressure coupling see cs_piso_t structure):

{
piso->nterup = 3;
}

Special fields and activate some automatic post-processings

For example, to force the presence of a boundary temperature field, which may be useful for postprocessing:

To add boundary values for all scalars, the following code can be added:

{
int n_fields = cs_field_n_fields();
for (int f_id = 0; f_id < n_fields; f_id++) {
}
}

Enforce existence of 'tplus' and 'tstar' fields, so that a Nusselt number may be computed using the post_boundary_nusselt subroutine. When postprocessing this quantity is activated, those fields are present, but if we need to compute them in the cs_user_extra_operations user subroutine without postprocessing them, forcing the definition of these fields to save the values computed for the boundary layer is necessary.

You can activate the post-processing of the Q-criterion on the whole domain mesh with:

Save contribution of slope test for variables in special fields. These fields are automatically created, with postprocessing output enabled, if the matching variable is convected, does not use a pure upwind scheme, and has a slope test (the slope_test_upwind_id key value for a given variable's field is automatically set to the matching postprocessing field's id, or -1 if not applicable).

{
int n_fields = cs_field_n_fields();
for (int f_id = 0; f_id < n_fields; f_id++) {
cs_field_set_key_int(f, cs_field_key_id("slope_test_upwind_id"), 0);
}
}

You can activate the post-processing of clipping on turbulent quantities on the whole domain mesh with:

Input-output related examples (usipes)

Basic options

Frequency of log output.

Change a property's label (here for density, first checking if it is variable). A field's name cannot be changed, but its label, used for logging and postprocessing output, may be redefined.

if (CS_F_(cp) != NULL)

Postprocessing output

Activate or deactivate probes output.

Probes for Radiative Transfer (Luminance and radiative density flux vector).

cs_field_t *f = cs_field_by_name_try('luminance');
if (f != NULL)
cs_field_key_id("post_vis"),
f = cs_field_by_name_try('radiative_flux');
if (f != NULL)
cs_field_key_id("post_vis"),

Base model for CDO/HHO schemes

Activation of CDO/HHO schemes

Two modes are available CS_DOMAIN_CDO_MODE_ONLY or CS_DOMAIN_CDO_MODE_WITH_FV

CDO/HHO schemes can be activated within this function as follows:

Domain boundary with CDO/HHO schemes

Several types of domain boundaries can be defined. There are gathered in cs_domain_boundary_type_t The definition of the domain boundaries for CDO/HHO schemes can be specified as follows:

{
cs_boundary_t *bdy = domain->boundaries;
/* Choose a boundary by default */
/* Add a new boundary
* >> cs_domain_add_boundary(boundary, type_of_boundary, zone_name);
*
* zone_name is either a predefined one or user-defined one
* type_of_boundary is one of the following keywords:
* CS_BOUNDARY_WALL,
* CS_BOUNDARY_INLET,
* CS_BOUNDARY_OUTLET,
* CS_BOUNDARY_SYMMETRY
*/
}

Output with CDO/HHO schemes

The management of the level and frequency of details written by the code can be specified for CDO/HHO schemes as follows:

{
-1, // restart frequency
10, // output log frequency
2); // verbosity (-1: no, 0, ...)
}

Time step with CDO/HHO schemes

The management of the time step with CDO/HHO schemes can be specified as follows:

{
/*
If there is an inconsistency between the max. number of iteration in
time and the final physical time, the first condition encountered stops
the calculation.
*/
100, // nt_max or -1 (automatic)
-1.); // t_max or < 0. (automatic)
/* Define the value of the time step
>> cs_domain_def_time_step_by_value(domain, dt_val);
>> cs_domain_def_time_step_by_func(domain, dt_func);
The second way to define the time step enable complex definitions.
*/
}

Wall distance with CDO/HHO schemes

The computation of the wall distance with CDO schemes is performed as follows:

Add a user-defined equation with CDO/HHO schemes

The add of a user-defined equation solved by CDO/HHO schemes is specified as follows:

{
/* Add a new user equation:
Set the default boundary condition among:
CS_PARAM_BC_HMG_DIRICHLET or
CS_PARAM_BC_HMG_NEUMANN
By default, initial values are set to zero (or the value given by the
restart file in case of restart).
*/
cs_equation_add_user("AdvDiff.Upw", // equation name
"Pot.Upw", // associated variable field name
1, // dimension of the unknown
CS_PARAM_BC_HMG_DIRICHLET); // default boundary
cs_equation_add_user("AdvDiff.SG", // equation name
"Pot.SG", // associated variable field name
1, // dimension of the unknown
CS_PARAM_BC_HMG_DIRICHLET); // default boundary
}

Add user-defined properties with CDO/HHO schemes

The add of a new user-defined property with CDO/HHO schemes is specified as follows:

{
cs_property_add("conductivity", /* property name */
CS_PROPERTY_ANISO); /* type of material property */
cs_property_add("rho.cp", /* property name */
CS_PROPERTY_ISO); /* type of material property */
}

If you want to compute the Fourier number related to a given property in an unsteady simulation

{
/* Retrieve a property named "conductivity" */
cs_property_t *pty = cs_property_by_name("conductivity");
/* Activate the computation of the Fourier number for this property */
}

Add user-defined advection field with CDO/HHO schemes

The definition of an advection field allows one to handle flows with a frozen velocity field or the transport of scalar quantities without solving the Navier-Stokes system. The add of a new user-defined advection field with CDO/HHO schemes is specified as follows:

{
/* Add a user-defined advection field named "adv_field" */
}

If you need to activate options related to advection fields, you can also specify

{
/* Retrieve an advection field named "adv_field" */
/* Compute the Courant number (if unsteady simulation) */
/* Set other advanced options: for instance, define an interpolation of the
advection field at vertices */
/* Both options in one call */
}

Settings related to the groundwater flow module with CDO/HHO schemes

Activation of the groundwater flow module. The second argument is either 0 if no option is needed or a list of flags among CS_GWF_GRAVITATION, CS_GWF_RICHARDS_UNSTEADY, CS_GWF_SOIL_PROPERTY_UNSTEADY, CS_GWF_SOIL_ALL_SATURATED

{
cs_gwf_activate(CS_PROPERTY_ISO, 0); // no flag to set
}

Add soils (settings of the soil is performed in cs_user_gwf_setup). Soils have to be added before adding tracers.

Add tracer equations which correspond to a transport equation using the darcean flux as the advection field. This call implies the creation of a new equation related to a new variable field name given as arguments. Advanced tracer equations can be defined using cs_gwf_add_tracer_user.

{
cs_gwf_add_tracer("Tracer_01","C1");
}

General options (\ref cs_user_parameters)

CDO/HHO schemes

Modifiy the numerical parameters related to a given equation.

{
/* The modification of the space discretization should be apply first */
/* Modify other parameters than the space discretization */
/* Linear algebra settings */
#if defined(HAVE_PETSC)
#else
#endif
}

Finalize the set-up for CDO/HHO schemes

Set up properties with CDO/HHO schemes

When a property has been added, the second step is to define this property. According to the type of property (isotropic, orthotropic or anisotropic) definitions differ. Here are two examples:

{
cs_property_t *conductivity = cs_property_by_name("conductivity");
cs_real_33_t tensor = {{1.0, 0.5, 0.0}, {0.5, 1.0, 0.5}, {0.0, 0.5, 1.0}};
cs_property_def_aniso_by_value(conductivity, // property structure
"cells", // name of the volume zone
tensor); // values of the property
cs_property_t *rhocp = cs_property_by_name("rho.cp");
cs_real_t iso_val = 2.0;
cs_property_def_iso_by_value(rhocp, // property structure
"cells", // name of the volume zone
iso_val); // value of the property
}

Set up user-defined advection field with CDO/HHO schemes

When an advection field has been added, the second step is to define this advection field. Here are is an example of definition using an anlytic function and the activation of optional features:

{
cs_advection_field_def_by_analytic(adv, _define_adv_field, NULL);
/* Enable also the defintion of the advection field at mesh vertices */
/* Activate the post-processing of the related Courant number */
}

Set up the boundary conditions for an equation

{
"boundary_faces", // zone name
_define_bcs, // pointer to the function
NULL); // input structure
}

Add terms to an equation

Add terms like diffusion term, advection term, unsteady term, reaction terms or source terms.

{
cs_property_t *rhocp = cs_property_by_name("rho.cp");
cs_property_t *conductivity = cs_property_by_name("conductivity");
/* Activate the unsteady term */
cs_equation_add_time(eqp, rhocp);
/* Activate the diffusion term */
cs_equation_add_diffusion(eqp, conductivity);
/* Activate advection effect */
/* Simple definition with cs_equation_add_source_term_by_val
where the value of the source term is given by m^3
*/
cs_real_t st_val = -0.1;
"cells",
&st_val);
}

Linear solver related options

By default, Code_Saturne will use a multigrid algorithm for pressure and iterative solver for other variables. For a given case, checking the setup file resulting from a first calculation will provide more info.

Available solvers include a variety of iterative linear solvers, described in more detail at cs_sles_it_create, and a multigrid solver, whose definition and settings are described at cs_multigrid_create, cs_multigrid_set_coarsening_options, and cs_multigrid_set_solver_options.

Simple options may be set using the GUI, but for more advanced settings are described in this section. It is also recommended to read the documentation of cs_sles.c (which is a solver definition "container"), cs_sles_it.c (iterative solvers, with available types cs_sles_it_type_t), and cs_multigrid.c (which are actual solver implementations). The API provided is extensible, so it is possible for a user to define other solvers or link to external solver libraries using this system, without requiring any modification to non-user source files.

The examples which follow illustrate mostly simple setting changes which may be useful.

Example: distance to wall

By default, the wall distance (active only with turbulence models which require it) is computed with a preconditionned conjugate gradient. The following example shows how to use a multigrid solver for this quantity (useful especially if computed repeatedly, such as for ALE).

Example: user variable

The following example shows how to set the linear solver for a given user variable field so as to use a BiCGStab solver with polynomial preconditioning of degree 1.

cs_field_t *cvar_user_1 = cs_field_by_name_try("user_1");
if (cvar_user_1 != NULL) {
cs_sles_it_define(cvar_user_1->id,
NULL, /* name passed is NULL if field_id > -1 */
1, /* polynomial precond. degree (default 0) */
10000); /* n_max_iter */
}

Changing the verbosity

By default, a linear solver uses the same verbosity as its matching variable, and is not verbose for non-variable quantities. The verbosity may be specifically set for linear system resolution, as shown in the following example:

{
cs_sles_t *sles_p = cs_sles_find_or_add(CS_F_(p)->id, NULL);
}

Example: error visualization

The following example shows how to activate local error visualization output (here for velocity and pressure).

Example: advanced multigrid settings

The following example shows how to set advanced settings for the multigrid solver used for the pressure solution.

{
NULL,
3, /* aggregation_limit (default 3) */
0, /* coarsening_type (default 0) */
10, /* n_max_levels (default 25) */
30, /* min_g_cells (default 30) */
0.95, /* P0P1 relaxation (default 0.95) */
20); /* postprocessing (default 0) */
(mg,
CS_SLES_JACOBI, /* descent smoother type (default: CS_SLES_PCG) */
CS_SLES_JACOBI, /* ascent smoother type (default: CS_SLES_PCG) */
CS_SLES_PCG, /* coarse solver type (default: CS_SLES_PCG) */
50, /* n max cycles (default 100) */
5, /* n max iter for descent (default 2) */
5, /* n max iter for asscent (default 10) */
1000, /* n max iter coarse solver (default 10000) */
0, /* polynomial precond. degree descent (default 0) */
0, /* polynomial precond. degree ascent (default 0) */
1, /* polynomial precond. degree coarse (default 0) */
-1.0, /* precision multiplier descent (< 0 forces max iters) */
-1.0, /* precision multiplier ascent (< 0 forces max iters) */
0.1); /* requested precision multiplier coarse (default 1) */
}

Example: multigrid preconditioner settings

The following example shows how to use multigrid as a preconditioner for a conjugate gradient solver (for the pressure correction), and set advanced settings for that multigrid preconditioner.

{
NULL,
-1,
10000);
assert(strcmp(cs_sles_pc_get_type(cs_sles_it_get_pc(c)), "multigrid") == 0);
(mg,
CS_SLES_P_GAUSS_SEIDEL, /* descent smoother (CS_SLES_P_SYM_GAUSS_SEIDEL) */
CS_SLES_P_GAUSS_SEIDEL, /* ascent smoother (CS_SLES_P_SYM_GAUSS_SEIDEL) */
CS_SLES_PCG, /* coarse solver (CS_SLES_P_GAUSS_SEIDEL) */
1, /* n max cycles (default 1) */
1, /* n max iter for descent (default 1) */
1, /* n max iter for asscent (default 1) */
500, /* n max iter coarse solver (default 1) */
0, /* polynomial precond. degree descent (default) */
0, /* polynomial precond. degree ascent (default) */
0, /* polynomial precond. degree coarse (default 0) */
-1.0, /* precision multiplier descent (< 0 forces max iters) */
-1.0, /* precision multiplier ascent (< 0 forces max iters) */
1.0); /* requested precision multiplier coarse (default 1) */
}

Multigrid parallel settings

In parallel, grids may optionally be merged across neigboring ranks when their local size becomes too small. This tends to deepen the grid hierarchy, as some parallel rank boundaries are removed. Depending on the architecture and processor/network performance ratios, this may increase or decrease performance.

cs_grid_set_merge_options(4, /* # of ranks merged at a time */
300, /* mean # of cells under which we merge */
500, /* global # of cells under which we merge */
1); /* # of ranks under which we do not merge */

Example: DOM radiation settings

For DOM radiation models, 1 solver is assigned for each direction this allows using a specific ordering for each direction for the default Block Gauss-Seidel solver.

The example below shows how to set a non-default linear solver for DOM radiation. Here, we assume a quadrature with 32 directions is used (if more solvers than directions are specified, the extra definitions will be unused, but this causes no further issues).

{
for (int i = 0; i < 32; i++) {
char name[16];
sprintf(name, "radiation_%03d", i+1);
name,
0, /* poly_degree */
1000); /* n_max_iter */
}
}

Plotting solver convergence

The following example shows how to activate convergence plotting for built-in iterative or multigrid solvers.

{
const cs_field_t *f = CS_F_(p);
cs_sles_t *sles_p = cs_sles_find_or_add(f->id, NULL);
bool use_iteration = true; /* use iteration or wall clock time for axis */
if (strcmp(cs_sles_get_type(sles_p), "cs_sles_it_t") == 0) {
cs_sles_it_set_plot_options(c, f->name, use_iteration);
}
else if (strcmp(cs_sles_get_type(sles_p), "cs_multigrid_t") == 0) {
cs_multigrid_set_plot_options(c, f->name, use_iteration);
}
}

Plots will appear as CSV (comma-separated value) files in the monitoring subdirectory.

Using PETSc

The following example shows how to setup a solver to use the PETSc library, if the code was built with PETSc support.

General options (those passed to PETSc through command line options) may be defined directly in cs_user_linear_solvers, for example:

{
/* Initialization must be called before setting options;
it does not need to be called before calling
cs_sles_petsc_define(), as this is handled automatically. */
PETSC_COMM_WORLD = cs_glob_mpi_comm;
PetscInitializeNoArguments();
/* See the PETSc documentation for the options database */
#if PETSC_VERSION_GE(3,7,0)
PetscOptionsSetValue(NULL, "-ksp_type", "cg");
PetscOptionsSetValue(NULL, "-pc_type", "jacobi");
#else
PetscOptionsSetValue("-ksp_type", "cg");
PetscOptionsSetValue("-pc_type", "jacobi");
#endif
}

A specific system may be set up to use PETsc, as is shown here for the pressure variable:

{
NULL,
MATSHELL,
_petsc_p_setup_hook,
NULL);
}

The basic matrix format to be used by PETSc is defined at this stage, using a PETSc MatType string (see PETSc documentation). Further options may be defined in a setup hook function, as follows:

static void
_petsc_p_setup_hook(const void *context,
KSP ksp)
{
PC pc;
KSPSetType(ksp, KSPCG); /* Preconditioned Conjugate Gradient */
KSPSetNormType(ksp, KSP_NORM_UNPRECONDITIONED); /* Try to have "true" norm */
KSPGetPC(ksp, &pc);
PCSetType(pc, PCJACOBI); /* Jacobi (diagonal) preconditioning */
}

If no additional settings are required, the matching parameter in cs_sles_petsc_define may be set to NULL.

Using PETSc with GAMG

To use PETSc's GAMG (geometric-algebraic multigrid) preconditioner, the following general options should be set:

{
/* Initialization must be called before setting options;
it does not need to be called before calling
cs_sles_petsc_define(), as this is handled automatically. */
PETSC_COMM_WORLD = cs_glob_mpi_comm;
PetscInitializeNoArguments();
/* See the PETSc documentation for the options database */
#if PETSC_VERSION_GE(3,7,0)
PetscOptionsSetValue(NULL, "-ksp_type", "cg");
PetscOptionsSetValue(NULL, "-pc_type", "gamg");
PetscOptionsSetValue(NULL, "-pc_gamg_agg_nsmooths", "1");
PetscOptionsSetValue(NULL, "-mg_levels_ksp_type", "richardson");
PetscOptionsSetValue(NULL, "-mg_levels_pc_type", "sor");
PetscOptionsSetValue(NULL, "-mg_levels_ksp_max_it", "1");
PetscOptionsSetValue(NULL, "-pc_gamg_threshold", "0.02");
PetscOptionsSetValue(NULL, "-pc_gamg_reuse_interpolation", "TRUE");
PetscOptionsSetValue(NULL, "-pc_gamg_square_graph", "4");
#else
PetscOptionsSetValue("-ksp_type", "cg");
PetscOptionsSetValue("-pc_type", "gamg");
PetscOptionsSetValue("-pc_gamg_agg_nsmooths", "1");
PetscOptionsSetValue("-mg_levels_ksp_type", "richardson");
PetscOptionsSetValue("-mg_levels_pc_type", "sor");
PetscOptionsSetValue("-mg_levels_ksp_max_it", "1");
PetscOptionsSetValue("-pc_gamg_threshold", "0.02");
PetscOptionsSetValue("-pc_gamg_reuse_interpolation", "TRUE");
PetscOptionsSetValue("-pc_gamg_square_graph", "4");
#endif
}

Setting GAMG-preconditioned PCG for the pressure may be done as in the previous option, ensuring here that a matrix structure native to PETSc is used (SEQAIJ, MPIAIJ, or some other AIJ variant), so all required matrix operations are available:

{
NULL,
MATMPIAIJ,
_petsc_p_setup_hook_gamg,
NULL);
}

With the associated setup hook function:

static void
_petsc_p_setup_hook_gamg(const void *context,
KSP ksp)
{
PC pc;
KSPSetType(ksp, KSPCG); /* Preconditioned Conjugate Gradient */
KSPGetPC(ksp, &pc);
PCSetType(pc, PCGAMG); /* GAMG (geometric-algebraic multigrid)
preconditioning */
}

Using PETSc with HYPRE

To use HYPRE's Boomer AMG as a PETSc preconditioner, the following general options should be set:

{
/* Initialization must be called before setting options;
it does not need to be called before calling
cs_sles_petsc_define(), as this is handled automatically. */
PETSC_COMM_WORLD = cs_glob_mpi_comm;
PetscInitializeNoArguments();
/* See the PETSc documentation for the options database */
#if PETSC_VERSION_GE(3,7,0)
PetscOptionsSetValue(NULL, "-ksp_type", "cg");
PetscOptionsSetValue(NULL, "-pc_type", "hypre");
PetscOptionsSetValue(NULL, "-pc_hypre_type","boomeramg");
PetscOptionsSetValue(NULL, "-pc_hypre_boomeramg_coarsen_type","HMIS");
PetscOptionsSetValue(NULL, "-pc_hypre_boomeramg_interp_type","ext+i-cc");
PetscOptionsSetValue(NULL, "-pc_hypre_boomeramg_agg_nl","2");
PetscOptionsSetValue(NULL, "-pc_hypre_boomeramg_P_max","4");
PetscOptionsSetValue(NULL, "-pc_hypre_boomeramg_strong_threshold","0.5");
PetscOptionsSetValue(NULL, "-pc_hypre_boomeramg_no_CF","");
#else
PetscOptionsSetValue("-ksp_type", "cg");
PetscOptionsSetValue("-pc_type", "hypre");
PetscOptionsSetValue("-pc_hypre_type","boomeramg");
PetscOptionsSetValue("-pc_hypre_boomeramg_coarsen_type","HMIS");
PetscOptionsSetValue("-pc_hypre_boomeramg_interp_type","ext+i-cc");
PetscOptionsSetValue("-pc_hypre_boomeramg_agg_nl","2");
PetscOptionsSetValue("-pc_hypre_boomeramg_P_max","4");
PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold","0.5");
PetscOptionsSetValue("-pc_hypre_boomeramg_no_CF","");
#endif
}

Setting BoomerAMG-preconditioned PCG for the pressure may be done as in the previous option, ensuring here that a matrix structure native to PETSc is used (SEQAIJ, MPIAIJ, or some other AIJ variant), so all required matrix operations are available:

{
NULL,
MATMPIAIJ,
_petsc_p_setup_hook_bamg,
NULL);
}

With the associated setup hook function:

static void
_petsc_p_setup_hook_bamg(const void *context,
KSP ksp)
{
PC pc;
KSPSetType(ksp, KSPCG); /* Preconditioned Conjugate Gradient */
KSPGetPC(ksp, &pc);
PCSetType(pc, PCHYPRE); /* HYPRE BoomerAMG preconditioning */
}

Additional PETSc features

Many additional features are possible with PETSc; for example, the following setup hook also outputs a view of the matrix, depending on an environment variable, CS_USER_PETSC_MAT_VIEW, which may take values DEFAULT, DRAW_WORLD, or DRAW:

static void
_petsc_p_setup_hook_view(const void *context,
KSP ksp)
{
const char *p = getenv("CS_USER_PETSC_MAT_VIEW");
if (p != NULL) {
/* Get system and preconditioner matrixes */
Mat a, pa;
KSPGetOperators(ksp, &a, &pa);
/* Output matrix in several ways depending on
CS_USER_PETSC_MAT_VIEW environment variable */
if (strcmp(p, "DEFAULT") == 0)
MatView(a, PETSC_VIEWER_DEFAULT);
else if (strcmp(p, "DRAW_WORLD") == 0)
MatView(a, PETSC_VIEWER_DRAW_WORLD);
else if (strcmp(p, "DRAW") == 0) {
PetscViewer viewer;
PetscDraw draw;
PetscViewerDrawOpen(PETSC_COMM_WORLD, NULL, "PETSc View",
0, 0, 600, 600, &viewer);
PetscViewerDrawGetDraw(viewer, 0, &draw);
PetscViewerDrawSetPause(viewer, -1);
MatView(a, viewer);
PetscDrawPause(draw);
PetscViewerDestroy(&viewer);
}
}
}

Time moment related options

Code_Saturne allows the calculation of temporal means or variances, either of expressions evaluated through a user function, or of expressions of the type $<f_1*f_2...*f_n>$. The variables may be fields or field components. This is done calling either through the GUI, or in the user function cs_user_time_moments. Each temporal mean is declared using either cs_time_moment_define_by_func, or cs_time_moment_define_by_field_ids.

For each time moment, a starting time step or value is defined. If the starting time step number is negative, the time value is used instead.

The moment values are stored as fields, and updated at each time step, using recursive formulas. Before the matching moment computation starting time step, a moment's values are uniformly 0. For visualization an interpretation reasons, only fields of dimension 1, 3, 6, or 9 (scalars, vectors, or tensors of rank 2) are allowed, so moment definitions not matching this constraint should be split.

To count defined moments, use the cs_time_moment_n_moments function, whether from Fortran or C. To access the matching fields, use time_moment_field_id in Fortran, or cs_time_moment_get_field in C.

Examples

Example 1

In the following example, we define a moment for the mean velocity. All components are used (component -1 means all components), so the moment is a vector.

{
/* Moment <U> calculated starting from time step 1000. */
/* resulting field is a vector */
int moment_f_id[] = {CS_F_(vel)->id};
int moment_c_id[] = {-1};
int n_fields = 1;
n_fields,
moment_f_id,
moment_c_id,
1000, /* nt_start */
-1, /* t_start */
NULL);
}

In the next example, we define the variance of the vector velocity. All components are used again (component -1 means all components), so the moment is a tensor.

{
/* Second order moment <UU>-<U><U> calculated starting from time step 1000. */
/* resulting field is a tensor */
int moment_f_id[] = {CS_F_(vel)->id};
int moment_c_id[] = {-1};
int n_fields = 1;
n_fields,
moment_f_id,
moment_c_id,
1000, /* nt_start */
-1, /* t_start */
NULL);
}

Example 2

In the next example, we multiply the expression by the density. As the density is of dimension 1, and the velocity of dimension 3, the resulting moment is of dimension 3.

{
/* Moment <rho U> (vector) calculated starting from time step 1000. */
/* resulting field is a vector */
int moment_f_id[] = {CS_F_(rho)->id, CS_F_(vel)->id};
int moment_c_id[] = {-1, -1};
int n_fields = 2;
n_fields,
moment_f_id,
moment_c_id,
1000, /* nt_start */
-1, /* t_start */
NULL);
}

Example 3

In the next example, we define a product of several field components, all of dimension 1, as we consider only the x and y components of the velocity; for the density, we cas use either component 0 or -1 (all), since the field is scalar.

This moment's computation is also restarted at each time step.

{
int moment_f_id[] = {CS_F_(rho)->id, CS_F_(vel)->id, CS_F_(vel)->id};
int moment_c_id[] = {-1, 0, 1};
int n_fields = 3;
n_fields,
moment_f_id,
moment_c_id,
-1, /* nt_start */
20.0, /* t_start */
NULL);

Example 4

This next example illustrates the use of user-defined functions to evaluate expressions. Here, we compute the moment of the sum ot two variables (which obviously cannot be expressed as a product), so we first need to define an appropriate function, matching the signature of a cs_time_moment_data_t function. We can name that function as we choose, so naming for clarity is recommmended. Note that in this case, the input argument is not used. This argument may be useful to pass data to the function, or distinguish between calls to a same function.

Note also that we compute both means and variances here.

static void
_simple_data_sum(const void *input,
cs_real_t *vals)
{
const int location_id = CS_MESH_LOCATION_CELLS;
const cs_lnum_t n_elts = cs_mesh_location_get_n_elts(location_id)[0];
const cs_real_t *s1 = cs_field_by_name("species_1")->val;
const cs_real_t *s2 = cs_field_by_name("species_2")->val;
for (cs_lnum_t i = 0; i < n_elts; i++) {
vals[i] = s1[i] + s2[i];
}
}

In cs_user_time_moments, we can now assign that function to a moments definition:

{
const char *sum_comp_name[] = {"species_sum_mean", "species_sum_variance"};
for (int i = 0; i < 2; i++) {
1, /* field dimension */
_simple_data_sum, /* data_func */
NULL, /* data_input */
NULL, /* w_data_func */
NULL, /* w_data_input */
m_type[i],
1000, /* nt_start */
-1, /* t_start */
NULL);
}
}

Example 5

This next example illustrates the use of another user-defined function to evaluate expressions. Here, we compute the moment of the thermal flux at the boundary. We also that we compute both means and variances here.

static void
_boundary_thermal_flux(const void *input,
cs_real_t *vals)
{
const int location_id = CS_MESH_LOCATION_BOUNDARY_FACES;
const cs_lnum_t n_elts = cs_mesh_location_get_n_elts(location_id)[0];
cs_post_boundary_flux(f->name, n_elts, NULL, vals);
}

In cs_user_time_moments, we assign that function to a moments definition:

{
const char *sum_comp_name[] = {"thermal_flux_mean", "thermal_flux_variance"};
for (int i = 0; i < 2; i++) {
1, /* field dimension */
_boundary_thermal_flux, /* data_func */
NULL, /* data_input */
NULL, /* w_data_func */
NULL, /* w_data_input */
m_type[i],
1, /* nt_start */
-1, /* t_start */
NULL);
}
}

Example 6

In this last example, we compute components of the mean velocity in the case of a rotating mesh. As the mesh orientation changes at each time step, it is necessary to compensate for this rotation when computing the mean, relative to a given mesh position. When using the matching moment, it will also be necessary to account for the mesh position.

Here, the same function will be called for each component, so an input array is defined, with a different key (here a simple character) used for each call.

static void
_velocity_moment_data(const void *input,
cs_real_t *vals)
{
const char key = *((const char *)input);
const int location_id = CS_MESH_LOCATION_CELLS;
const cs_lnum_t n_elts = cs_mesh_location_get_n_elts(location_id)[0];
const cs_real_3_t *vel = (const cs_real_3_t *)(CS_F_(vel)->val);
const cs_real_3_t *restrict cell_cen
double omgnrm = fabs(rot->omega);
/* Axial, tangential and radial unit vectors */
cs_real_3_t e_ax = {rot->axis[0], rot->axis[1], rot->axis[2]};
for (cs_lnum_t i = 0; i < n_elts; i++) {
cs_rotation_velocity(rot, cell_cen[i], e_th);
double xnrm = sqrt(cs_math_3_square_norm(e_th));
e_th[0] /= xnrm;
e_th[1] /= xnrm;
e_th[2] /= xnrm;
cs_rotation_coriolis_v(rot, -1., e_th, e_r);
xnrm = sqrt(cs_math_3_square_norm(e_r));
e_r[0] /= xnrm;
e_r[1] /= xnrm;
e_r[2] /= xnrm;
/* Radius */
cs_real_t xr = cs_math_3_dot_product(cell_cen[i], e_r);
/* Axial, tangential and radial components of velocity */
cs_real_t xva = vel[i][0]*e_ax[0] + vel[i][1]*e_ax[1] + vel[i][2]*e_ax[2];
cs_real_t xvt = vel[i][0]*e_th[0] + vel[i][1]*e_th[1] + vel[i][2]*e_th[2];
cs_real_t xvr = vel[i][0]*e_r[0] + vel[i][1]*e_r[1] + vel[i][2]*e_r[2];
/* Entrainment velocity is removed */
xvt -= omgnrm*xr;
/* Store value */
if (key == 'r')
vals[i] = xvr;
else if (key == 't')
vals[i] = xvt;
else if (key == 'a')
vals[i] = xva;
}
}

Note that the input arrays must be accessible when updating moments at each time step, so the array of inputs is declared static in cs_user_time_moments. Fo more complex inputs, we would have an array of inputs here; in this simple case, we could pass a simple call id as the input, casting from point to integer.

{
const char *vel_comp_name[] = {"Wr_moy", "Wt,moy", "Wa_moy"};
/* Data input must be "static" so it can be used in later calls */
static char vel_comp_input[3] = {'r', 't', 'a'};
for (int comp_id = 0; comp_id < 3; comp_id++) {
cs_time_moment_define_by_func(vel_comp_name[comp_id],
1,
_velocity_moment_data, /* data_func */
&(vel_comp_input[comp_id]), /* data_input */
NULL, /* w_data_func */
NULL, /* w_data_input */
74000, /* nt_start */
-1, /* t_start */
NULL);
}
}

To activate means for all variables:

for (int f_id = 0; f_id < cs_field_n_fields(); f_id++) {
if (f->type & CS_FIELD_VARIABLE) {
int moment_f_id[] = {f_id};
int moment_c_id[] = {-1};
int n_fields = 1;
char *extension = "_mean";
char *mean_name;
BFT_MALLOC(mean_name, strlen(f->name) + 1 + 5, char);
strcpy(mean_name, f->name); /* copy field name into the new var */
strcat(mean_name, extension); /* add the extension */
n_fields,
moment_f_id,
moment_c_id,
10, /* nt_start */
-1, /* t_start */
NULL);
}
}

Fan modelling options

Code_Saturne allows modelling of some circular fans as volume regions, defined by simple geometric characteristics, and modeled as explicit momentum source terms in those regions.

Fan pressure characteristic curves are defined as a 2nd order polynomial, and a torque may also be specified. For correct results, it is important that the mesh match the fan dimensions and placement (thickness, hub, blades, and total radius).

The following example shows how a fan may be defined:

{
const cs_real_t inlet_axis_coords[3] = {0, 0, 0};
const cs_real_t outlet_axis_coords[3] = {0.1, 0, 0};
const cs_real_t fan_radius = 0.7;
const cs_real_t blades_radius = 0.5;
const cs_real_t hub_radius = 0.1;
const cs_real_t pressure_curve_coeffs[3] = {0.6, -0.1, -0.05};
const cs_real_t axial_torque = 0.01;
cs_fan_define(3, /* fan (mesh) dimension (2D or 3D) */
inlet_axis_coords,
outlet_axis_coords,
fan_radius,
blades_radius,
hub_radius,
pressure_curve_coeffs,
axial_torque);
}
cs_thermal_model_t::itherm
int itherm
Definition: cs_thermal_model.h:76
cs_field_set_key_int_bits
int cs_field_set_key_int_bits(cs_field_t *f, int key_id, int mask)
Set integer bits matching a mask to 1 for a given key for a field.
Definition: cs_field.c:3045
CS_PARAM_BC_HMG_DIRICHLET
Definition: cs_param.h:304
CS_DRIFT_SCALAR_CENTRIFUGALFORCE
Definition: cs_parameters.h:164
f_id
void const int * f_id
Definition: cs_gui.h:292
cs_get_glob_piso
cs_piso_t * cs_get_glob_piso(void)
Provide acces to cs_glob_piso.
Definition: cs_parameters.c:878
cs_time_moment_type_t
cs_time_moment_type_t
Moment type.
Definition: cs_time_moment.h:56
input
static int input(void)
CS_MULTIGRID_V_CYCLE
Definition: cs_multigrid.h:59
CS_EQKEY_SPACE_SCHEME
Definition: cs_equation_param.h:668
CS_BOUNDARY_INLET
Definition: cs_boundary.h:55
cs_domain_t
Structure storing the main features of the computational domain and pointers to the main geometrical ...
Definition: cs_domain.h:87
CS_PARAM_BC_DIRICHLET
Definition: cs_param.h:305
cs_field_set_key_struct
int cs_field_set_key_struct(cs_field_t *f, int key_id, void *s)
Assign a simple structure for a given key to a field.
Definition: cs_field.c:3336
CS_MESH_LOCATION_BOUNDARY_FACES
Definition: cs_mesh_location.h:65
cs_rotation_velocity
static void cs_rotation_velocity(const cs_rotation_t *r, const cs_real_t coords[3], cs_real_t vr[3])
Compute velocity relative to a fixed frame at a given point.
Definition: cs_rotation.h:127
cs_property_add
cs_property_t * cs_property_add(const char *name, cs_property_type_t type)
Create and initialize a new property structure.
Definition: cs_property.c:293
cs_sles_petsc_define
cs_sles_petsc_t * cs_sles_petsc_define(int f_id, const char *name, MatType matrix_type, cs_sles_petsc_setup_hook_t *setup_hook, void *context)
Define and associate a PETSc linear system solver for a given field or equation name.
Definition: cs_sles_petsc.c:536
CS_EQKEY_PRECOND
Definition: cs_equation_param.h:665
CS_TIME_MOMENT_VARIANCE
Definition: cs_time_moment.h:59
CS_BOUNDARY_OUTLET
Definition: cs_boundary.h:56
cs_xdef_t
Structure storing medata for defining a quantity in a very flexible way.
Definition: cs_xdef.h:126
cs_sles_set_verbosity
void cs_sles_set_verbosity(cs_sles_t *sles, int verbosity)
Set the verbosity for a given linear equation solver.
Definition: cs_sles.c:1345
cs_var_cal_opt_t::ischcv
int ischcv
Definition: cs_parameters.h:70
cs_sles_it_transfer_pc
void cs_sles_it_transfer_pc(cs_sles_it_t *context, cs_sles_pc_t **pc)
Assign a preconditioner to an iterative sparse linear equation solver, transfering its ownership to t...
Definition: cs_sles_it.c:5407
cs_fuel_incl::a
double precision, save a
Definition: cs_fuel_incl.f90:146
cs_glob_post_util_flag
int cs_glob_post_util_flag[]
cs_sles_get_type
const char * cs_sles_get_type(cs_sles_t *sles)
Return type name of solver context.
Definition: cs_sles.c:1443
cs_glob_mesh_quantities
cs_mesh_quantities_t * cs_glob_mesh_quantities
cs_math_3_square_norm
static cs_real_t cs_math_3_square_norm(const cs_real_t v[3])
Compute the square norm of a vector of 3 real values.
Definition: cs_math.h:388
rho
Definition: cs_field_pointer.h:103
cs_sles_it_get_pc
cs_sles_pc_t * cs_sles_it_get_pc(cs_sles_it_t *context)
Return a preconditioner context for an iterative sparse linear equation solver.
Definition: cs_sles_it.c:5378
cs_rotation_coriolis_v
static void cs_rotation_coriolis_v(const cs_rotation_t *r, cs_real_t c, const cs_real_t v[3], cs_real_t vr[3])
Compute a vector Coriolis term.
Definition: cs_rotation.h:173
cs_time_moment_define_by_func
int cs_time_moment_define_by_func(const char *name, int location_id, int dim, cs_time_moment_data_t *data_func, const void *data_input, cs_time_moment_data_t *w_data_func, void *w_data_input, cs_time_moment_type_t type, int nt_start, double t_start, cs_time_moment_restart_t restart_mode, const char *restart_name)
Define a moment whose data values will be computed using a specified function.
Definition: cs_time_moment.c:1577
cs_sles_pc_get_context
void * cs_sles_pc_get_context(cs_sles_pc_t *pc)
Return pointer to preconditioner context structure pointer.
Definition: cs_sles_pc.c:850
cs_equation_param_by_name
cs_equation_param_t * cs_equation_param_by_name(const char *eqname)
Return the cs_equation_param_t structure associated to a cs_equation_t structure thanks to the equati...
Definition: cs_equation.c:414
restrict
#define restrict
Definition: cs_defs.h:127
cs_stokes_model_t::arak
double arak
Definition: cs_stokes_model.h:64
cs_sles_pc_t
struct _cs_sles_pc_t cs_sles_pc_t
Definition: cs_sles_pc.h:66
cs_real_3_t
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:315
cs_property_set_option
void cs_property_set_option(cs_property_t *pty, cs_property_key_t key)
Set optional parameters related to a cs_property_t structure.
Definition: cs_property.c:383
cs_property_def_aniso_by_value
cs_xdef_t * cs_property_def_aniso_by_value(cs_property_t *pty, const char *zname, cs_real_t tens[3][3])
Define an anisotropic cs_property_t structure by value for entities related to a volume zone.
Definition: cs_property.c:643
cs_fan_define
void cs_fan_define(int fan_dim, const cs_real_t inlet_axis_coords[3], const cs_real_t outlet_axis_coords[3], cs_real_t fan_radius, cs_real_t blades_radius, cs_real_t hub_radius, const cs_real_t curve_coeffs[3], cs_real_t axial_torque)
Fan definition (added to the ones previously defined)
Definition: cs_fan.c:226
cs_get_glob_turb_rans_model
cs_turb_rans_model_t * cs_get_glob_turb_rans_model(void)
Provide acces to cs_glob_turb_rans_model.
Definition: cs_turbulence_model.c:1239
cs_field_t::id
int id
Definition: cs_field.h:128
cs_sles_it_define
cs_sles_it_t * cs_sles_it_define(int f_id, const char *name, cs_sles_it_type_t solver_type, int poly_degree, int n_max_iter)
Define and associate an iterative sparse linear system solver for a given field or equation name.
Definition: cs_sles_it.c:4632
cs_boundary_t
Structure storing information related to the "physical" boundaries that one want to set on the comput...
Definition: cs_boundary.h:77
cs_sles_it_t
struct _cs_sles_it_t cs_sles_it_t
Definition: cs_sles_it.h:79
cs_sles_it_set_plot_options
void cs_sles_it_set_plot_options(cs_sles_it_t *context, const char *base_name, bool use_iteration)
Set plotting options for an iterative sparse linear equation solver.
Definition: cs_sles_it.c:5743
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.
Definition: cs_equation_param.c:1266
cs_glob_domain
cs_domain_t * cs_glob_domain
Definition: cs_domain.c:85
cs_field_by_name
cs_field_t * cs_field_by_name(const char *name)
Return a pointer to a field based on its name.
Definition: cs_field.c:2331
cs_domain_t::boundaries
cs_boundary_t * boundaries
Definition: cs_domain.h:101
cs_glob_rotation
cs_rotation_t * cs_glob_rotation
cs_glob_ale
int cs_glob_ale
Definition: cs_ale.c:79
cs_boundary_add
void cs_boundary_add(cs_boundary_t *bdy, cs_boundary_type_t type, const char *zone_name)
Add a new boundary type for a given boundary zone.
Definition: cs_boundary.c:383
cs_multigrid_set_plot_options
void cs_multigrid_set_plot_options(cs_multigrid_t *mg, const char *base_name, bool use_iteration)
Set plotting options for multigrid.
Definition: cs_multigrid.c:4051
cs_vof_parameters_t
VOF model parameters. Void fraction variable tracks fluid 2.
Definition: cs_vof.h:81
rij
Definition: cs_field_pointer.h:79
cs_get_glob_time_step_options
cs_time_step_options_t * cs_get_glob_time_step_options(void)
Provide acces to cs_glob_time_step_options.
Definition: cs_time_step.c:414
cs_field_t::name
const char * name
Definition: cs_field.h:126
CS_ALE_LEGACY
Definition: cs_ale.h:56
cs_property_by_name
cs_property_t * cs_property_by_name(const char *name)
Find the related property definition from its name.
Definition: cs_property.c:336
cs_real_t
double cs_real_t
Floating-point value.
Definition: cs_defs.h:302
CS_EQKEY_SOLVER_FAMILY
Definition: cs_equation_param.h:667
cs_field_t::type
int type
Definition: cs_field.h:129
CS_VOF_ENABLED
#define CS_VOF_ENABLED
Definition: cs_vof.h:60
CS_DRIFT_SCALAR_ON
Definition: cs_parameters.h:159
cs_var_cal_opt_t::isstpc
int isstpc
Definition: cs_parameters.h:72
cs_thermal_model_t
Thermal model descriptor.
Definition: cs_thermal_model.h:74
cs_field_by_name_try
cs_field_t * cs_field_by_name_try(const char *name)
Return a pointer to a field based on its name if present.
Definition: cs_field.c:2357
cs_var_cal_opt_t::blend_st
double blend_st
Definition: cs_parameters.h:87
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.
Definition: cs_equation_param.c:2450
cs_time_moment_define_by_field_ids
int cs_time_moment_define_by_field_ids(const char *name, int n_fields, const int field_id[], const int component_id[], cs_time_moment_type_t type, int nt_start, double t_start, cs_time_moment_restart_t restart_mode, const char *restart_name)
Define a moment of a product of existing fields components.
Definition: cs_time_moment.c:1520
cs_gwf_add_tracer
cs_gwf_tracer_t * cs_gwf_add_tracer(const char *eq_name, const char *var_name)
Add a new equation related to the groundwater flow module This equation is a particular type of unste...
Definition: cs_gwf.c:857
CS_ADVKEY_DEFINE_AT_VERTICES
Definition: cs_advection_field.h:89
CS_EQKEY_HODGE_DIFF_COEF
Definition: cs_equation_param.h:657
cs_sles_pc_get_type
const char * cs_sles_pc_get_type(cs_sles_pc_t *pc)
Return type name of preconditioner context.
Definition: cs_sles_pc.c:804
CS_EQKEY_VERBOSITY
Definition: cs_equation_param.h:671
cs_multigrid_set_solver_options
void cs_multigrid_set_solver_options(cs_multigrid_t *mg, cs_sles_it_type_t descent_smoother_type, cs_sles_it_type_t ascent_smoother_type, cs_sles_it_type_t coarse_solver_type, int n_max_cycles, int n_max_iter_descent, int n_max_iter_ascent, int n_max_iter_coarse, int poly_degree_descent, int poly_degree_ascent, int poly_degree_coarse, double precision_mult_descent, double precision_mult_ascent, double precision_mult_coarse)
Set multigrid parameters for associated iterative solvers.
Definition: cs_multigrid.c:3150
cs_parameters_add_boundary_temperature
cs_field_t * cs_parameters_add_boundary_temperature(void)
Define a boundary values field for temperature, if applicable.
Definition: cs_parameters.c:1451
cs_advection_field_set_option
void cs_advection_field_set_option(cs_adv_field_t *adv, cs_advection_field_key_t key)
Set optional parameters related to a cs_adv_field_t structure.
Definition: cs_advection_field.c:568
cs_gwf_activate
cs_gwf_t * cs_gwf_activate(cs_property_type_t pty_type, cs_flag_t flag)
Initialize the module dedicated to groundwater flows.
Definition: cs_gwf.c:605
cs_get_glob_turb_model
cs_turb_model_t * cs_get_glob_turb_model(void)
Provide access to global turbulence model structure cs_glob_turb_model.
Definition: cs_turbulence_model.c:1167
cs_multigrid_t
struct _cs_multigrid_t cs_multigrid_t
Definition: cs_multigrid.h:67
CS_DOMAIN_CDO_MODE_ONLY
#define CS_DOMAIN_CDO_MODE_ONLY
Definition: cs_domain.h:59
cs_parameters_add_boundary_values
cs_field_t * cs_parameters_add_boundary_values(cs_field_t *f)
Define a boundary values field for a variable field.
Definition: cs_parameters.c:1346
cs_get_glob_thermal_model
cs_thermal_model_t * cs_get_glob_thermal_model(void)
Definition: cs_thermal_model.c:233
CS_EQKEY_ADV_SCHEME
Definition: cs_equation_param.h:646
cs_equation_param_t
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition: cs_equation_param.h:148
cs_field_t::val
cs_real_t * val
Definition: cs_field.h:145
cs_sles_set_post_output
void cs_sles_set_post_output(cs_sles_t *sles, int writer_id)
Activate postprocessing output for a given linear equation solver.
Definition: cs_sles.c:1383
CS_TIME_MOMENT_MEAN
Definition: cs_time_moment.h:58
cs_domain_def_time_step_by_value
void cs_domain_def_time_step_by_value(cs_domain_t *domain, double dt)
Define the value of the time step.
Definition: cs_domain_setup.c:408
eps
Definition: cs_field_pointer.h:71
cs_multigrid_set_coarsening_options
void cs_multigrid_set_coarsening_options(cs_multigrid_t *mg, int aggregation_limit, int coarsening_type, int n_max_levels, cs_gnum_t min_g_rows, double p0p1_relax, int postprocess)
Set multigrid coarsening parameters.
Definition: cs_multigrid.c:3097
cs_glob_mpi_comm
MPI_Comm cs_glob_mpi_comm
Definition: cs_defs.c:181
cs_mesh_quantities_t::cell_cen
cs_real_t * cell_cen
Definition: cs_mesh_quantities.h:92
cs_property_def_iso_by_value
cs_xdef_t * cs_property_def_iso_by_value(cs_property_t *pty, const char *zname, double val)
Define an isotropic cs_property_t structure by value for entities related to a volume zone.
Definition: cs_property.c:555
cs_thermal_model_field
cs_field_t * cs_thermal_model_field(void)
Definition: cs_thermal_model.c:205
cs_turb_model_t
Turbulence model general options descriptor.
Definition: cs_turbulence_model.h:75
CS_TIME_MOMENT_RESTART_AUTO
Definition: cs_time_moment.h:68
CS_ADVKEY_POST_COURANT
Definition: cs_advection_field.h:91
CS_DRIFT_SCALAR_ZERO_BNDY_FLUX
Definition: cs_parameters.h:166
cs_multigrid_pc_create
cs_sles_pc_t * cs_multigrid_pc_create(cs_multigrid_type_t mg_type)
Create a multigrid preconditioner.
Definition: cs_multigrid.c:3804
cs_turb_model_t::iturb
int iturb
Definition: cs_turbulence_model.h:77
CS_EQKEY_ITSOL_EPS
Definition: cs_equation_param.h:661
cs_field_get_key_struct
const void * cs_field_get_key_struct(const cs_field_t *f, const int key_id, void *s)
Return a structure for a given key associated with a field.
Definition: cs_field.c:3386
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...
Definition: cs_equation_param.c:2251
cs_boundary_set_default
void cs_boundary_set_default(cs_boundary_t *boundaries, cs_boundary_type_t type)
Set the default boundary related to the given cs_boundary_t structure.
Definition: cs_boundary.c:310
CS_PROPERTY_ANISO
Definition: cs_property.h:97
BFT_MALLOC
#define BFT_MALLOC(_ptr, _ni, _type)
Allocate memory for _ni elements of type _type.
Definition: bft_mem.h:62
cs_turb_rans_model_t
RANS turbulence model descriptor.
Definition: cs_turbulence_model.h:117
cs_field_set_key_int
int cs_field_set_key_int(cs_field_t *f, int key_id, int value)
Assign a integer value for a given key to a field.
Definition: cs_field.c:2930
cs_sles_get_context
void * cs_sles_get_context(cs_sles_t *sles)
Return pointer to solver context structure pointer.
Definition: cs_sles.c:1463
CS_SLES_PCG
Definition: cs_sles_it.h:59
cs_rotation_t
Subdomain rotation description.
Definition: cs_rotation.h:46
CS_F_
#define CS_F_(e)
Macro used to return a field pointer by its enumerated value.
Definition: cs_field_pointer.h:51
cs_piso_t
PISO options descriptor.
Definition: cs_parameters.h:205
cs_time_step_options_t
time step options descriptor
Definition: cs_time_step.h:80
CS_DRIFT_SCALAR_TURBOPHORESIS
Definition: cs_parameters.h:162
CS_EQKEY_ITSOL_RESNORM_TYPE
Definition: cs_equation_param.h:663
cs_post_boundary_flux
void cs_post_boundary_flux(const char *scalar_name, cs_lnum_t n_loc_b_faces, const cs_lnum_t b_face_ids[], cs_real_t b_face_flux[])
Compute scalar flux on a specific boundary region.
Definition: cs_post_util.c:844
CS_SLES_JACOBI
Definition: cs_sles_it.h:62
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.
Definition: cs_equation_param.c:2393
cs_property_t
Definition: cs_property.h:104
CS_VOF_MERKLE_MASS_TRANSFER
#define CS_VOF_MERKLE_MASS_TRANSFER
Definition: cs_vof.h:66
CS_SLES_BICGSTAB2
Definition: cs_sles_it.h:64
CS_PTYKEY_POST_FOURIER
Definition: cs_property.h:72
cs_domain_set_output_param
void cs_domain_set_output_param(cs_domain_t *domain, int nt_interval, int nt_list, int verbosity)
Set auxiliary parameters related to the way output is done.
Definition: cs_domain_setup.c:282
CS_POST_MONITOR
#define CS_POST_MONITOR
Definition: cs_post.h:62
cs_lnum_t
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:298
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.
Definition: cs_equation_param.c:2371
cs_field_set_key_double
int cs_field_set_key_double(cs_field_t *f, int key_id, double value)
Assign a floating point value for a given key to a field.
Definition: cs_field.c:3106
cs_glob_log_frequency
int cs_glob_log_frequency
cs_field_by_id
cs_field_t * cs_field_by_id(int id)
Return a pointer to a field based on its id.
Definition: cs_field.c:2307
cs_multigrid_define
cs_multigrid_t * cs_multigrid_define(int f_id, const char *name, cs_multigrid_type_t mg_type)
Define and associate a multigrid sparse linear system solver for a given field or equation name.
Definition: cs_multigrid.c:2841
cs_stokes_model_t
Stokes equation model descriptor.
Definition: cs_stokes_model.h:51
cs_equation_add_user
cs_equation_t * cs_equation_add_user(const char *eqname, const char *varname, int dim, cs_param_bc_type_t default_bc)
Add a new user equation structure and set a first set of parameters.
Definition: cs_equation.c:1078
cs_walldistance_activate
void cs_walldistance_activate(void)
Activate the future computation of the wall distance.
Definition: cs_walldistance.c:417
CS_EQKEY_ITSOL_MAX_ITER
Definition: cs_equation_param.h:662
CS_SLES_P_GAUSS_SEIDEL
Definition: cs_sles_it.h:66
cs_parameters_add_variable_variance
void cs_parameters_add_variable_variance(const char *name, const char *variable_name)
Define a user variable which is a variance of another variable.
Definition: cs_parameters.c:1123
cs_sles_t
struct _cs_sles_t cs_sles_t
Definition: cs_sles.h:68
cs_math_3_dot_product
static cs_real_t cs_math_3_dot_product(const cs_real_t u[3], const cs_real_t v[3])
Compute the dot product of two vectors of 3 real values.
Definition: cs_math.h:326
cs_get_glob_stokes_model
cs_stokes_model_t * cs_get_glob_stokes_model(void)
Definition: cs_stokes_model.c:361
cs_domain_set_cdo_mode
void cs_domain_set_cdo_mode(cs_domain_t *domain, int mode)
Set the global variable storing the mode of activation to apply to CDO/HHO schemes.
Definition: cs_domain.c:255
cs_domain_set_time_param
void cs_domain_set_time_param(cs_domain_t *domain, int nt_max, double t_max)
Set parameters for unsteady computations: the max number of time steps or the final physical time of ...
Definition: cs_domain_setup.c:309
cs_advection_field_by_name
cs_adv_field_t * cs_advection_field_by_name(const char *name)
Search in the array of advection field structures which one has the name given in argument.
Definition: cs_advection_field.c:259
cs_real_33_t
cs_real_t cs_real_33_t[3][3]
3x3 matrix of floating-point values
Definition: cs_defs.h:321
cs_vof_parameters_t::vof_model
unsigned vof_model
Definition: cs_vof.h:83
vel
void const cs_int_t *const const cs_int_t *const const cs_int_t *const const cs_int_t *const const cs_int_t *const const cs_int_t *const const cs_int_t *const const cs_int_t *const const cs_real_t *const const cs_real_t *const const cs_real_t const cs_real_t const cs_real_3_t vel[]
Definition: cs_divergence.h:64
CS_POST_UTIL_Q_CRITERION
Definition: cs_post_util.h:48
cs_var_cal_opt_t
structure containing the variable calculation options.
Definition: cs_parameters.h:60
cs_time_step_options_t::idtvar
int idtvar
Definition: cs_time_step.h:86
t
Definition: cs_field_pointer.h:98
cs_field_key_id
int cs_field_key_id(const char *name)
Return an id associated with a given key name.
Definition: cs_field.c:2490
CS_DRIFT_SCALAR_ADD_DRIFT_FLUX
Definition: cs_parameters.h:160
cs_var_cal_opt_t::iwarni
int iwarni
Definition: cs_parameters.h:61
cs_get_glob_vof_parameters
cs_vof_parameters_t * cs_get_glob_vof_parameters(void)
Definition: cs_vof.c:366
cp
Definition: cs_field_pointer.h:106
cs_parameters_add_variable
void cs_parameters_add_variable(const char *name, int dim)
Define a user variable.
Definition: cs_parameters.c:1082
CS_EQKEY_HODGE_DIFF_ALGO
Definition: cs_equation_param.h:656
CS_PROPERTY_ISO
Definition: cs_property.h:95
cs_piso_t::nterup
int nterup
Definition: cs_parameters.h:207
cs_grid_set_merge_options
void cs_grid_set_merge_options(int rank_stride, int cells_mean_threshold, cs_gnum_t cells_glob_threshold, int min_ranks)
Set global multigrid parameters for parallel grid merging behavior.
Definition: cs_grid.c:6388
CS_FIELD_VARIABLE
#define CS_FIELD_VARIABLE
Definition: cs_field.h:63
cs_advection_field_def_by_analytic
void cs_advection_field_def_by_analytic(cs_adv_field_t *adv, cs_analytic_func_t *func, void *input)
Define a cs_adv_field_t structure thanks to an analytic function.
Definition: cs_advection_field.c:637
CS_TIME_MOMENT_RESTART_RESET
Definition: cs_time_moment.h:67
cs_turb_rans_model_t::irijco
int irijco
Definition: cs_turbulence_model.h:155
CS_NVD_SUPERBEE
Definition: cs_convection_diffusion.h:90
CS_POST_WRITER_DEFAULT
#define CS_POST_WRITER_DEFAULT
Definition: cs_post.h:68
cs_adv_field_t
Definition: cs_advection_field.h:149
cs_parameters_add_property
void cs_parameters_add_property(const char *name, int dim, int location_id)
Define a user property.
Definition: cs_parameters.c:1160
cs_gwf_soil_add
cs_gwf_soil_t * cs_gwf_soil_add(const char *z_name, cs_gwf_soil_hydraulic_model_t model)
Create and add a new cs_gwf_soil_t structure. A first initialization of all members by default is per...
Definition: cs_gwf_soil.c:381
cs_field_t
Field descriptor.
Definition: cs_field.h:124
CS_TURB_K_EPSILON_LIN_PROD
Definition: cs_turbulence_model.h:56
CS_EQKEY_ITSOL
Definition: cs_equation_param.h:660
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.
Definition: cs_equation_param.c:2344
CS_GWF_SOIL_SATURATED
Definition: cs_gwf_soil.h:104
CS_BOUNDARY_WALL
Definition: cs_boundary.h:53
CS_MESH_LOCATION_CELLS
Definition: cs_mesh_location.h:63
CS_DRIFT_SCALAR_THERMOPHORESIS
Definition: cs_parameters.h:161
cs_advection_field_add_user
cs_adv_field_t * cs_advection_field_add_user(const char *name)
Add and initialize a new user-defined advection field structure.
Definition: cs_advection_field.c:310
cs_mesh_location_get_n_elts
const cs_lnum_t * cs_mesh_location_get_n_elts(int id)
Get a mesh location's number of elements.
Definition: cs_mesh_location.c:823
cs_field_n_fields
int cs_field_n_fields(void)
Return the number of defined fields.
Definition: cs_field.c:1527
cs_field_set_key_str
int cs_field_set_key_str(cs_field_t *f, int key_id, const char *str)
Assign a character string for a given key to a field.
Definition: cs_field.c:3219
p
Definition: cs_field_pointer.h:67
cs_sles_find_or_add
cs_sles_t * cs_sles_find_or_add(int f_id, const char *name)
Return pointer to linear system object, based on matching field id or system name.
Definition: cs_sles.c:1167