7.0
general documentation
cs_equation_param_t Struct Reference

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

#include <cs_equation_param.h>

+ Collaboration diagram for cs_equation_param_t:

Data Fields

General settings
char *restrict name
 
cs_equation_type_t type
 
int dim
 
union {
   int   verbosity
 
   int   iwarni
 
}; 
 
cs_flag_t flag
 
cs_flag_t process_flag
 
cs_param_space_scheme_t space_scheme
 
cs_param_dof_reduction_t dof_reduction
 
int space_poly_degree
 
Legacy Settings
int iconv
 
int istat
 
int idircl
 
int ndircl
 
int idiff
 
int idifft
 
int idften
 
int iswdyn
 
int ischcv
 
int ibdtso
 
int isstpc
 
int nswrgr
 
int nswrsm
 
int imrgra
 
int imligr
 
int ircflu
 
int iwgrec
 
int icoupl
 
double thetav
 
double blencv
 
double blend_st
 
double epsilo
 
double epsrsm
 
double epsrgr
 
double climgr
 
double extrag
 
double relaxv
 
Settings for the boundary conditions
cs_param_bc_type_t default_bc
 
int n_bc_defs
 
cs_xdef_t ** bc_defs
 
cs_param_bc_enforce_t default_enforcement
 
cs_real_t strong_pena_bc_coeff
 
cs_real_t weak_pena_bc_coeff
 
Numerical settings for the time-dependent parameters
int n_ic_defs
 
cs_xdef_t ** ic_defs
 
bool do_lumping
 
cs_hodge_param_t time_hodgep
 
cs_property_ttime_property
 
cs_param_time_scheme_t time_scheme
 
cs_real_t theta
 
Numerical settings for the diffusion (Laplacian div-grad) term
cs_hodge_param_t diffusion_hodgep
 
cs_property_tdiffusion_property
 
Numerical settings for the curl-curl term
cs_hodge_param_t curlcurl_hodgep
 
cs_property_tcurlcurl_property
 
Numerical settings for the grad-div term
cs_hodge_param_t graddiv_hodgep
 
cs_property_tgraddiv_property
 
Numerical settings for the advection term
cs_param_advection_form_t adv_formulation
 
cs_param_advection_scheme_t adv_scheme
 
cs_param_advection_strategy_t adv_strategy
 
cs_param_advection_extrapol_t adv_extrapol
 
cs_real_t upwind_portion
 
cs_adv_field_tadv_field
 
cs_property_tadv_scaling_property
 
Numerical settings for the reaction term

The contribution of a reaction term to the algebraic system is either at the left-hand and/or right-hand side according to the choice of time scheme

cs_hodge_param_t reaction_hodgep
 
int n_reaction_terms
 
cs_property_t ** reaction_properties
 
Definition of the related source terms

The contribution of a source term to the algebraic system is always at right-hand side whatever is the choice of time scheme

int n_source_terms
 
cs_xdef_t ** source_terms
 
Definition of the related volume mass injection

The contribution of a volume mass injection to the algebraic system is always at right-hand side whatever is the choice of time scheme, and is defined in a manner similar to boundary conditions. For variables whose injection value matches the "ambient" value, no term needs to be added.

int n_volume_mass_injections
 
cs_xdef_t ** volume_mass_injections
 
Enforcement of values inside the computational domain

This is different from the enforcement of boundary conditions but rely on the same algebraic manipulation. Two mechanisms are possible.

1) CELL SELECTION: defined a selection of cells and then automatically built the related selection of degrees of freedom and assigned a value to each selected degrees of freedom

2) DOF SELECTION: defined a selection of degrees of freedom (DoFs) and assign a values to a selection of degrees of freedom inside the domain

cs_flag_t enforcement_type
 
cs_real_tenforcement_ref_value
 
cs_lnum_t n_enforced_cells
 
cs_lnum_tenforced_cell_ids
 
cs_real_tenforced_cell_values
 
cs_lnum_t n_enforced_dofs
 
cs_lnum_tenforced_dof_ids
 
cs_real_tenforced_dof_values
 
Settings related to the resolution of the algebraic system
cs_param_sles_tsles_param
 
Settings related to the optimization of the performance
cs_param_assemble_omp_strategy_t omp_assembly_choice
 

Detailed Description

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

Field Documentation

◆ @10

union { ... }

◆ adv_extrapol

adv_extrapol

Extrapolation used to estimate the advection field (please refer to cs_param_advection_extrapol_t)

◆ adv_field

adv_field

Pointer to the cs_adv_field_t structure associated to the advection term

◆ adv_formulation

adv_formulation

Type of formulation (conservative, non-conservative...) for the advective term

◆ adv_scaling_property

adv_scaling_property

May be set to NULL even if the advection term is activated. The value of this property in each cell is multiplicative coefficient in front of the advection term (boundary terms are also considered) This is useful to treat the thermal module using the variable temperature instead of the enthalpy for instance or in the solidification module.

◆ adv_scheme

adv_scheme

Numerical scheme used for the discretization of the advection term

◆ adv_strategy

adv_strategy

Strategy used to handle the advection term (please refer to cs_param_advection_strategy_t)

◆ bc_defs

bc_defs

Pointers to the definitions of the boundary conditions

◆ blencv

blencv

Proportion of second-order convective scheme (0 corresponds to an upwind first-order scheme); in case of LES calculation, a second-order scheme is recommended and activated by default (blencv = 1).
Relevant where iconv = 1.

◆ blend_st

blend_st

Proportion of second-order convective scheme (0 corresponds to an upwind first-order scheme) after the slope test is activated; in case of LES calculation, a second-order scheme is recommended and activated by default (blend_st = 1).
Relevant whereiconv = 1.

◆ climgr

climgr

For least squares gradients, factor of gradient limitation (high value means little limitation).
Relevant for all the variables using least-squares gradientsfor which imligr > CS_GRADIENT_LIMIT_NONE.

◆ curlcurl_hodgep

curlcurl_hodgep

Set of parameters for the discrete Hodge operator related to the curl-curl operator

◆ curlcurl_property

curlcurl_property

Pointer to the property related to the curl-curl term

◆ default_bc

default_bc

Default boundary condition related to this equation. Valid choices:

◆ default_enforcement

default_enforcement

Type of strategy to enforce an essential boundary conditions (Dirichlet for instance) when no predefined strategy is prescribed. See cs_param_bc_enforce_t for more details.

◆ diffusion_hodgep

diffusion_hodgep

Set of parameters for the discrete Hodge operator related to the diffusion

◆ diffusion_property

diffusion_property

Pointer to the property related to the diffusion term

◆ dim

int dim

Dimension of the unknown

◆ do_lumping

do_lumping

Set to true or false. Activate several actions: Perform a mass lumping of the matrices related to the time and reaction discretization. All source terms are evaluated using a barycentric quadrature.

◆ dof_reduction

cs_param_dof_reduction_t dof_reduction

How is defined DoF

◆ enforced_cell_ids

enforced_cell_ids

List of selected cell ids related to an enforcement

◆ enforced_cell_values

cs_real_t* enforced_cell_values

◆ enforced_dof_ids

enforced_dof_ids

List of related DoF ids

◆ enforced_dof_values

enforced_dof_values

List of related values to enforce

◆ enforcement_ref_value

enforcement_ref_value

Reference value to use. Avod to allocate an array with the same value for all selected entities

◆ enforcement_type

enforcement_type

Flag specifying which kind of enforcement to perform

◆ epsilo

epsilo

Relative precision for the solution of the linear system. The default is epsilo = $ 10^-8 $ . When there are enough iterations on the reconstruction of the right-hand side of the equation, the value may be increased (by default, in case of second-order in time, with nswrsm = 5 or 10, epsilo is increased to $ 10^-5 $.

◆ epsrgr

epsrgr

Relative precision for the iterative gradient reconstruction. (when imrgra = 0).

◆ epsrsm

epsrsm

Relative precision on the reconstruction of the right hand-side. The default is epsrsm = $ 10^-8 $. When there are not enough iterations on the reconstruction of the right-hand side of the equation, the value may be increased (by default, in case of second-order in time, with nswrsm = 5 or 10, epsrsm is increased to $ 10^-5 $ ).

◆ extrag

extrag

Removed option, has no effect.

◆ flag

flag

Flag to know if unsteady or diffusion or convection or reaction or source terms are activated or not

◆ graddiv_hodgep

graddiv_hodgep

Set of parameters for the discrete Hodge operator related to the grad-div operator

◆ graddiv_property

graddiv_property

Pointer to the property related to the grad-div term

◆ ibdtso

ibdtso

Backward differential scheme in time order.

◆ ic_defs

ic_defs

List of pointers to the definition of the inititial condition

◆ iconv

iconv

Indicate if the convection is taken into account (1) or not (0). By default, 0 for the pressure or f in v2f model, 1 for the other unknowns.

◆ icoupl

icoupl

Internal coupling indicator

  • -1: not coupled (default)
  • 1: coupled

◆ idften

idften

Type of diffusivity flag (sum of mask constants defining if diffusivity is isotropic, anisotropic, ... Masks are defined in Transported scalars parameters).

◆ idiff

idiff

Indifcate if diffusion is taken into account (1) or not (0).

◆ idifft

idifft

When diffusion is taken into account (idiff = 1), indicate if the turbulent diffusion is taken into account (idifft = 1) or not (0).

◆ idircl

idircl

Indicate whether the diagonal of the matrix should be slightly shifted if there is no Dirichlet boundary condition and if istat = 0 (0: false / 1: true). Indeed, in such a case, the matrix for the general advection/diffusion equation is singular. A slight shift in the diagonal will make it invertible again. By default, idircl is set to 1 for all the unknowns, except $\overline{f}$ in the v2f model (whose equation already contain another diagonal term).

Remarks
this code is defined automatically based on the presence of Dirichlet BCs.

◆ imligr

imligr

Type of gradient limiter

  • -1 (CS_GRADIENT_LIMIT_NONE): no limitation
  • 0 (CS_GRADIENT_LIMIT_CELL): based on the neighbors
  • 1 (CS_GRADIENT_LIMIT_FACE): superior order
    imligr is applied only to least-squares gradients. In the case of the Green-Gauss gradient with least-squares based face gradients, applied to the least-squares step.

◆ imrgra

imrgra

Type of gradient reconstruction

  • 0: iterative reconstruction of the non-orthogonalities
  • 1: least squares method based on the first neighbor cells (those which share a face with the treated cell)
  • 2, 3: least squares method using the extended neighborhood
  • 4: Green-Gauss based using the least squares method (first neighbors) to compute face values
  • 5, 6: Green-Gauss based using the least squares method with and extended neighborhood to compute face values
    If the computation fails due to mesh quality aspects, it is usually effective to use imrgra = 3, 5, or 6.

◆ ircflu

ircflu

Indicate whether the convective and diffusive fluxes at the faces should be reconstructed:

  • 0: no reconstruction
  • 1: reconstruction Deactivating the reconstruction of the fluxes can have a stabilizing effect on the calculation. It is sometimes useful with the $ k-\epsilon $ model, if the mesh is strongly non-orthogonal in the near-wall region, where the gradients of k and $ \epsilon $ are strong. In such a case, setting ircflu = 0 will probably help (switching to a first order convective scheme, blencv = 0, for k and $ \epsilon $ might also help in that case).

◆ ischcv

ischcv

Indicate the type of second-order convective scheme

  • 0: Second Order Linear Upwind
  • 1: Centered
  • 2: Second Order with upwind-gradient reconstruction (SOLU)
  • 3: Blending between SOLU and Centered scheme
  • 4: NVD/TVD Scheme 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

◆ isstpc

isstpc

Indicate whether a slope test should be used to switch from a second-order to an upwind convective scheme under certain conditions, to ensure stability.

  • 0: slope test activated
  • 1: slope test deactivated
  • 2: continuous limiter ensuring boundedness (beta limiter) The use of the slope test stabilises the calculation but may reduce the spatial convergence order.

◆ istat

istat

Indicate whether unsteady terms are present (1) or not (0) in the matrices. By default, 0 for the pressure or f in v2f model, 1 for the other unknowns.

◆ iswdyn

iswdyn

Dynamic relaxation type:

  • 0 no dynamic relaxation
  • 1 dynamic relaxation depending on $ \delta \varia^k $
  • 2 dynamic relaxation depending on $ \delta \varia^k $ and $ \delta \varia^{k-1} $.

◆ iwarni

iwarni
Deprecated:
use verbosity instead (iwarni is an alias to verbosity)

◆ iwgrec

iwgrec

Gradient calculation weighting

  • 0: standard
  • 1: weighted

◆ n_bc_defs

n_bc_defs

Number of boundary conditions which are defined for this equation

◆ n_enforced_cells

n_enforced_cells

Number of selected cells related to an enforcement

◆ n_enforced_dofs

n_enforced_dofs

Number of degrees of freedom (DoFs) to enforce

◆ n_ic_defs

n_ic_defs

Number of definitions for setting the intial condition

◆ n_reaction_terms

n_reaction_terms

Number of reaction terms to consider.

◆ n_source_terms

n_source_terms

Number of source terms to consider.

◆ n_volume_mass_injections

n_volume_mass_injections

Number of volume injections to consider.

◆ name

char* restrict name

name of the equation

◆ ndircl

ndircl

Number of Dirichlet BCs

◆ nswrgr

nswrgr

Iteration limit for the iterative gradient reconstruction (imrgra = 0). If imrgra = 0 and nswrgr <= 1, gradients are not reconstructed.

◆ nswrsm

nswrsm

Iteration limit for the reconstruction of the right-hand sides of the equations with a first-order scheme in time (standard case), the default values are 2 for pressure and 1 for the other variables. With a second-order scheme in time (ischtp = 2) or LES, the default values are 5 for pressure and 10 for the other variables.

◆ omp_assembly_choice

omp_assembly_choice

When OpenMP is active, choice of parallel reduction for the assembly

◆ process_flag

process_flag

Flag to determine if predefined post-treatments such as Peclet, are requested

◆ reaction_hodgep

reaction_hodgep

Set of parameters for the discrete Hodge operator related to the reaction terms

◆ reaction_properties

reaction_properties

List of properties associated to each reaction term

◆ relaxv

relaxv

Relaxation coefficient for the associated variable. This relaxation parameter is only useful for the pressure with the unsteady algorithm (so as to improve the convergence in case of meshes of insufficient quality or of some turbulent models (k-epsilon, v2f, k-omega) and ikecou = 0; if ikecou = 1, relaxv is ignored.
Default values are 0.7 for turbulent variables and 1. for pressure. relaxv also stores the value of the relaxation coefficient when using the steady algorithm, deduced from the value of relxst (defaulting to relaxv = 1. - relxst).
Used only for the pressure and for turbulent variables ( $ k-\epsilon $, v2f or $ k-\omega $ models without coupling) with the unsteady algorithm. Always used with the steady algorithm.

◆ sles_param

sles_param

Set of parameters to specify how to to solve the algebraic system

  • iterative solver
  • preconditioner
  • tolerance...

◆ source_terms

source_terms

List of definition of each source term

◆ space_poly_degree

space_poly_degree

Maximum degree of the polynomial basis

◆ space_scheme

Space discretization scheme

◆ strong_pena_bc_coeff

strong_pena_bc_coeff

Value of the penalization coefficient used to enforce the Dirichlet boundary conditions when CS_PARAM_BC_ENFORCE_PENALIZED is set. This value should be sufficiently large in order to neglect off-diagonal terms.

◆ theta

theta

Value of the coefficient for a theta scheme (between 0 and 1)

◆ thetav

thetav

Value of $ \theta $ used to express at the second order the terms of convection, diffusion and the source terms which are linear functions of the solved variable, according to the formula $ \phi^{n+\theta} = (1-\theta) \phi^n + \theta \phi^{n+1} $. Generally, only the values 1 and 0.5 are used. The user is not allowed to modify this variable.

  • 1: first-order
  • 0.5: second-order
    For the pressure, thetav is always 1. For the other variables, thetav = 0.5 is used when the second-order time scheme is activated (ischtp = 2, standard for LES calculations), otherwise thetav = 1.

◆ time_hodgep

time_hodgep

Set of parameters for the discrete Hodge operator related to the unsteady term

◆ time_property

time_property

Pointer to the cs_property_t structure related to the unsteady term

◆ time_scheme

time_scheme

Type of numerical scheme used for the time discretization

◆ type

type of equation: predefined...

◆ upwind_portion

upwind_portion

Value between 0. and 1. (0: centered scheme, 1: pure upwind scheme) Introduce a constant portion of upwinding in a centered scheme Only useful if the advection scheme is set to CS_PARAM_ADVECTION_SCHEME_HYBRID_CENTERED_UPWIND

◆ verbosity

verbosity

Verbosity for the resolution (0 or 1 for a reasonable log size, 2 or more for troubleshooting).

◆ volume_mass_injections

volume_mass_injections

List of definitions of injection values

◆ weak_pena_bc_coeff

weak_pena_bc_coeff

Value of the penalization coefficient used to enforce the Dirichlet boundary condition when CS_PARAM_BC_ENFORCE_WEAK_NITSCHE or CS_PARAM_BC_ENFORCE_WEAK_SYM is set. This two latter strategies have a lesser influence on the conditioning number of the linear system than the choice CS_PARAM_BC_ENFORCE_PENALIZED


The documentation for this struct was generated from the following file: