8.0
general documentation
Input of calculation parameters (Fortran subroutines in cs_user_parameters.f90)

Introduction

User subroutines for input of calculation parameters (Fortran modules). 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 routines are present in the file, each destined to defined specific parameters.

To modify the default value of parameters which do not appear in the examples provided, code should be placed as follows:

  • usipsu for numerical and physical options
  • usipes for input-output related options
  • usppmo for specific physics options
  • usipph for additional input of parameters
  • usati1 for calculation options for the atmospheric module
  • cs_user_combustion for calculation options for the combustion module
  • uscfx1 and uscfx2 for non-standard options for the compressible module
  • uscti1 for the definition of cooling tower model and exchange zones
  • user_darcy_ini1 for calculation options for the Darcy module

As a convention, "specific physics" defers to the following modules only: pulverized coal, gas combustion, electric arcs.

In addition, specific routines are provided for the definition of some "specific physics" options. These routines are described at the end of this file and will be activated when the corresponding option is selected in the usppmo routine.

General options (usipsu)

The following code block presents options available in the usipsu subroutine.

! Calculation options (optcal)
! ============================
! Error estimators for Navier-Stokes (non-frozen velocity field)
! We recommend running a calculation restart on a few time steps
! with the activation of the most interesting of those.
! (=2 to activate, =0 to deactivate).
iescal(iescor) = 2 ! div(rho u) -Gamma
iescal(iestot) = 2 ! resolution precision for the momentum

Postprocessing options (usipes)

The usipes routine can be found in the cs_user_parameters.f90 file, and is used for setting of postprocessing options, as it is called once all fields are defined.

Specific physic activation (usppmo)

The usppmo routine can be found in the cs_user_parameters.f90 file.

!===============================================================================
! 1. Choice for a specific physics
!===============================================================================
! --- cod3p: Diffusion flame with complete fast chemistry (3 points)
! ==========
! if = -1 module not activated
! if = 0 adiabatic model
! if = 1 extended model with enthalpy source term
ippmod(icod3p) = -1
! --- coebu: Eddy-Break Up pre-mixed flame
! ==========
! if = -1 module not activated
! if = 0 reference Spalding model
! (adiabatic, homogeneous mixture fraction)
! if = 1 extended model with enthalpy source term
! (homogeneous mixture fraction : perfect premix)
! if = 2 extended model with mixture fraction transport
! (adiabatic, no variance of mixture fraction)
! if = 3 extended model with enthalpy and mixture fraction transport
! (dilution, thermal losses, etc.)
ippmod(icoebu) = -1
! --- colwc: Libby-Williams pre-mixed flame
! ==========
! if = -1 module not activated
! if = 0 reference two-peak model with adiabatic condition
! if = 1 extended two-peak model with enthapy source terms
! if = 2 extended three-peak model, adiabatic
! if = 3 extended three-peak model with enthalpy source terms
! if = 4 extended four-peak model, adiabatic
! if = 5 extended four-peak model with enthalpy source terms
ippmod(icolwc) = -1
! --- Soot model
! =================
! if = -1 module not activated
! if = 0 constant soot yield
! if = 1 2 equations model of Moss et al.
isoot = 0
xsoot = 0.1d0 ! (if isoot = 0 and only if the soot yield is not
! defined in the thermochemistry data file)
rosoot = 2000.d0 ! kg/m3
! --- cfuel: Heavy fuel oil combustion
! ==========
! Progressive evaporation (temperature gap)
! Char residue
! Sulphur tracking
! if = -1 module not activated
! if = 0 module activated
ippmod(icfuel) = -1
! --- coal :
! ==========
!
! Pulverized coal combustion
! Description of granulometry
! Assumption of diffusion flame around particles
! (extension of 3-point fast chemistry "D3P")
! Between a mixture of gaseous fuels (volatiles matters, CO from char
! oxydation)
! and a mixture of oxidisers (air and water vapor)
! Enthalpy for both mix and solid phase are solved
!
! if = -1 module not activated
! if = 0 module activated
! if = 1 with drying
ippmod(iccoal) = -1
! Activate the drift: 0 (no activation),
! 1 (transported particle velocity)
! 2 (limit drop particle velocity)
i_comb_drift = 1
! --- cpl3c: Pulverized coal with Lagrangian reciprocal approach
! ==========
! Not recently tested... at least outdated, may be obsolete
! if = -1 module not activated
! if = 0 module activated
! if = 1 with drying (NOT functional)
ippmod(icpl3c) = -1
! --- eljou: Joule effect
! ==========
! if = -1 module not activated
! if = 1 real potential
! if = 2 complex potential
! if = 3 real potential + Transfo
! if = 4 complex potential + Transfo
ippmod(ieljou) = -1
! --- elarc: Electric arcs
! ==========
! if = -1 module not activated
! if = 1 electric potential
! if = 2 electric potential and vector potential (hence 3D modelling)
ippmod(ielarc) = -1
! --- aeros: Cooling towers
! ==========
! if = -1 module not activated
! if = 0 no model (NOT functional)
! if = 1 Poppe's model
! if = 2 Merkel's model
ippmod(iaeros) = -1
! Radiative transfer module (iirayo)
!--------------------------
! if = 0: not activated (Default)
! if = 1: DOM
! if = 2: approximation P1 method
iirayo = 1
! --- richards model
! ==========
! if = -1 module not activated
! if = 1 module activated
ippmod(idarcy) = -1
!===============================================================================
! 2. Specific options related to herebefore modules
!===============================================================================
! These options are defined here at the moment, this might change in the future
! --- Enthalpy-Temperature conversion law (for gas combustion modelling)
! if = 0 user-specified
! if = 1 tabulated by JANAF (default)
indjon = 1
!===============================================================================
! 2. Data file related to modules above
!===============================================================================
! Combustion
if ( ippmod(icod3p).ge.0 &
.or. ippmod(icoebu).ge.0 .or. ippmod(icolwc).ge.0) then
if (indjon.eq.1) then
ficfpp = 'dp_C3P'
else
ficfpp = 'dp_C3PSJ'
endif
endif
! Fuel combustion
if (ippmod(icfuel).ge.0) then
ficfpp = 'dp_FUE'
endif
! Specific condensation modelling (used with ippmod(igmix) >= 0)
! wall condensation model
! -1: not activated
! 0: condensation source terms
icondb = -1
! internal condensation model
! -1: not activated
! 0: condensation source terms with metal structures
icondv = -1

Additional input of parameters (usipph)

The usipph routine can be found in the cs_user_parameters.f90 file.

Calculation options for the atmospheric module (usati1)

The usati1 routine can be found in the cs_user_parameters.f90 file.

! -----------------------------------------------------------------------------
! Atmospheric imbrication on large scale meteo (atimbr module)
! -----------------------------------------------------------------------------
!
! --------------------------------------------------------------
! activation flag
! --------------------------------------------------------------
imbrication_flag = .false.
imbrication_verbose = .false.
! ------------------------------------------------------------------------------
! flags for activating the cressman interpolation for the boundary conditions
! ------------------------------------------------------------------------------
cressman_u = .true.
cressman_v = .true.
cressman_tke = .true.
cressman_eps = .true.
cressman_theta = .true.
cressman_qw = .true.
cressman_nc = .true.
! --------------------------------------------------------------
! numerical parameters for the cressman interpolation formulas
! --------------------------------------------------------------
horizontal_influence_radius = 8500.d0
vertical_influence_radius = 100.d0
! --------------------------------------------------------------
! ifilechemistry: choice to read (=1,2,3,4, according to the scheme)
! or not (0) a concentration profile file
! if ichemistry>0 ifilechemistry is automaticaly set to ichemistry
ifilechemistry = 0
! Change the name of the chemistry concentration profile
call atmo_set_chem_conc_file_name('chem_01_01_2001')
! Change the name of the aerosol concentration profile
call atmo_set_aero_conc_file_name('aero_01_01_2001')
! isepchemistry: splitted (=1) or semi-coupled (=2, pu-sun)
! resolution of chemistry.
! Splitted (=1) mandatory for aerosols.
! Semi-coupled (=2) by default.
isepchemistry = 1
! dtchemmax: maximal time step (s) for chemistry resolution
dtchemmax = 10.0d0
! computation / storage of downward and upward infrared radiative fluxes
irdu = 1

Calculation options for the combustion module (cs_user_combustion)

The cs_user_combustion routine can be found in the cs_user_parameters.f90 file.

!===============================================================================
! 1. Additional Calculation Options
!===============================================================================
! --- Kinetic model for CO <=> CO2
! if = 0 unused (maximal conversion in turbulent model)
! if = 1 transport of CO2 mass fraction
! if = 2 transport of CO mass fraction (coal and fuel only)
ieqco2 = 0
! --- Density Relaxation
! RHO(n+1) = SRROM * RHO(n) + (1-SRROM) * RHO(n+1)
srrom = 0.8d0
!===============================================================================
! 2. Physical Constants
!===============================================================================
! diftl0: Dynamic Diffusion Coefficient (kg/(m s))
diftl0 = 4.25d-5
! -----------------------------------------------------------------------------
! 2.1 For 3 points combusution model ONLY
! -----------------------------------------------------------------------------
! Reference temperature for fuel and oxydant (K)
tinfue = 436.d0
tinoxy = 353.d0
! -----------------------------------------------------------------------------
! 2.2 For EBU-model ONLY
! -----------------------------------------------------------------------------
! cebu: EBU-model constant
cebu = 2.5d0
! -----------------------------------------------------------------------------
! 2.3 For Libby-Williams model ONLY
! -----------------------------------------------------------------------------
! Reference velocity
vref = 60.d0
! Reference length scale
lref = 0.1d0
! Activation Temperature
ta = 0.2d5
! Cross-over Temperature (combustion of propane)
tstar= 0.12d4

Non-standard options for the compressible module (uscfx1)

The uscfx1 routine can be found in the cs_user_parameters.f90 file.

!===============================================================================
! 1. Properties options
!===============================================================================
! --> Molecular thermal conductivity
! constant : ifcvsl = -1
! variable : ifcvsl = 0
ifcvsl = -1
call field_set_key_int(ivarfl(isca(itempk)), kivisl, ifcvsl)
! --> Volumetric molecular viscosity
! iviscv = -1 : uniform in space and constant in time
! = 0 : variable in space and time
iviscv = -1

Non-standard options for the compressible module (uscfx2)

The uscfx2 routine can be found in the cs_user_parameters.f90 file.

!===============================================================================
! 1. Physical properties
!===============================================================================
! --> Molecular viscosity
! constant : ivivar = 0
! variable : ivivar = 1
ivivar = 0
! --> Reference molecular thermal conductivity
! visls0 = lambda0 (molecular thermal conductivity, W/(m K))
! Reference conductivity must be strictly positive
! (set a realistic value here even if conductivity is variable)
call field_set_key_double(ivarfl(isca(itempk)), kvisl0, 3.d-2)
! If the molecular thermal conductivity is variable, its values
! must be provided in the user subroutine 'cs_user_physical_properties'
! --> Volumetric molecular viscosity
! Reference volumetric molecular viscosity
! viscv0 = kappa0 (volumetric molecular viscosity, kg/(m s))
viscv0 = 0.d0
! If the volumetric molecular viscosity is variable, its values
! must be provided in the user subroutine 'cs_user_physical_properties'
! --> Molar mass of the gas (kg/mol)
! For example with dry air, xmasml is around 28.8d-3 kg/mol
xmasmr = 0.028966
! --> Hydrostatic equilibrium
! Specify if the hydrostatic equilibrium must be accounted for
! (yes = 1 , no = 0)
icfgrp = 1

Definition of cooling tower model and exchange zones (uscti1)

The uscti1 routine can be found in the cs_user_parameters.f90 file.

Calculation options for the Darcy module (user_darcy_ini1)

The user_darcy_ini1 routine can be found in the cs_user_parameters.f90 file.

darcy_anisotropic_permeability = 0 ! permeability : 0 isotrop, 1 anisotrop
darcy_anisotropic_dispersion = 0 ! dispersion : 0 isotrop, 1 anisotrop
darcy_unsteady = 0 ! 0 steady flow, 1 unsteady flow
darcy_convergence_criterion = 0 ! convergence criterion of Newton scheme:
! 0, over pressure, 1, over velocity