My Project
programmer's documentation
Turbulence model source terms.

Additional right_hand side source terms for turbulence models

Turbulence source terms may be modified using the cs_user_source_terms user-defined function.

Local variables

/* mesh quantities */
const cs_lnum_t n_cells = cs_glob_mesh->n_cells;

Remaining initialization

Get the density array in cpro_rom

const cs_real_t *cpro_rom = CS_F_(rho)->val;

Get the array of the current turbulent variable and its name

/* Define pointer f to the current turbulent variable field */
const cs_var_cal_opt_t *var_cal_opt
cs_field_key_id("var_cal_opt"));

Example

Example of arbitrary additional source term for turbulence models (Source term on the TKE "k" here).

Source term for $\varia:$

\[ \rho \norm{\vol{\celli}} \frac{d(\varia)}{dt} = ... - \rho \norm{\vol{\celli}} \cdot ff - \rho \frac{ \varia}{\tau}\]

with $ ff $ = 3.0 an $ \tau $ = 4.0

Body

Note
The turbulence variable names are:
  • 'k' and 'epsilon' for the k-epsilon models
  • 'rij' and 'epsilon' for the Rij-epsilon LRR and S SG
  • 'rij', 'epsilon' and 'alpha' for the EBRSM
  • 'k', 'epsilon', 'phi' and 'f_bar' for the phi-model
  • 'k', 'epsilon', 'phi' and 'alpha' for the Bl-v2-k model
  • 'k' and 'omega' for the k-omega turbulence model
  • 'nu_tilda' for the Spalart Allmaras model

Calculation of the explicit and implicit source terms

if (f == CS_F_(k)) {
if (var_cal_opt->iwarni >= 1)
bft_printf(" User source terms for turbulence variable %s\n",
const cs_real_t ff = 3.0;
const cs_real_t tau = 4.0;
for (cs_lnum_t i = 0; i < n_cells; i++) {
/* explicit source terms */
st_exp[i] = -cpro_rom[i] * cell_f_vol[i] * ff;
/* implicit source terms */
st_imp[i] = -cpro_rom[i] * cell_f_vol[i] / tau;
}
}
cs_mesh_t::n_cells
cs_lnum_t n_cells
Definition: cs_mesh.h:73
f_id
void const int * f_id
Definition: cs_gui.h:292
cs_glob_mesh_quantities
cs_mesh_quantities_t * cs_glob_mesh_quantities
mesh::cell_f_vol
double precision, dimension(:), pointer cell_f_vol
Definition: mesh.f90:156
rho
Definition: cs_field_pointer.h:103
cs_real_t
double cs_real_t
Floating-point value.
Definition: cs_defs.h:302
cs_glob_mesh
cs_mesh_t * cs_glob_mesh
cs_field_get_key_struct_const_ptr
const void * cs_field_get_key_struct_const_ptr(const cs_field_t *f, int key_id)
Return a read-only pointer to a simple structure for a given key to a field.
Definition: cs_field.c:3533
CS_F_
#define CS_F_(e)
Macro used to return a field pointer by its enumerated value.
Definition: cs_field_pointer.h:51
cs_lnum_t
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:298
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_mesh_quantities_t::cell_vol
cs_real_t * cell_vol
Definition: cs_mesh_quantities.h:93
cs_var_cal_opt_t
structure containing the variable calculation options.
Definition: cs_parameters.h:60
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_var_cal_opt_t::iwarni
int iwarni
Definition: cs_parameters.h:61
cs_field_t
Field descriptor.
Definition: cs_field.h:124
coincl::ff
double precision, dimension(nmaxfm), save ff
Definition: coincl.f90:85
k
Definition: cs_field_pointer.h:70
bft_printf
int bft_printf(const char *const format,...)
Replacement for printf() with modifiable behavior.
Definition: bft_printf.c:140
cs_field_get_label
const char * cs_field_get_label(const cs_field_t *f)
Return a label associated with a field.
Definition: cs_field.c:4257