8.3
general documentation
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
How to access and manage variables and properties using the cs_field API?

Variables and properties can be accessed using the cs_field API.

First, let us specify a few conventions used in code_saturne:

  • The dynamic variables designation includes velocity, pressure, and variables related to the turbulence model (k, ε, Rij, ω, ϕ, $ \overline{f} $, α, νt).
  • The designation “scalar” refers to (usually scalar) variables which are solution of an advection equation, apart from the dynamic variables listed above; for instance the temperature, scalars which may be passive or not, user-defined or model-defined.

    The mean value of the square of the fluctuations of a “scalar” is a “scalar”, too.

    Scalars may be divided into two groups: user-defined scalars and model-defined (sometimes referred to as “specific physics”) scalars.

  • Depending on the thermal model used; the solved thermal scalar variable may be the temperature, the enthalpy, or in the case of the compressible module, the total energy. When the temperature is not the solved thermal variable, it is usually available as a property.
  • When a user scalar represents the mean of the square of the fluctuations of another scalar (i.e. the average $ \overline{\varphi^\prime\varphi^\prime} $ for a fluctuating scalar ϕ). This can be made either via the interface or by declaring that scalar using cs_parameters_add_variable_variance in cs_user_parameters.c.
Accessing variables and properties in C:
  • The most general way of accessing a field is to use its name:

    cs_real_t *cvar_pr = cs_field_by_name("pressure")->val;
    then:
    cvar_pr[cell_id]

  • The main variables and properties can be accessed using the CS_F_ macro:

  • Indexed variables (such as user scalars) and indexed properties are accessed as:
    CS_FI_(name,ii-1)->val[cell_id].

  • The thermal scalar can be accessed using the cs_thermal_model_field function:
    field_t *tf = cs_thermal_model_field();
    then:
    tf->val[cell_id].
    if there is no thermal scalar, cs_thermal_model_field returns NULL.

  • In any case, all variables can be accessed using the function cs_field_by_name :
    cs_real_t *cvar_pr = cs_field_by_name("pressure")->val

  • A field's scalar_id key-word should be ≥ 0 if it is a scalar.

  • When a field is user-defined (rather than model-defined), its type flag should also match the CS_FIELD_USER mask:
    if (f->type & CS_FIELD_USER) ...
Remarks
Note that indexes in C begin at 0.

Cross-reference tables are available for the variables and properties of the standard solver and the specific physics features:

Variables

The C variables names are defined in cs_field_pointer.h.
Note that dt is accessible by its name (using the cs_field_by_name family of functions), but not through a permanent index..

C code Description
CS_F_(p)->val Pressure
CS_F_(vel)->val Velocity
CS_F_(void_f)->val Void fraction for Volume of Fluid model
CS_F_(k)->val Turbulent kinetic energy $ k $
CS_F_(eps)->val Turbulent dissipation $ \varepsilon $
CS_F_(rij)->val Reynolds stress tensor $ R $ (symmetric $ R_{ij} $ components, in the following order for each cell: xx, yy, zz, xy, yz, xz)
CS_F_(phi)->val $ \phi $ for $ \phi-f_b $ model
CS_F_(f_bar)->val $ f_b $ for $ \phi-f_b $ model
CS_F_(alp_bl)->val $ \alpha $ for $ Bl-v^2-k $
or EBRSM model
CS_F_(omg)->val $ \omega $ for $ k-\omega $ SST model
CS_F_(nusa)->val $ \widetilde{\nu}_T $ for Spalart-Allmaras
CS_F_(mesh_u)->val Mesh velocity
CS_F_(h)->val Enthalpy
CS_F_(t)->val Temperature

Properties

These properties are defined in cs_field_pointer.h.

C code Description
CS_F_(dt)->val Local time step
CS_F_(mu)->val Molecular viscosity
CS_F_(mu_t)->val Turbulent dynamic viscosity
CS_F_(cp)->val Specific heat
CS_F_(rho)->val Density (at cells)
CS_F_(rho_b)->val[face_id] Density (at boundary faces)
cs_real_t *cpro_smago = cs_field_by_name("smagorinsky_constant^2")->val Field id of the anisotropic turbulent viscosity
cs_real_t *cpro_cour = cs_field_by_name("courant_number")->val Courant number
cs_real_t *cpro_four = cs_field_by_name("fourier_number")->val Fourier number
cs_real_t *cpro_prtot = cs_field_by_name("total_pressure")->val Total pressure at cell centers
cs_real_t *cpro_visma_s = cs_field_by_name("mesh_viscosity")->val Mesh velocity viscosity (scalar) for the ALE module
cs_real_t *cpro_visma_v = cs_field_by_name("mesh_viscosity")->val Mesh velocity viscosity (vector) for the ALE module
cs_real_t *cpro_tsrho = cs_field_by_name("dila_st")->val Global dilatation source terms
cs_real_t *cpro_beta = cs_field_by_name("thermal_expansion")->val Thermal expansion coefficient
CS_F_(poro)->val Porosity
CS_F_(t_poro)->val Tensorial porosity
cs_real_t *bpro_forbr = cs_field_by_name("boundary_forces")->val Field id of the stresses at boundary
cs_real_t *bpro_yplus = cs_field_by_name("yplus")->val Field id of $y^+$ at boundary
cs_real_t *dttens = cs_field_by_name("dttens")->val Field id for the dttens tensor
CS_F_(t_b)->val Boundary temperature

Some physical properties such as specific heat or dynamic diffusivity are often constant (depending on the model or user parameters). In that case, these properties are stored as a simple real numbers rather than in a field over all mesh cells.

  • This is the case for the specific heat Cp.
    • If Cp is constant, it is based on the reference value cs_glob_fluid_properties->cp0.
      When this is the case, cs_glob_fluid_properties->icp should remain set to 0.
      When icp is initialized to 1 (by the GUI, or in cs_user_parameters.c, it is automatically reset to the id of the cell-based property field referenced in the above table.
  • This is the same for the dynamic diffusivity K of each scalar.
    • When K is constant, its value is based on the field's reference dynamic diffusivity, accessible through the scalar field's diffusivity_ref keyword.

      When it is variable, the matching field can be specified and accessed using the base scalar field's diffusivity_id key (accessible using cs_field_key_id("diffusivity_id")).
      For example, for a scalar field f:
      int k_f_id = cs_field_get_key_int(f, cs_field_key_id("diffusivity_id"));
      cs_field_t *kf = cs_field_by_id(k_f_id);

Reference physical values

Reference physical values represent either the fluid properties if they are constant, either simple mean values for the initialization if properties are variable and defined in cs_user_physical_properties. Reference values can be set in the GUI or in cs_user_parameters.c.

C code Description
cs_glob_fluid_properties->p0 Reference total pressure
cs_glob_fluid_properties->ro0 Reference density
cs_glob_fluid_properties->viscl0 Reference molecular dynamic viscosity
cs_glob_fluid_properties->cp0 Reference specific heat

Specific physical models

Atmospheric

Defined in cs_field_pointer.h. Note, esp_id = iesp -1.

C code Description
CS_F_(pot_t)->val Potential temperature
CS_F_(ym_w)->val Total water content
CS_F_(ntdrp)->val Total number of droplets
CS_FI_(chemistry,esp_id)->val Chemistry species (indexed)

Coal combustion

Defined in cs_field_pointer.h. Note, esp_id = iesp -1.

C code Description
CS_FI_(np,esp_id)->val Particles per kg for coal class
CS_FI_(xch,esp_id)->val Reactive coal mass fraction for coal class
CS_FI_(xck,esp_id)->val Coke mass fraction for coal class
CS_FI_(xwt,esp_id)->val Water mass fraction for coal class
CS_FI_(h2,esp_id)->val Mass enthalpy for coal class (permeatic case)
CS_FI_(f1m,esp_id)->val Mean value light volatiles for coal class
CS_FI_(f2m,esp_id)->val Mean value heavy volatiles for coal class
CS_F_(f4m)->val Oxydant 2 mass fraction
CS_F_(f5m)->val Oxydant 3 mass fraction
CS_F_(f6m)->val Water from coal drying mass fraction
CS_F_(f7m)->val Carbon from coal oxidyzed by O2 mass fraction
CS_F_(f8m)->val Carbon from coal gasified by CO2 mass fraction
CS_F_(f9m)->val Carbon from coal gasified by H2O mass fraction
CS_F_(fvp2m)->val f1f2 variance
CS_F_(yco2)->val CO2 fraction
CS_F_(yhcn)->val HCN fraction
CS_F_(yno)->val NO fraction
CS_F_(ynh3)->val NH3 enthalpy
CS_F_(hox)->val Ox enthalpy

Compressible

Defined in cs_field_pointer.h.

C code Description
CS_F_(e_tot)->val Total energy
CS_F_(t_kelvin)->val Temperature, in Kelvin

Electric arcs

Defined in cs_field_pointer.h. Note, esp_id = iesp -1.

C code Description
CS_F_(potr)->val Electric potential, real part
CS_F_(poti)->val Electric potential, imaginary part
CS_F_(potva)->val Vector potential
CS_FI_(ycoel,esp_id)->val Constituent mass fraction

Gas combustion

Defined in cs_field_pointer.h.

C code Description
CS_F_(fm)->val Mixture fraction
CS_F_(fp2m)->val Mixture fraction variance
CS_F_(fsm)->val Soot mass fraction
CS_F_(npm)->val Soot precursor number
CS_F_(ygfm)->val Fresh gas fraction
CS_F_(yfm)->val Mass fraction
CS_F_(yfp2m)->val Mass fraction variance
CS_F_(coyfp)->val Mass fraction covariance

Radiative transfer

Defined in cs_field_pointer.h. Note, esp_id = iesp -1.

C code Description
CS_F_(rad_energy)->val Radiative luminance
CS_F_(rad_q)->val Radiative flux
CS_FI_(rad_ets,esp_id)->val Radiative flux explicit source term
CS_FI_(rad_its,esp_id)->val Radiative flux implicit source term
CS_FI_(rad_abs,esp_id)->val Radiative absorption
CS_FI_(rad_emi,esp_id)->val Radiative emission
CS_FI_(rad_cak,esp_id)->val Radiative absorption coefficient
CS_F_(qinci)->val Radiative incident radiative flux density
CS_F_(xlam)->val Wall thermal conductivity
CS_F_(epa)->val Wall thickness
CS_F_(emissivity)->val Wall emissivity
CS_F_(fnet)->val Boundary radiative flux
CS_F_(fconv)->val Boundary radiative convective flux
CS_F_(hconv)->val Radiative exchange coefficient

Eulerian-Eulerian multiphase flows

Defined in cs_field_pointer.h.

C code Description
CS_FI_(yf_ncond,inc)->val Non-condensable gas mass fractions
CS_FI_(qp,ip)->val Particles turbulent kinetic energy Q2
CS_FI_(qfp,ip)->val Covariance of the turbulent Q12
CS_FI_(qfpxx,ip)->val XX component of qfp
CS_FI_(qfpxy,ip)->val XY component of qfp
CS_FI_(qfpxz,ip)->val XZ component of qfp
CS_FI_(qfpyx,ip)->val YX component of qfp
CS_FI_(qfpyy,ip)->val YY component of qfp
CS_FI_(qfpyz,ip)->val YZ component of qfp
CS_FI_(qfpzx,ip)->val ZX component of qfp
CS_FI_(qfpzy,ip)->val ZY component of qfp
CS_FI_(qfpzz,ip)->val ZZ component of qfp
CS_FI_(gamma,ip)->val Interfacial mass transfer
CS_FI_(ia,ip)->val Interfacial area
CS_FI_(x2,ip)->val Droplets x2
CS_FI_(d32,ip)->val Droplets Sauter mean diameter
CS_FI_(drag,ipcpl)->val Drag between phases
CS_FI_(ad_mass,ip)->val Added mass
CS_FI_(th_diff,ip)->val Thermal dynamic diffusivity ( $kg.m^{-1}.s^{-1}$)
CS_FI_(th_diff_t,ip)->val Turbulent thermal dynamic diffusivity ( $kg.m^{-1}.s^{-1}$)
CS_FI_(drho_dp,ip)->val dRho over dP
CS_FI_(drho_dh,ip)->val dRho over dH
CS_FI_(tau12_t,ip)->val Turbulent tau12 for particles
CS_FI_(lift,ip)->val Particles lift
CS_FI_(disp_t,ip)->val Particles turbulent dispersion
CS_FI_(drift_vel,ip)->val Particles drift velocity

Liste of reserved advanded field names.

A list can be found in liste of predefined fields.