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 */
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:
Define a boundary values field for a variable field.
Definition cs_parameters.c:1346
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.
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).
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.
Structure storing the main features of the computational domain and pointers to the main geometrical ...
Definition cs_domain.h:87
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:
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" */
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.
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition cs_equation_param.h:148
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:
Define an anisotropic cs_property_t structure by value for entities related to a volume zone.
Definition cs_property.c:643
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:
Structure storing medata for defining a quantity in a very flexible way.
Definition cs_xdef.h:126
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.
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).
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.
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:
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)
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.
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.
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.
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
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).
Define and associate a PETSc linear system solver for a given field or equation name.
Definition cs_sles_petsc.c:536
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:
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:
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:
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:
Code_Saturne allows the calculation of temporal means or variances, either of expressions evaluated through a user function, or of expressions of the type . 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.
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. */
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.
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. */
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. */
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.
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.
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.
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.
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.
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: