8.1
general documentation
cs_prototypes.h File Reference
#include "cs_base.h"
#include "cs_domain.h"
#include "cs_field.h"
#include "cs_mesh.h"
#include "cs_mesh_quantities.h"
#include "cs_mesh_bad_cells.h"
#include "cs_probe.h"
#include "cs_time_control.h"
#include "cs_volume_zone.h"
+ Include dependency graph for cs_prototypes.h:

Go to the source code of this file.

Functions

void caltri (void)
 
void cpthp1 (const int *mode, cs_real_t *eh, cs_real_t *xesp, cs_real_t *f1mc, cs_real_t *f2mc, cs_real_t *tp)
 
void csinit (const int *irgpar, const int *nrgpar)
 
void findpt (const cs_lnum_t *ncelet, const cs_lnum_t *ncel, const cs_real_t *xyzcen, const cs_real_t *xx, const cs_real_t *yy, const cs_real_t *zz, cs_lnum_t *node, int *ndrang)
 
void initi1 (void)
 
int cs_add_model_field_indexes (int f_id)
 
void cs_add_model_thermal_field_indexes (int f_id)
 
void cs_atmo_chem_source_terms (int iscal, cs_real_t st_exp[], cs_real_t st_imp[])
 
void cs_at_source_term_for_inlet (cs_real_3_t st_exp[])
 Additional right-hand side source terms for momentum equation in case of free inlet. More...
 
void cs_cavitation_compute_source_term (const cs_real_t pressure[], const cs_real_t voidf[])
 Compute the vaporization source term $ \Gamma_V \left(\alpha, p\right) = m^+ + m^- $ using the Merkle model: More...
 
void cs_compressible_convective_mass_flux (int iterns, cs_real_t dt[], cs_real_t vela[][3])
 Solves the continuity equation in pressure formulation and then updates the density and the mass flux. More...
 
void cs_coal_bt2h (cs_lnum_t n_faces, const cs_lnum_t face_ids[], const cs_real_t t[], cs_real_t h[])
 Convert temperature to enthalpy at boundary for coal combustion. More...
 
void cs_coal_thfieldconv1 (int location_id, const cs_real_t eh[], cs_real_t tp[])
 Calculation of the gas temperature Function with the gas enthalpy and concentrations. More...
 
void cs_drift_convective_flux (int f_id, const cs_real_t dt[], cs_real_t imasfl[], cs_real_t bmasf[], cs_real_t divflu[])
 Compute the modified convective flux for scalars with a drift. More...
 
cs_real_tcs_get_cavitation_dgdp_st (void)
 Return pointer to cavitation "dgdpca" array. More...
 
cs_real_tcs_get_cavitation_gam (void)
 Return pointer to cavitation "gamcav" array. More...
 
void cs_hydrostatic_pressure_prediction (cs_real_t grdphd[][3], int iterns)
 Computes a hydrostatic pressure $ P_{hydro} $ solving an a priori simplified momentum equation: More...
 
void cs_lagr_status (int *model_flag, int *restart_flag, int *frozen_flag)
 
void cs_field_map_and_init_bcs (void)
 Initialize all field BC coefficients. More...
 
void cs_mass_flux_prediction (cs_real_t dt[])
 Update the convective mass flux before the Navier Stokes equations (prediction and correction steps). More...
 
void cs_physical_model_source_terms (int iscal, cs_real_t st_imp[], cs_real_t st_exp[])
 
void cs_prediction_mass_flux (cs_lnum_t ncesmp, cs_lnum_t icetsm[], cs_real_t dt[])
 Update the convective mass flux before the Navier Stokes equations (prediction and correction steps). More...
 
void cs_fortran_resize_aux_arrays (void)
 resize some Fortran auxiliary arrays More...
 
void cs_sat_coupling_exchange_at_cells (int f_id, cs_real_t st_exp[], cs_real_t st_imp[])
 
void cs_secondary_viscosity (cs_real_t secvif[], cs_real_t secvib[])
 Computes the secondary viscosity contribution $\kappa -\dfrac{2}{3} \mu$ in order to compute: More...
 
void cs_turbomachinery_update (void)
 Sync turbomachinery module components to global c turbomachinery structure. More...
 
void cs_wall_condensation_volume_exchange_surf_at_cells (cs_real_t *surf)
 Return condensing volume structures surface at each cell. More...
 
void cs_user_1d_wall_thermal (int iappel, int isuit1)
 
void cs_user_wall_condensation (int nvar, int nscal, int iappel)
 Source terms associated at the boundary faces and the neighboring cells with surface condensation. More...
 
void cs_user_boundary_conditions_setup (cs_domain_t *domain)
 Set boundary conditions to be applied. More...
 
void cs_user_boundary_conditions (cs_domain_t *domain, int bc_type[])
 User definition of boundary conditions. More...
 
void cs_user_boundary_conditions_ale (cs_domain_t *domain, int bc_type[], int ale_bc_type[], int impale[])
 User definition of boundary conditions for ALE. More...
 
void cs_user_extra_operations_initialize (cs_domain_t *domain)
 Initialize variables. More...
 
void cs_user_extra_operations (cs_domain_t *domain)
 This function is called at the end of each time step. More...
 
void cs_user_extra_operations_finalize (cs_domain_t *domain)
 This function is called at the end of the calculation. More...
 
void cs_user_fsi_structure_define (int is_restart, int n_structs, int *plot, cs_time_control_t *plot_time_control, cs_real_t *aexxst, cs_real_t *bexxst, cs_real_t *cfopre, cs_real_t xstr0[][3], cs_real_t vstr0[][3], cs_real_t xstreq[][3])
 Definition of internal mobile structures and corresponding initial conditions (initial displacement and velocity ). More...
 
void cs_user_fsi_structure_values (int n_structs, const cs_time_step_t *ts, const cs_real_t xstreq[][3], const cs_real_t xstr[][3], const cs_real_t vstr[][3], cs_real_t xmstru[][3][3], cs_real_t xcstru[][3][3], cs_real_t xkstru[][3][3], cs_real_t forstr[][3], cs_real_t dtstr[])
 Time-based settings for internal mobile structures. More...
 
void cs_user_fsi_structure_num (cs_domain_t *domain, int structure_num[])
 Define structure numbers for faces associated with internal or external (code_aster) structures. More...
 
void cs_user_head_losses (const cs_zone_t *zone, cs_real_t cku[][6])
 Compute GUI-defined head losses for a given volume zone. More...
 
void cs_user_initialization (cs_domain_t *domain)
 This function is called one time step to initialize problem. More...
 
void cs_user_internal_coupling (void)
 Define internal coupling options. More...
 
void cs_user_internal_coupling_add_volumes (cs_mesh_t *mesh)
 Define volumes as internal coupling zones. More...
 
void cs_user_internal_coupling_from_disjoint_meshes (cs_mesh_t *mesh)
 Define volumesi from separated meshes as internal coupling zones. More...
 
void cs_user_physical_properties (cs_domain_t *domain)
 This function is called each time step to define physical properties. More...
 
void cs_user_physical_properties_h_to_t (cs_domain_t *domain, const cs_zone_t *z, bool z_local, const cs_real_t h[restrict], cs_real_t t[restrict])
 User definition of enthalpy to temperature conversion. More...
 
void cs_user_physical_properties_t_to_h (cs_domain_t *domain, const cs_zone_t *z, bool z_local, const cs_real_t t[restrict], cs_real_t h[restrict])
 User definition of temperature to enthalpy conversion. More...
 
void cs_user_physical_properties_td_pressure (cs_real_t *td_p)
 User function to define a custom law for the thermodynamic pressure. More...
 
void cs_user_physical_properties_turb_viscosity (cs_domain_t *domain)
 User modification of the turbulence viscosity. More...
 
void cs_user_source_terms (cs_domain_t *domain, int f_id, cs_real_t *st_exp, cs_real_t *st_imp)
 Additional user-defined source terms for variable equations (momentum, scalars, turbulence...). More...
 
void cs_user_porosity (cs_domain_t *domain)
 Compute the porosity (volume factor $ \epsilon $ when the porosity model is activated. (cs_glob_porous_model > 0). More...
 
void cs_user_join (void)
 Define mesh joinings. More...
 
void cs_user_linear_solvers (void)
 Define linear solver options. More...
 
void cs_user_finalize_setup (cs_domain_t *domain)
 Define or modify output user parameters. For CDO schemes, specify the elements such as properties, advection fields, user-defined equations and modules which have been previously added. More...
 
void cs_user_mesh_bad_cells_tag (cs_mesh_t *mesh, cs_mesh_quantities_t *mesh_quantities)
 Tag bad cells within the mesh based on user-defined geometric criteria. More...
 
void cs_user_mesh_restart_mode (void)
 Force preprocessing behavior in case of restart. More...
 
void cs_user_mesh_input (void)
 Define mesh files to read and optional associated transformations. More...
 
void cs_user_mesh_modify (cs_mesh_t *mesh)
 Modify geometry and mesh. More...
 
void cs_user_mesh_boundary (cs_mesh_t *mesh)
 Insert boundaries into a mesh. More...
 
void cs_user_mesh_smoothe (cs_mesh_t *mesh)
 Mesh smoothing. More...
 
void cs_user_mesh_save (cs_mesh_t *mesh)
 Enable or disable mesh saving. More...
 
void cs_user_mesh_warping (void)
 Set options for cutting of warped faces. More...
 
void cs_user_mesh_modify_partial (cs_mesh_t *mesh, cs_mesh_quantities_t *mesh_quantities)
 Apply partial modifications to the mesh after the preprocessing and initial postprocessing mesh building stage. More...
 
void cs_user_mesh_cartesian_define (void)
 Define a cartesian mesh. More...
 
void cs_user_model (void)
 Select physical model options, including user fields. More...
 
void cs_user_numbering (void)
 Define advanced mesh numbering options. More...
 
void cs_user_parallel_io (void)
 Define parallel IO settings. More...
 
void cs_user_partition (void)
 Define advanced partitioning options. More...
 
void cs_user_matrix_tuning (void)
 Define sparse matrix tuning options. More...
 
void cs_user_parameters (cs_domain_t *domain)
 Define or modify general numerical and physical user parameters. More...
 
void cs_user_radiative_transfer_parameters (void)
 User function for input of radiative transfer module options. More...
 
void cs_user_radiative_transfer_bcs (cs_domain_t *domain, const int bc_type[], int isothp[], cs_real_t *tmin, cs_real_t *tmax, cs_real_t *tx, const cs_real_t dt[], const cs_real_t thwall[], const cs_real_t qincid[], cs_real_t hfcnvp[], cs_real_t flcnvp[], cs_real_t xlamp[], cs_real_t epap[], cs_real_t epsp[], cs_real_t textp[])
 User definition of radiative transfer boundary conditions. More...
 
void cs_user_periodicity (void)
 Define periodic faces. More...
 
void cs_user_postprocess_writers (void)
 Define post-processing writers. More...
 
void cs_user_postprocess_probes (void)
 Define monitoring probes and profiles. More...
 
void cs_user_postprocess_meshes (void)
 Define post-processing meshes. More...
 
void cs_user_postprocess_values (const char *mesh_name, int mesh_id, int cat_id, cs_probe_set_t *probes, cs_lnum_t n_cells, cs_lnum_t n_i_faces, cs_lnum_t n_b_faces, cs_lnum_t n_vertices, const cs_lnum_t cell_list[], const cs_lnum_t i_face_list[], const cs_lnum_t b_face_list[], const cs_lnum_t vertex_list[], const cs_time_step_t *ts)
 User function for output of values on a post-processing mesh. More...
 
void cs_user_postprocess_activate (int nt_max_abs, int nt_cur_abs, double t_cur_abs)
 
void cs_user_rad_transfer_absorption (const int bc_type[], cs_real_t ck[])
 Absorption coefficient for radiative module. More...
 
void cs_user_rad_transfer_net_flux (const int itypfb[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t twall[], const cs_real_t qincid[], const cs_real_t xlam[], const cs_real_t epa[], const cs_real_t eps[], const cs_real_t ck[], cs_real_t net_flux[])
 Compute the net radiation flux. More...
 
int cs_user_solver_set (void)
 Set user solver. More...
 
void cs_user_solver (const cs_mesh_t *mesh, const cs_mesh_quantities_t *mesh_quantities)
 Main call to user solver. More...
 
void cs_user_saturne_coupling (void)
 Define couplings with other instances of code_saturne. More...
 
void cs_user_syrthes_coupling (void)
 Define couplings with SYRTHES code. More...
 
void cs_user_syrthes_coupling_volume_h (int coupling_id, const char *syrthes_name, cs_lnum_t n_elts, const cs_lnum_t elt_ids[], cs_real_t h_vol[])
 Compute a volume exchange coefficient for SYRTHES couplings. More...
 
void cs_user_cathare_coupling (void)
 Define couplings with SYRTHES code. More...
 
void cs_user_time_moments (void)
 Define time moments. More...
 
void cs_user_time_table (void)
 Define time tables using C API. This function is called at the begin of the simulation only. More...
 
void cs_user_turbomachinery (void)
 Define rotor/stator model. More...
 
void cs_user_turbomachinery_rotor (void)
 Define rotor axes, associated cells, and rotor/stator faces. More...
 
void cs_user_turbomachinery_set_rotation_velocity (void)
 Define rotation velocity of rotor. More...
 
void cs_user_zones (void)
 Define volume and surface zones. More...
 
void cs_user_scaling_elec (const cs_mesh_t *mesh, const cs_mesh_quantities_t *mesh_quantities, cs_real_t *dt)
 Define scaling parameter for electric model. More...
 
void cs_user_hgn_thermo_relax_time (const cs_mesh_t *mesh, const cs_real_t *alpha_eq, const cs_real_t *y_eq, const cs_real_t *z_eq, const cs_real_t *ei, const cs_real_t *v, cs_real_t *relax_tau)
 Computation of the relaxation time-scale. More...
 
void cs_user_paramedmem_define_couplings (void)
 Define ParaMEDMEM coupling(s) More...
 
void cs_user_paramedmem_define_meshes (void)
 Define coupled meshes. More...
 
void cs_user_paramedmem_define_fields (void)
 Define fields to couple with ParaMEDMEM. More...
 

Function Documentation

◆ caltri()

void caltri ( void  )

◆ cpthp1()

void cpthp1 ( const int *  mode,
cs_real_t eh,
cs_real_t xesp,
cs_real_t f1mc,
cs_real_t f2mc,
cs_real_t tp 
)

◆ cs_add_model_field_indexes()

int cs_add_model_field_indexes ( int  f_id)

◆ cs_add_model_thermal_field_indexes()

void cs_add_model_thermal_field_indexes ( int  f_id)

◆ cs_at_source_term_for_inlet()

void cs_at_source_term_for_inlet ( cs_real_3_t  st_exp[])

Additional right-hand side source terms for momentum equation in case of free inlet.

Parameters
[in]st_expexplicit part of the momentum source term

◆ cs_atmo_chem_source_terms()

void cs_atmo_chem_source_terms ( int  iscal,
cs_real_t  st_exp[],
cs_real_t  st_imp[] 
)

◆ cs_cavitation_compute_source_term()

void cs_cavitation_compute_source_term ( const cs_real_t  pressure[],
const cs_real_t  voidf[] 
)

Compute the vaporization source term $ \Gamma_V \left(\alpha, p\right) = m^+ + m^- $ using the Merkle model:

\[ m^+ = -\dfrac{C_{prod} \rho_l \min \left( p-p_V,0 \right)\alpha(1-\alpha)} {0.5\rho_lu_\infty^2t_\infty}, \]

\[ m^- = -\dfrac{C_{dest} \rho_v \max \left( p-p_V,0 \right)\alpha(1-\alpha)} {0.5\rho_lu_\infty^2t_\infty}, \]

with $ C_{prod}, C_{dest} $ empirical constants, $ t_\infty=l_\infty/u_\infty $ a reference time scale and $p_V$ the reference saturation pressure. $ l_\infty $, $ u_\infty $ and $p_V$ may be provided by the user (user function). Note that the r.h.s. of the void fraction transport equation is $ \Gamma_V/\rho_v $.

Parameters
[in]pressurepressure array
[in]voidfvoid fraction array

◆ cs_coal_bt2h()

void cs_coal_bt2h ( cs_lnum_t  n_faces,
const cs_lnum_t  face_ids[],
const cs_real_t  t[],
cs_real_t  h[] 
)

Convert temperature to enthalpy at boundary for coal combustion.

Parameters
[in]n_facesnumber of faces in list
[in]face_idslist of boundary faces at which conversion is requested (0-based numbering)
[in]t_btemperature at boundary
[out]h_benthalpy at boundary

◆ cs_coal_thfieldconv1()

void cs_coal_thfieldconv1 ( int  location_id,
const cs_real_t  eh[],
cs_real_t  tp[] 
)

Calculation of the gas temperature Function with the gas enthalpy and concentrations.

Parameters
[in]location_idmesh location id (cells or boundary faces)
[in]ehgas enthalpy ( $ j . kg $ of gaseous mixture)
[in,out]tpgas temperature (in kelvin)

◆ cs_compressible_convective_mass_flux()

void cs_compressible_convective_mass_flux ( int  iterns,
cs_real_t  dt[],
cs_real_t  vela[][3] 
)

Solves the continuity equation in pressure formulation and then updates the density and the mass flux.

Parameters
[in]iternsNavier-Stokes iteration number
[in]dttime step (per cell)
[in]velavelocity value at time step beginning

◆ cs_drift_convective_flux()

void cs_drift_convective_flux ( int  f_id,
const cs_real_t  dt[],
cs_real_t  imasfl[],
cs_real_t  bmasf[],
cs_real_t  divflu[] 
)

Compute the modified convective flux for scalars with a drift.

Parameters
[in]iflidindex of the current drift scalar field
[in]dttime step (per cell)
[in,out]imasflscalar mass flux at interior face centers
[in,out]bmasflscalar mass flux at boundary face centers
[in,out]divfludivergence of drift flux

◆ cs_field_map_and_init_bcs()

void cs_field_map_and_init_bcs ( void  )

Initialize all field BC coefficients.

◆ cs_fortran_resize_aux_arrays()

void cs_fortran_resize_aux_arrays ( void  )

resize some Fortran auxiliary arrays

◆ cs_get_cavitation_dgdp_st()

cs_real_t* cs_get_cavitation_dgdp_st ( void  )

Return pointer to cavitation "dgdpca" array.

Returns
pointer to "dgdpca" array.

◆ cs_get_cavitation_gam()

cs_real_t* cs_get_cavitation_gam ( void  )

Return pointer to cavitation "gamcav" array.

Returns
pointer to "gamcav" array.

◆ cs_hydrostatic_pressure_prediction()

void cs_hydrostatic_pressure_prediction ( cs_real_t  grdphd[][3],
int  iterns 
)

Computes a hydrostatic pressure $ P_{hydro} $ solving an a priori simplified momentum equation:

Parameters
[out]grdphdthe a priori hydrostatic pressure gradient $ \partial _x (P_{hydro}) $
[in]iternsNavier-Stokes iteration number

◆ cs_lagr_status()

void cs_lagr_status ( int *  model_flag,
int *  restart_flag,
int *  frozen_flag 
)

◆ cs_mass_flux_prediction()

void cs_mass_flux_prediction ( cs_real_t  dt[])

Update the convective mass flux before the Navier Stokes equations (prediction and correction steps).

This function computes a potential $ \varia $ solving the equation:

\[ D \left( \Delta t, \varia \right) = \divs \left( \rho \vect{u}^n\right) - \Gamma^n + \dfrac{\rho^n - \rho^{n-1}}{\Delta t} \]

This potential is then used to update the mass flux as follows:

\[ \dot{m}^{n+\frac{1}{2}}_\ij = \dot{m}^{n}_\ij - \Delta t \grad_\fij \varia \cdot \vect{S}_\ij \]

Parameters
[in]dttime step (per cell)

◆ cs_physical_model_source_terms()

void cs_physical_model_source_terms ( int  iscal,
cs_real_t  st_imp[],
cs_real_t  st_exp[] 
)

◆ cs_prediction_mass_flux()

void cs_prediction_mass_flux ( cs_lnum_t  ncesmp,
cs_lnum_t  icetsm[],
cs_real_t  dt[] 
)

Update the convective mass flux before the Navier Stokes equations (prediction and correction steps).

This function computes a potential $ \varia $ solving the equation:

\[ D \left( \Delta t, \varia \right) = \divs \left( \rho \vect{u}^n\right) - \Gamma^n + \dfrac{\rho^n - \rho^{n-1}}{\Delta t} \]

This potential is then used to update the mass flux as follows:

\[ \dot{m}^{n+\frac{1}{2}}_\ij = \dot{m}^{n}_\ij - \Delta t \grad_\fij \varia \cdot \vect{S}_\ij \]

Parameters
[in]ncesmpnumber of cells with mass source term
[in]icetsmindex of cells with mass source term
[in]dttime step (per cell)

◆ cs_sat_coupling_exchange_at_cells()

void cs_sat_coupling_exchange_at_cells ( int  f_id,
cs_real_t  st_exp[],
cs_real_t  st_imp[] 
)

◆ cs_secondary_viscosity()

void cs_secondary_viscosity ( cs_real_t  secvif[],
cs_real_t  secvib[] 
)

Computes the secondary viscosity contribution $\kappa -\dfrac{2}{3} \mu$ in order to compute:

\[ \grad\left( (\kappa -\dfrac{2}{3} \mu) \trace( \gradt(\vect{u})) \right) \]

with:

  • $ \mu = \mu_{laminar} + \mu_{turbulent} $
  • $ \kappa $ is the volume viscosity (generally zero)
Remarks
In LES, the tensor $\overline{\left(\vect{u}-\overline{\vect{u}}\right)\otimes\left(\vect{u} *-\overline{\vect{u}}\right)}$ is modeled by $\mu_t \overline{\tens{S}}$ and not by $\mu_t\overline{\tens{S}}-\dfrac{2}{3}\mu_t \trace\left(\overline{\tens{S}}\right)\tens{1}+\dfrac{2}{3}k\tens{1}$ so that no term $\mu_t \dive \left(\overline{\vect{u}}\right)$ is needed.

Please refer to the visecv section of the theory guide for more informations.

Parameters
[in,out]secviflambda*surface at interior faces
[in,out]secviblambda*surface at boundary faces

◆ cs_turbomachinery_update()

void cs_turbomachinery_update ( void  )

Sync turbomachinery module components to global c turbomachinery structure.

◆ cs_user_1d_wall_thermal()

void cs_user_1d_wall_thermal ( int  iappel,
int  isuit1 
)

◆ cs_user_boundary_conditions()

void cs_user_boundary_conditions ( cs_domain_t domain,
int  bc_type[] 
)

User definition of boundary conditions.

Parameters
[in,out]domainpointer to a cs_domain_t structure
[in,out]bc_typeboundary face types

The icodcl and rcodcl arrays are pre-initialized based on default and GUI-defined definitions, and may be modified here.

For a given variable field f, and a given face "face_id", these arrays may be used as follows:

  • Boundary condition type code given at: f->bc_coeffs->icodcl[face_id]
  • Dirichlet value defined at: f->bc_coeffs->rcodcl1[face_id]
  • Interior exchange coefficient (infinite if no exchange) at: f->bc_coeffs->rcodcl2[face_id]
  • Flux density defined at: f->bc_coeffs->rcodcl3[face_id]

For vector or tensor fields, these arrays are not interleaved, so for a given face "face_id" and field component "comp_id", access is as follows (where n_b_faces is domain->mesh->n_b_faces):

f->bc_coeffs->rcodcl1[n_b_faces*comp_id + face_id]
f->bc_coeffs->rcodcl2[n_b_faces*comp_id + face_id]
f->bc_coeffs->rcodcl3[n_b_faces*comp_id + face_id]

Only the icodcl code values from the first component are used in the case of vector or tensor fields, so the icodcl values can be defined as for a scalar.

◆ cs_user_boundary_conditions_ale()

void cs_user_boundary_conditions_ale ( cs_domain_t domain,
int  bc_type[],
int  ale_bc_type[],
int  impale[] 
)

User definition of boundary conditions for ALE.

See Boundary conditions for ALE (Arbitrary Lagrangian Eulerian) for additional details.

Parameters
[in,out]domainpointer to a cs_domain_t structure
[in,out]bc_typeboundary face types
[in,out]ale_bc_typeboundary face types for mesh velocity (see cs_boundary_ale_subtype_bits_t)
[in]impaleindicator for prescribed node displacement (0 or 1)

The icodcl and rcodcl arrays are pre-initialized based on default and GUI-defined definitions, and may be modified here.

For a given variable field f, and a given face "face_id", these arrays may be used as follows:

  • Boundary condition type code given at: f->bc_coeffs->icodcl[face_id]
  • Dirichlet value defined at: f->bc_coeffs->rcodcl1[face_id]
  • Interior exchange coefficient (infinite if no exchange) at: f->bc_coeffs->rcodcl2[face_id]
  • Flux density defined at: f->bc_coeffs->rcodcl3[face_id]

For vector or tensor fields, these arrays are not interleaved, so for a given face "face_id" and field component "comp_id", access is as follows (where n_b_faces is domain->mesh->n_b_faces):

f->bc_coeffs->rcodcl1[n_b_faces*comp_id + face_id]
f->bc_coeffs->rcodcl2[n_b_faces*comp_id + face_id]
f->bc_coeffs->rcodcl3[n_b_faces*comp_id + face_id]

Only the icodcl code values from the first component are used in the case of vector or tensor fields, so the icodcl values can be defined as for a scalar.

◆ cs_user_boundary_conditions_setup()

void cs_user_boundary_conditions_setup ( cs_domain_t domain)

Set boundary conditions to be applied.

This function is called just before cs_user_finalize_setup, and boundary conditions can be defined in either of those functions, depending on whichever is considered more readable or practical for a given use.

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_cathare_coupling()

void cs_user_cathare_coupling ( void  )

Define couplings with SYRTHES code.

This is done by calling the cs_syr_coupling_define function for each coupling to add.

◆ cs_user_extra_operations()

void cs_user_extra_operations ( cs_domain_t domain)

This function is called at the end of each time step.

It has a very general purpose, although it is recommended to handle mainly postprocessing or data-extraction type operations.

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_extra_operations_finalize()

void cs_user_extra_operations_finalize ( cs_domain_t domain)

This function is called at the end of the calculation.

It has a very general purpose, although it is recommended to handle mainly postprocessing or data-extraction type operations.

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_extra_operations_initialize()

void cs_user_extra_operations_initialize ( cs_domain_t domain)

Initialize variables.

This function is called at beginning of the computation (restart or not) before the time step loop.

This is intended to initialize or modify (when restarted) variable and time step values.

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_finalize_setup()

void cs_user_finalize_setup ( cs_domain_t domain)

Define or modify output user parameters. For CDO schemes, specify the elements such as properties, advection fields, user-defined equations and modules which have been previously added.

Define or modify output user parameters. For CDO schemes, specify the elements such as properties, advection fields, user-defined equations and modules which have been previously added.

For CDO schemes, this function concludes the setup of properties, equations, source terms...

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_fsi_structure_define()

void cs_user_fsi_structure_define ( int  is_restart,
int  n_structs,
int *  plot,
cs_time_control_t plot_time_control,
cs_real_t aexxst,
cs_real_t bexxst,
cs_real_t cfopre,
cs_real_t  xstr0[][3],
cs_real_t  vstr0[][3],
cs_real_t  xstreq[][3] 
)

Definition of internal mobile structures and corresponding initial conditions (initial displacement and velocity ).

Parameters
[in]is_restartindicate if computation is restarted
[in]n_structsnumber of mobile structures
[in,out]plot;monitoring format mask 0: no plot 1: plot to text (.dat) format 2: plot to .csv format 3: plot to both formats
[in,out]plot_time_controlplot time output frequency control
[in,out]aexxstcoefficient for predicted displacement
[in,out]bexxstcoefficient for predicted displacement
[in,out]cfoprecoefficient for predicted force
[in,out]xstr0initial displacement per structure
[in,out]vstr0initial velocity per structure
[in,out]xstreqdisplacement of initial mesh relative to structures position at equilibrium
[in,out]plot_time_controltime control associated to plotting

◆ cs_user_fsi_structure_num()

void cs_user_fsi_structure_num ( cs_domain_t domain,
int  structure_num[] 
)

Define structure numbers for faces associated with internal or external (code_aster) structures.

Structure numbers associated to a given face have the following values:

  • -i where coupled to i-th (1-to n) external (code_aster) structure.
  • 0 where not coupled with an internal or external structure.
  • i where coupled to i-th (1-to n) internal (mass-spring) structure.
Parameters
[in,out]domainpointer to a cs_domain_t structure
[in,out]structure_numstructure id associated to each face

◆ cs_user_fsi_structure_values()

void cs_user_fsi_structure_values ( int  n_structs,
const cs_time_step_t ts,
const cs_real_t  xstreq[][3],
const cs_real_t  xstr[][3],
const cs_real_t  vstr[][3],
cs_real_t  xmstru[][3][3],
cs_real_t  xcstru[][3][3],
cs_real_t  xkstru[][3][3],
cs_real_t  forstr[][3],
cs_real_t  dtstr[] 
)

Time-based settings for internal mobile structures.

Parameters
[in]n_structsnumber of mobile structures
[in]tstime step structure
[in]xstreqdisplacement of initial mesh rel. to equilibrium
[in]xstrstructural displacement
[in]vstrstructural velocity
[in,out]xmstrumatrix of structural mass
[in,out]xcstrumatrix of structural friction
[in,out]xkstrumatrix of structural stiffness
[in,out]forstrforces acting on structures (take forces)
[in,out]dtstrstructural time step

◆ cs_user_head_losses()

void cs_user_head_losses ( const cs_zone_t zone,
cs_real_t  cku[][6] 
)

Compute GUI-defined head losses for a given volume zone.

Head loss tensor coefficients for each cell are organized as follows: cku11, cku22, cku33, cku12, cku13, cku23.

Parameters
[in]zonepointer to zone structure
[in,out]ckuhead loss coefficients

Compute GUI-defined head losses for a given volume zone.

Head loss tensor coefficients for each cell are organized as follows: ck11, ck22, ck33, ck12, ck13, ck23.

Coefficients are set to zero (then computed based on definitions provided through the GUI if this is the case) before calling this function, so setting values to zero is usually not necessary, unless we want to fully overwrite a GUI-based definition.

Diagonal coefficients must be positive; the calculation may diverge if this is not the case.

Parameters
[in]zonepointer to zone structure
[in,out]ckuhead loss coefficients

◆ cs_user_hgn_thermo_relax_time()

void cs_user_hgn_thermo_relax_time ( const cs_mesh_t mesh,
const cs_real_t alpha_eq,
const cs_real_t y_eq,
const cs_real_t z_eq,
const cs_real_t ei,
const cs_real_t v,
cs_real_t relax_tau 
)

Computation of the relaxation time-scale.

This function computes the value of the relaxation time-scale (for the return to equilibrium).

Parameters
[in]meshpointer to mesh
[in]alpha_eqequilibrium volume fraction
[in]y_eqequilibrium mass fraction
[in]z_eqequilibrium energy fraction
[in]eispecific internal energy
[in]vspecific volume
[in]relax_taurelaxation time scale towards equilibrium

◆ cs_user_initialization()

void cs_user_initialization ( cs_domain_t domain)

This function is called one time step to initialize problem.

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_internal_coupling()

void cs_user_internal_coupling ( void  )

Define internal coupling options.

Options are usually defined using cs_internal_coupling_add_entity.

◆ cs_user_internal_coupling_add_volumes()

void cs_user_internal_coupling_add_volumes ( cs_mesh_t mesh)

Define volumes as internal coupling zones.

These zones will be separated from the rest of the domain using automatically defined thin walls.

Parameters
[in,out]meshpointer to a cs_mesh_t structure

These zones will be separated from the rest of the domain using automatically defined thin walls.

Deprecated:
move contents tocs_user_internal_coupling instead.
Parameters
[in,out]meshpointer to a cs_mesh_t structure

◆ cs_user_internal_coupling_from_disjoint_meshes()

void cs_user_internal_coupling_from_disjoint_meshes ( cs_mesh_t mesh)

Define volumesi from separated meshes as internal coupling zones.

These zones must be disjoint and the face selection criteria must be specified.

Parameters
[in,out]meshpointer to a cs_mesh_t structure

Define volumesi from separated meshes as internal coupling zones.

These zones must be disjoint and the face selection criteria must be specified.

Deprecated:
move contents tocs_user_internal_coupling instead.
Parameters
[in,out]meshpointer to a cs_mesh_t structure

◆ cs_user_join()

void cs_user_join ( void  )

Define mesh joinings.

◆ cs_user_linear_solvers()

void cs_user_linear_solvers ( void  )

Define linear solver options.

This function is called at the setup stage, once user and most model-based fields are defined.

Available native iterative linear solvers include conjugate gradient, Jacobi, BiCGStab, BiCGStab2, and GMRES. For symmetric linear systems, an algebraic multigrid solver is available (and recommended).

External solvers may also be setup using this function, the cs_sles_t mechanism allowing such through user-define functions.

◆ cs_user_matrix_tuning()

void cs_user_matrix_tuning ( void  )

Define sparse matrix tuning options.

◆ cs_user_mesh_bad_cells_tag()

void cs_user_mesh_bad_cells_tag ( cs_mesh_t mesh,
cs_mesh_quantities_t mesh_quantities 
)

Tag bad cells within the mesh based on user-defined geometric criteria.

Parameters
[in,out]meshpointer to a cs_mesh_t structure
[in,out]mesh_quantitiespointer to a cs_mesh_quantities_t structure

◆ cs_user_mesh_boundary()

void cs_user_mesh_boundary ( cs_mesh_t mesh)

Insert boundaries into a mesh.

Parameters
[in,out]meshpointer to a cs_mesh_t structure

◆ cs_user_mesh_cartesian_define()

void cs_user_mesh_cartesian_define ( void  )

Define a cartesian mesh.

◆ cs_user_mesh_input()

void cs_user_mesh_input ( void  )

Define mesh files to read and optional associated transformations.

◆ cs_user_mesh_modify()

void cs_user_mesh_modify ( cs_mesh_t mesh)

Modify geometry and mesh.

Parameters
[in,out]meshpointer to a cs_mesh_t structure

◆ cs_user_mesh_modify_partial()

void cs_user_mesh_modify_partial ( cs_mesh_t mesh,
cs_mesh_quantities_t mesh_quantities 
)

Apply partial modifications to the mesh after the preprocessing and initial postprocessing mesh building stage.

Parameters
[in,out]meshpointer to a cs_mesh_t structure
[in,out]mesh_quantitiespointer to a cs_mesh_quantities_t structure

Apply partial modifications to the mesh after the preprocessing and initial postprocessing mesh building stage.

Parameters
[in,out]meshpointer to a cs_mesh_t structure
[in,out]mesh_quantitiespointer to a cs_mesh_quantities_t structure

◆ cs_user_mesh_restart_mode()

void cs_user_mesh_restart_mode ( void  )

Force preprocessing behavior in case of restart.

By default, in case of restart, if a "restart/mesh_input.csm" file is present, it will be read and proprocessing will be skipped.

This behavior may be changed in the GUI (in the "Mesh" section, unchecking "Use unmodified checkpoint mesh in case of restart"), or by calling cs_preprocessor_data_set_restart_mode in this function.

◆ cs_user_mesh_save()

void cs_user_mesh_save ( cs_mesh_t mesh)

Enable or disable mesh saving.

By default, mesh is saved when modified.

Parameters
[in,out]meshpointer to a cs_mesh_t structure

◆ cs_user_mesh_smoothe()

void cs_user_mesh_smoothe ( cs_mesh_t mesh)

Mesh smoothing.

Parameters
[in,out]meshpointer to a cs_mesh_t structure

◆ cs_user_mesh_warping()

void cs_user_mesh_warping ( void  )

Set options for cutting of warped faces.

◆ cs_user_model()

void cs_user_model ( void  )

Select physical model options, including user fields.

This function is called at the earliest stages of the data setup, so field ids are not available yet.

◆ cs_user_numbering()

void cs_user_numbering ( void  )

Define advanced mesh numbering options.

◆ cs_user_parallel_io()

void cs_user_parallel_io ( void  )

Define parallel IO settings.

◆ cs_user_paramedmem_define_couplings()

void cs_user_paramedmem_define_couplings ( void  )

Define ParaMEDMEM coupling(s)

◆ cs_user_paramedmem_define_fields()

void cs_user_paramedmem_define_fields ( void  )

Define fields to couple with ParaMEDMEM.

◆ cs_user_paramedmem_define_meshes()

void cs_user_paramedmem_define_meshes ( void  )

Define coupled meshes.

◆ cs_user_parameters()

void cs_user_parameters ( cs_domain_t domain)

Define or modify general numerical and physical user parameters.

At the calling point of this function, most model-related most variables and other fields have been defined, so specific settings related to those fields may be set here.

At this stage, the mesh is not built or read yet, so associated data such as field values are not accessible yet, though pending mesh operations and some fields may have been defined.

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_partition()

void cs_user_partition ( void  )

Define advanced partitioning options.

◆ cs_user_periodicity()

void cs_user_periodicity ( void  )

Define periodic faces.

◆ cs_user_physical_properties()

void cs_user_physical_properties ( cs_domain_t domain)

This function is called each time step to define physical properties.

Parameters
[in,out]domainpointer to a cs_domain_t structure

This function is called each time step to define physical properties.

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_physical_properties_h_to_t()

void cs_user_physical_properties_h_to_t ( cs_domain_t domain,
const cs_zone_t z,
bool  z_local,
const cs_real_t  h[restrict],
cs_real_t  t[restrict] 
)

User definition of enthalpy to temperature conversion.

This allows overwriting the solver defaults if necessary.

This function may be called on a per-zone basis, so as to allow different conversion relations in zones representing solids or different fluids.

Parameters
[in,out]domainpointer to a cs_domain_t structure
[in]zzone (volume or boundary) applying to current call
[in]z_localif true, h and t arrays are defined in a compact (contiguous) manner for this zone only; if false, h and t are defined on the zone's parent location (usually all cells or boundary faces)
[in]henthalpy values
[in,out]ttemperature values

◆ cs_user_physical_properties_t_to_h()

void cs_user_physical_properties_t_to_h ( cs_domain_t domain,
const cs_zone_t z,
bool  z_local,
const cs_real_t  t[restrict],
cs_real_t  h[restrict] 
)

User definition of temperature to enthalpy conversion.

This allows overwriting the solver defaults if necessary.

This function may be called on a per-zone basis, so as to allow different conversion relations in zones representing solids or different fluids.

Parameters
[in,out]domainpointer to a cs_domain_t structure
[in]zzone (volume or boundary) applying to current call
[in]z_localif true, h and t arrays are defined in a compact (contiguous) manner for this zone only; if false, h and t are defined on the zone's parent location (usually all cells or boundary faces)
[in]htemperature values
[in,out]tenthalpy values

◆ cs_user_physical_properties_td_pressure()

void cs_user_physical_properties_td_pressure ( cs_real_t td_p)

User function to define a custom law for the thermodynamic pressure.

Allows to define a custom law for the constant uniform thermodynamic pressure (whenn cs_velocity_pressure_model_t::idilat = 3 or cs_fluid_properties_t::ipthrm = 1).

The density is then updated (in pthrbm.f90) as:

\[\rho^{n+1} =\rho^{n} \cdot \frac{P_{th}^{n+1}}{P_{th}^{n}}\]

.

Parameters
[in,out]td_pUpdated value of the thermodynamic pressure

◆ cs_user_physical_properties_turb_viscosity()

void cs_user_physical_properties_turb_viscosity ( cs_domain_t domain)

User modification of the turbulence viscosity.

Parameters
[in,out]domainpointer to a cs_domain_t structure

Turbulent viscosity $ \mu_T $ (kg/(m s)) can be modified. You can access the field by its name.

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_porosity()

void cs_user_porosity ( cs_domain_t domain)

Compute the porosity (volume factor $ \epsilon $ when the porosity model is activated. (cs_glob_porous_model > 0).

This function is called at the begin of the simulation only.

Parameters
[in,out]domainpointer to a cs_domain_t structure

Compute the porosity (volume factor $ \epsilon $ when the porosity model is activated. (cs_glob_porous_model > 0).

This function is called at the begin of the simulation only.

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_postprocess_activate()

void cs_user_postprocess_activate ( int  nt_max_abs,
int  nt_cur_abs,
double  t_cur_abs 
)

Override default frequency or calculation end based output.

This allows fine-grained control of activation or deactivation,

Parameters
[in]nt_max_absmaximum time step number
[in]nt_cur_abscurrent time step number
[in]t_cur_absabsolute time at the current time step

◆ cs_user_postprocess_meshes()

void cs_user_postprocess_meshes ( void  )

Define post-processing meshes.

The main post-processing meshes may be configured, and additional post-processing meshes may be defined as a subset of the main mesh's cells or faces (both interior and boundary).

◆ cs_user_postprocess_probes()

void cs_user_postprocess_probes ( void  )

Define monitoring probes and profiles.

Profiles are defined as sets of probes.

◆ cs_user_postprocess_values()

void cs_user_postprocess_values ( const char *  mesh_name,
int  mesh_id,
int  cat_id,
cs_probe_set_t probes,
cs_lnum_t  n_cells,
cs_lnum_t  n_i_faces,
cs_lnum_t  n_b_faces,
cs_lnum_t  n_vertices,
const cs_lnum_t  cell_list[],
const cs_lnum_t  i_face_list[],
const cs_lnum_t  b_face_list[],
const cs_lnum_t  vertex_list[],
const cs_time_step_t ts 
)

User function for output of values on a post-processing mesh.

Parameters
[in]mesh_namename of the output mesh for the current call
[in]mesh_idid of the output mesh for the current call
[in]cat_idcategory id of the output mesh for the current call
[in]probespointer to associated probe set structure if the mesh is a probe set, NULL otherwise
[in]n_cellslocal number of cells of post_mesh
[in]n_i_faceslocal number of interior faces of post_mesh
[in]n_b_faceslocal number of boundary faces of post_mesh
[in]n_verticeslocal number of vertices faces of post_mesh
[in]cell_listlist of cells (0 to n-1) of post-processing mesh
[in]i_face_listlist of interior faces (0 to n-1) of post-processing mesh
[in]b_face_listlist of boundary faces (0 to n-1) of post-processing mesh
[in]vertex_listlist of vertices (0 to n-1) of post-processing mesh
[in]tstime step status structure, or NULL

◆ cs_user_postprocess_writers()

void cs_user_postprocess_writers ( void  )

Define post-processing writers.

The default output format and frequency may be configured, and additional post-processing writers allowing outputs in different formats or with different format options and output frequency than the main writer may be defined.

◆ cs_user_rad_transfer_absorption()

void cs_user_rad_transfer_absorption ( const int  bc_type[],
cs_real_t  ck[] 
)

Absorption coefficient for radiative module.

It is necessary to define the value of the fluid's absorption coefficient Ck.

This value is defined automatically for specific physical models, such as gas and coal combustion, so this function is not used by these models.

For a transparent medium, the coefficient should be set to 0.

In the case of the P-1 model, we check that the optical length is at least of the order of 1.

Parameters
[in]bc_typeboundary face types
[out]ckmedium's absorption coefficient (zero if transparent)

◆ cs_user_rad_transfer_net_flux()

void cs_user_rad_transfer_net_flux ( const int  bc_type[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
const cs_real_t  twall[],
const cs_real_t  qincid[],
const cs_real_t  xlam[],
const cs_real_t  epa[],
const cs_real_t  eps[],
const cs_real_t  ck[],
cs_real_t  net_flux[] 
)

Compute the net radiation flux.

The density of net radiation flux must be calculated consistently with the boundary conditions of the intensity. The density of net flux is the balance between the radiative emiting part of a boundary face (and not the reflecting one) and the radiative absorbing part.

Parameters
[in]bc_typeboundary face types
[in]coefapboundary condition work array for the radiance (explicit part)
[in]coefbpboundary condition work array for the radiance (implicit part)
[in]cofafpboundary condition work array for the diffusion of the radiance (explicit part)
[in]cofbfpboundary condition work array for the diffusion of the radiance (implicit part)
[in]twallinside current wall temperature (K)
[in]qincidradiative incident flux (W/m2)
[in]xlamconductivity (W/m/K)
[in]epathickness (m)
[in]epsemissivity (>0)
[in]ckabsorption coefficient
[out]net_fluxnet flux (W/m2)

◆ cs_user_radiative_transfer_bcs()

void cs_user_radiative_transfer_bcs ( cs_domain_t domain,
const int  bc_type[],
int  isothp[],
cs_real_t tmin,
cs_real_t tmax,
cs_real_t tx,
const cs_real_t  dt[],
const cs_real_t  thwall[],
const cs_real_t  qincid[],
cs_real_t  hfcnvp[],
cs_real_t  flcnvp[],
cs_real_t  xlamp[],
cs_real_t  epap[],
cs_real_t  epsp[],
cs_real_t  textp[] 
)

User definition of radiative transfer boundary conditions.

See Examples of data settings for radiative transfers for examples.

Warning
the temperature unit here is the Kelvin

For each boundary face face_id, a specific output (logging and postprocessing) class id may be assigned. This allows realizing balance sheets by treating them separately for each zone. By default, the output class id is set to the general (input) zone id associated to a face.

To access output class ids (both for reading and modifying), use the cs_boundary_zone_face_class_id function. The zone id values are arbitrarily chosen by the user, but must be positive integers; very high numbers may also lead to higher memory consumption.

Wall characteristics

The following face characteristics must be set:

  • isothp(face_id) boundary face type
    • CS_BOUNDARY_RAD_WALL_GRAY: Gray wall with temperature based on fluid BCs
    • CS_BOUNDARY_RAD_WALL_GRAY_EXTERIOR_T: Gray wall with fixed outside temperature
    • CS_BOUNDARY_RAD_WALL_REFL_EXTERIOR_T: Reflecting wall with fixed outside temperature
    • CS_BOUNDARY_RAD_WALL_GRAY_COND_FLUX: Gray wall with fixed conduction flux
    • CS_BOUNDARY_RAD_WALL_REFL_COND_FLUX: Reflecting wall with fixed conduction flux

Depending on the value of isothp, other values may also need to be set:

  • rcodcl = conduction flux
  • epsp = emissivity
  • xlamp = conductivity (W/m/K)
  • epap = thickness (m)
  • textp = outside temperature (K)
Parameters
[in,out]domainpointer to a cs_domain_t structure
[in]bc_typeboundary face types
[in]isothpboundary face type for radiative transfer
[out]tminmin allowed value of the wall temperature
[out]tmaxmax allowed value of the wall temperature
[in]txrelaxation coefficient (0 < tx < 1)
[in]dttime step (per cell)
[in]thwallinside current wall temperature (K)
[in]qincidradiative incident flux (W/m2)
[in]hfcnvpconvective exchange coefficient (W/m2/K)
[in]flcnvpconvective flux (W/m2)
[out]xlampconductivity (W/m/K)
[out]epapthickness (m)
[out]epspemissivity (>0)
[out]textpoutside temperature (K)

◆ cs_user_radiative_transfer_parameters()

void cs_user_radiative_transfer_parameters ( void  )

User function for input of radiative transfer module options.

◆ cs_user_saturne_coupling()

void cs_user_saturne_coupling ( void  )

Define couplings with other instances of code_saturne.

This is done by calling the cs_sat_coupling_define function for each coupling to add.

◆ cs_user_scaling_elec()

void cs_user_scaling_elec ( const cs_mesh_t mesh,
const cs_mesh_quantities_t mesh_quantities,
cs_real_t dt 
)

Define scaling parameter for electric model.

Define scaling parameter for electric model.

Parameters
[in]meshpointer to a cs_mesh_t structure
[in,out]mesh_quantitiespointer to a cs_mesh_quantities_t structure
[in]dtpointer to a cs_real_t

These options allow defining the time step synchronization policy, as well as a time step multiplier.

◆ cs_user_solver()

void cs_user_solver ( const cs_mesh_t mesh,
const cs_mesh_quantities_t mesh_quantities 
)

Main call to user solver.

Parameters
[in]meshpointer to a cs_mesh_t structure
[in,out]mesh_quantitiespointer to a cs_mesh_quantities_t structure

◆ cs_user_solver_set()

int cs_user_solver_set ( void  )

Set user solver.

Returns
1 if user solver is called, 0 otherwise

◆ cs_user_source_terms()

void cs_user_source_terms ( cs_domain_t domain,
int  f_id,
cs_real_t st_exp,
cs_real_t st_imp 
)

Additional user-defined source terms for variable equations (momentum, scalars, turbulence...).

Parameters
[in,out]domainpointer to a cs_domain_t structure
[in]f_idfield id of the variable
[out]st_expexplicit source term
[out]st_impimplicit part of the source term

This function is called at each time step, for each relevant field. It is therefore necessary to test the value of the field id or name to separate the treatments of the different variables.

The additional source term is decomposed into an explicit part (st_exp) and an implicit part (st_imp) that must be provided here. The resulting equation solved by the code for a scalar f is:

\[ \rho*volume*\frac{df}{dt} + .... = st\_imp*f + st\_exp \]

Note that st_exp and st_imp are defined after the Finite Volume integration over the cells, so they include the "volume" term. More precisely:

  • st_exp is expressed in kg.[var]/s, where [var] is the unit of the variable. Its dimension is the one of the variable (3 for vectors)
  • st_imp is expressed in kg/s. Its dimension is 1 for scalars, 3x3 for vectors.

The st_exp and st_imp arrays are already initialized to 0 (or a value defined through the GUI or defined by a model) before entering the function. It is generally not useful to reset them here.

For stability reasons, code_saturne will not add -st_imp directly to the diagonal of the matrix, but Max(-st_imp,0). This way, the st_imp term is treated implicitely only if it strengthens the diagonal of the matrix. However, when using the second-order in time scheme, this limitation cannot be done anymore and -st_imp is added directly. The user should therefore check for the negativity of st_imp.

When using the second-order in time scheme, one should supply:

  • st_exp at time n
  • st_imp at time n+1/2
Warning

If the variable is a temperature, the resulting equation solved is:

rho*Cp*volume*dT/dt + .... = st_imp*T + st_exp

Note that st_exp and st_imp are defined after the Finite Volume integration over the cells, so they include the "volume" term. More precisely:

  • st_exp is expressed in W
  • st_imp is expressed in W/K
Steep source terms

In case of a complex, non-linear source term, say F(f), for variable f, the easiest method is to implement the source term explicitly.

df/dt = .... + F(f(n)) where f(n) is the value of f at time tn, the beginning of the time step.

This yields: st_exp = volume*F(f(n)) st_imp = 0

However, if the source term is potentially steep, this fully explicit method will probably generate instabilities. It is therefore wiser to partially implicit the term by writing:

df/dt = .... + dF/df*f(n+1) - dF/df*f(n) + F(f(n))

This yields: st_exp = volume*( F(f(n)) - dF/df*f(n) ) st_imp = volume*dF/df

Parameters
[in,out]domainpointer to a cs_domain_t structure
[in]f_idfield id of the variable
[out]st_expexplicit source term
[out]st_impimplicit part of the source term

◆ cs_user_syrthes_coupling()

void cs_user_syrthes_coupling ( void  )

Define couplings with SYRTHES code.

This is done by calling the cs_syr_coupling_define function for each coupling to add.

◆ cs_user_syrthes_coupling_volume_h()

void cs_user_syrthes_coupling_volume_h ( int  coupling_id,
const char *  syrthes_name,
cs_lnum_t  n_elts,
const cs_lnum_t  elt_ids[],
cs_real_t  h_vol[] 
)

Compute a volume exchange coefficient for SYRTHES couplings.

Parameters
[in]coupling_idSyrthes coupling id
[in]syrthes_namename of associated Syrthes instance
[in]n_eltsnumber of associated cells
[in]elt_idsassociated cell ids
[out]h_volassociated exchange coefficient (size: n_elts)

◆ cs_user_time_moments()

void cs_user_time_moments ( void  )

Define time moments.

This function is called at the setup stage, once user and most model-based fields are defined, and before fine control of field output options is defined.

◆ cs_user_time_table()

void cs_user_time_table ( void  )

Define time tables using C API. This function is called at the begin of the simulation only.

Parameters
[in,out]domainpointer to a cs_domain_t structure

◆ cs_user_turbomachinery()

void cs_user_turbomachinery ( void  )

Define rotor/stator model.

◆ cs_user_turbomachinery_rotor()

void cs_user_turbomachinery_rotor ( void  )

Define rotor axes, associated cells, and rotor/stator faces.

◆ cs_user_turbomachinery_set_rotation_velocity()

void cs_user_turbomachinery_set_rotation_velocity ( void  )

Define rotation velocity of rotor.

◆ cs_user_wall_condensation()

void cs_user_wall_condensation ( int  nvar,
int  nscal,
int  iappel 
)

Source terms associated at the boundary faces and the neighboring cells with surface condensation.

This function fills the condensation source terms for each variable at the cell center associated to the boundary faces identifed in the mesh. The fluid exchange coefficient is computed with a empirical law to be imposed at the boundary face where the condensation phenomenon occurs.

This user function is called which allows the setting of $ \gamma_{\mbox{cond}} $ the condensation source term.

This function fills the condensation source term array gamma_cond adding to the following equations:

  • The equation for mass conservation:

    \[ D\frac{rho}{dt} + divs \left( \rho \vect{u}^n\right) = \Gamma _{cond} \]

  • The equation for a variable $\Phi $:

    \[ D\frac{\phi}{dt} = ... + \Gamma _{cond}*(\Phi _i - \Phi) \]

discretized as below:

\[ \rho*\dfrac{\Phi^{n+1}-\Phi^{n}}/dt = ... + \Gamma _{cond}*(\Phi _i - \Phi^{n+1}) \]

Remarks
  • $ \Phi _i $ is the value of $ \Phi $ associated to the injected condensation rate.

    With 2 options are available:

    • the condensation rate is injected with the local value of variable $ \Phi = \Phi ^{n+1}$ in this case the $ \Phi $ variable is not modified.
    • the condensation rate is injected with a specific value for $ \Phi = \Phi _i $ the specified value given by the user.

Usage

The three stages in the code where this User subroutine is called (with

iappel = 1, 2 and 3

)

iappel = 1
  • Calculation of the number of cells where a mass source term is imposed: ncesmp Called once at the beginning of the calculation
iappel = 2
  • Identification of the cells where a mass source term is imposed: array icesmp(ncesmp) Called once at the beginning of the calculation
iappel = 3
  • Calculation of the values of the mass source term Called at each time step

specific variables to define with is user subroutine

  • ifbpcd(ieltcd): identification of the faces where a condensation source term is imposed.
  • itypcd(ieltcd,ivar): type of treatment for variable ivar in the ieltcd cell with condensation source term.
    • itypcd = 0 – * injection of ivar at local value
    • itypcd = 1 – * injection of ivar at user specified value.
  • spcond(ielscd,ipr): value of the injection condensation rate gamma_cond (kg/m3/s) in the ieltcd cell with condensation source term.
  • spcond(ieltcd,ivar): specified value for variable ivar associated to the injected condensation in the ieltcd cell with a condensation source term except for ivar=ipr.
Remarks
  • For each face where a condensation source terms is imposed ielscd in [1;nfbpcd]), ifbpcd(ielscd) is the global index number of the corresponding face (ifbpcd(ieltcd) in [1;ncel]).
  • if itypcd(ieltcd,ivar)=0, spcond(ielpcd,ivar) is not used.
  • if spcond(ieltcd,ipr)<0, mass is removed from the system, therefore Code_Saturna automatically considers f_i=f^(n+1), whatever the values of itypcd or smacel specified by the user
Examples of settings for boundary condensation mass source terms
Examples are available here.
Parameters
[in]nvartotal number of variables
[in]nscaltotal number of scalars
[in]iappelindicates which at which stage the routine is

◆ cs_user_zones()

void cs_user_zones ( void  )

Define volume and surface zones.

See Element selection criteria for details on selection criteria.

◆ cs_wall_condensation_volume_exchange_surf_at_cells()

void cs_wall_condensation_volume_exchange_surf_at_cells ( cs_real_t surf)

Return condensing volume structures surface at each cell.

Parameters
[out]surfarray of volume structures surface at each cell
[out]surfarray of volume structure surfaces at each cell

◆ csinit()

void csinit ( const int *  irgpar,
const int *  nrgpar 
)

◆ findpt()

void findpt ( const cs_lnum_t ncelet,
const cs_lnum_t ncel,
const cs_real_t xyzcen,
const cs_real_t xx,
const cs_real_t yy,
const cs_real_t zz,
cs_lnum_t node,
int *  ndrang 
)

◆ initi1()

void initi1 ( void  )