7.0
general documentation
cs_user_boundary_conditions.f90 File Reference

User subroutine which fills boundary conditions arrays (icodcl, rcodcl) for unknown variables. More...

Functions/Subroutines

subroutine cs_f_user_boundary_conditions (nvar, nscal, icodcl, itrifb, itypfb, izfppp, dt, rcodcl)
 

Detailed Description

User subroutine which fills boundary conditions arrays (icodcl, rcodcl) for unknown variables.

See cs_user_boundary_conditions.f90 for examples.

Introduction

Here one defines boundary conditions on a per-face basis.

Boundary faces may be selected using the getfbr subroutine.

getfbr(string, nelts, lstelt)
  • string is a user-supplied character string containing selection criteria;
  • nelts is set by the subroutine. It is an integer value corresponding to the number of boundary faces verifying the selection criteria;
  • lstelt is set by the subroutine. It is an integer array of size nelts containing the list of boundary faces verifying the selection criteria.

string may contain:

  • references to colors (ex.: 1, 8, 26, ...)
  • references to groups (ex.: inlet, group1, ...)
  • geometric criteria (ex. x < 0.1, y >= 0.25, ...)

These criteria may be combined using logical operators (and,or) and parentheses.

Example
1 and(group2 or group3) and y < 1
will select boundary faces of color 1, belonging to groups 'group2' or 'group3' and with face center coordinate y less than 1.

Operators priority, from highest to lowest: '( )' > 'not' > 'and' > 'or' > 'xor'

Similarly, interior faces and cells can be identified using the getfac and getcel subroutines (respectively). Their syntax are identical to getfbr syntax.

For a more thorough description of the criteria syntax, see the user guide.

Boundary condition types

Boundary conditions may be assigned in two ways.

For "standard" boundary conditions:

One defines a code in the itypfb array (of dimensions number of boundary faces). This code will then be used by a non-user subroutine to assign the following conditions. The available codes are:

  • ientre: Inlet
  • isolib: Free outlet
  • isymet: Symmetry
  • iparoi: Wall (smooth)
  • iparug: Rough wall

These integers are defined elsewhere (in paramx.f90 module). Their value is greater than or equal to 1 and less than or equal to ntypmx (value fixed in paramx.h)

In addition, some values must be defined:

  • Inlet (more precisely, inlet/outlet with prescribed flow, as the flow may be prescribed as an outflow):
    • Dirichlet conditions on variables other than pressure are mandatory if the flow is incoming, optional if the flow is outgoing (the code assigns zero flux if no Dirichlet is specified); thus, at face ifac, for the variable ivar: rcodcl(ifac, ivar, 1)
  • Smooth wall: (= impermeable solid, with smooth friction)
    • Velocity value for sliding wall if applicable:
      • rcodcl(ifac, iu, 1) = fluid velocity in the x direction
      • rcodcl(ifac, iv, 1) = fluid velocity in the y direction
      • rcodcl(ifac, iw, 1) = fluid velocity in the z direction
    • Specific code and prescribed temperature value at wall if applicable:
      • icodcl(ifac, ivar) = 5
      • rcodcl(ifac, ivar, 1) = prescribed temperature
    • Specific code and prescribed flux value at wall if applicable:
      • icodcl(ifac, ivar) = 3
      • rcodcl(ifac, ivar, 3) = prescribed flux
    Note that the default condition for scalars (other than k and epsilon) is homogeneous Neumann.
  • Rough wall: (= impermeable solid, with rough friction)
    • Velocity value for sliding wall if applicable:
      • rcodcl(ifac, iu, 1) = fluid velocity in the x direction
      • rcodcl(ifac, iv, 1) = fluid velocity in the y direction
      • rcodcl(ifac, iw, 1) = fluid velocity in the z direction
    • Value of the dynamic roughness height to specify in
      • boundary_roughness field
    • Value of the scalar roughness height (if required) to specify in
      • boundary_thermal_roughness
    • Specific code and prescribed temperature value at wall if applicable:
      • icodcl(ifac, ivar) = 6
      • rcodcl(ifac, ivar, 1) = prescribed temperature
    • Specific code and prescribed flux value at rough wall, if applicable:
      • icodcl(ifac, ivar) = 3
      • rcodcl(ifac, ivar, 3) = prescribed flux
    Note that the default condition for scalars (other than k and epsilon) is homogeneous Neumann.
  • Symmetry (= slip wall):
    • Nothing to specify
  • Free outlet (more precisely free inlet/outlet with prescribed pressure)
    • Nothing to prescribe for pressure and velocity. For scalars and turbulent values, a Dirichlet value may optionally be specified. The behavior is as follows:
      • pressure is always handled as a Dirichlet condition
      • if the mass flow is inflowing: one retains the velocity at infinity Dirichlet condition for scalars and turbulent values (or zero flux if the user has not specified a Dirichlet value)
      • if the mass flow is outflowing: one prescribes zero flux on the velocity, the scalars, and turbulent values
    Note that the pressure will be reset to p0 on the first free outlet face found.

For "non-standard" conditions:

Other than (inlet, free outlet, wall, symmetry), one defines

  • on one hand, for each face:
    • an admissible itypfb value (i.e. greater than or equal to 1 and less than or equal to ntypmx; see its value in paramx.h). The values predefined in paramx.h: ientre, isolib, isymet, iparoi, iparug are in this range, and it is preferable not to assign one of these integers to itypfb randomly or in an inconsiderate manner. To avoid this, one may use iindef if one wish to avoid checking values in paramx.h. iindef is an admissible value to which no predefined boundary condition is attached. Note that the itypfb array is reinitialized at each time step to the non-admissible value of 0. If one forgets to modify itypfb for a given face, the code will stop.
  • and on the other hand, for each face and each variable:
    • a code
      • icodcl(ifac, ivar)
    • three real values
      • rcodcl(ifac, ivar, 1)
      • rcodcl(ifac, ivar, 2)
      • rcodcl(ifac, ivar, 3)

The value of icodcl is taken from the following:

  • 1: Dirichlet (usable for any variable)
  • 3: Neumann (usable for any variable)
  • 4: Symmetry (usable only for the velocity and components of the Rij tensor)
  • 5: Smooth wall (usable for any variable except for pressure)
  • 6: Rough wall (usable for any variable except for pressure)
  • 9: Free outlet (usable only for velocity)
  • 13: Dirichlet for the advection operator and Neumann for the diffusion operator

The values of the 3 rcodcl components are:

  • rcodcl(ifac, ivar, 1):
    • Dirichlet for the variable if icodcl(ifac, ivar) = 1 or 13
    • Wall value (sliding velocity, temp) if icodcl(ifac, ivar) = 5
    The dimension of rcodcl(ifac, ivar, 1) is that of the resolved variable, for instance:
    • U (velocity in m/s),
    • T (temperature in degrees)
    • H (enthalpy in J/kg)
    • F (passive scalar in -)
  • rcodcl(ifac, ivar, 2): "exterior" exchange coefficient (between the prescribed value and the value at the domain boundary) rinfin = infinite by default
    • For velocities U, in kg/(m2 s): rcodcl(ifac, ivar, 2) = (viscl+visct) / d
    • For the pressure P, in s/m: rcodcl(ifac, ivar, 2) = dt / d
    • For temperatures T, in Watt/(m2 degres): rcodcl(ifac, ivar, 2) = Cp*(viscls+visct/turb_schmidt) / d
    • For enthalpies H, in kg /(m2 s): rcodcl(ifac, ivar, 2) = (viscls+visct/turb_schmidt) / d
    • For other scalars F in: rcodcl(ifac, ivar, 2) = (viscls+visct/turb_schmidt) / d (d has the dimension of a distance in m)
  • rcodcl(ifac, ivar, 3) if icodcl(ifac, ivar) = 3 or 13: Flux density (< 0 if gain, n outwards-facing normal)
    • For velocities U, in kg/(m s2) = J: rcodcl(ifac, ivar, 3) = -(viscl+visct) * (grad U).n
    • For pressure P, in kg/(m2 s): rcodcl(ifac, ivar, 3) = -dt * (grad P).n
    • For temperatures T, in Watt/m2: rcodcl(ifac, ivar, 3) = -Cp*(viscls+visct/turb_schmidt) * (grad T).n
    • For enthalpies H, in Watt/m2: rcodcl(ifac, ivar, 3) = -(viscls+visct/turb_schmidt) * (grad H).n
    • For other scalars F in: rcodcl(ifac, ivar, 3) = -(viscls+visct/turb_schmidt) * (grad F).n
  • rcodcl(ifac, ivar, 3) if icodcl(ifac, ivar) = 6: Roughness length scale for the rough wall law
    • For velocities U, dynamic roughness boundary_roughness field contains the dynamic roughness length scale
    • For other scalars, thermal roughness boundary_thermal_roughness field contains the scalar roughness length scale

Note that if the user assigns a value to itypfb equal to ientre, isolib, isymet, iparoi, or iparug and does not modify icodcl (zero value by default), itypfb will define the boundary condition type.

To the contrary, if the user prescribes icodcl(ifac, ivar) (nonzero), the values assigned to rcodcl will be used for the considered face and variable (if rcodcl values are not set, the default values will be used for the face and variable, so:

  • rcodcl(ifac, ivar, 1) = 0.d0
  • rcodcl(ifac, ivar, 2) = rinfin
  • rcodcl(ifac, ivar, 3) = 0.d0)

Especially, one may have for example:

  • set itypfb(ifac) = iparoi which prescribes default wall conditions for all variables at face ifac,
  • and define IN ADDITION for variable ivar on this face specific conditions by specifying icodcl(ifac, ivar) and the 3 rcodcl values.

The user may also assign to itypfb a value not equal to ientre, isolib, isymet, iparoi, iparug, iindef but greater than or equal to 1 and less than or equal to ntypmx (see values in param.h) to distinguish groups or colors in other subroutines which are specific to the case and in which itypfb is accessible. In this case though it will be necessary to prescribe boundary conditions by assigning values to icodcl and to the 3 rcodcl fields (as the value of itypfb will not be predefined in the code).

Boundary condition types for enthalpy

For enthalpy, prescribed values in rcodcl(ifac, ivar, 1) may be defined using the temperature instead. In this case, icodcl(ifac, ivar) must be replaced by -icodcl(ifac, ivar) to mark the face for automatic conversion.

Boundary condition types for enthalpy

For compressible flows, only predefined boundary conditions may be assigned among: iparoi, isymet, iesicf, isspcf, isopcf, iephcf, ieqhcf

  • iparoi : standard wall
  • isymet : standard symmetry
  • iesicf, isspcf, isopcf, iephcf, ieqhcf : inlet/outlet

For inlets/outlets, we can prescribe a value for turbulence and passive scalars in rcodcl(.,.,1) for the case in which the mass flux is incoming. If this is not done, a zero flux condition is applied.

  • iesicf: prescribed inlet/outlet (for example supersonic inlet) the user prescribes the velocity and all thermodynamic variables
  • isspcf: supersonic outlet the user does not prescribe anything
  • isopcf: subsonic outlet with prescribed pressure the user presribes the pressure
  • iephcf: mixed inlet with prescribed total pressure and enthalpy the user prescribes the total pressure and total enthalpy
  • ieqhcf: subsonic inlet with prescribed mass and enthalpy flow to be implemented

Consistency rules

A few consistency rules between icodcl codes for variables with non-standard boundary conditions:

  • Codes for velocity components must be identical
  • Codes for Rij components must be identical
  • If code (velocity or Rij) = 4 one must have code (velocity and Rij) = 4
  • If code (velocity or turbulence) = 5 one must have code (velocity and turbulence) = 5
  • If code (velocity or turbulence) = 6 one must have code (velocity and turbulence) = 6
  • If scalar code (except pressure or fluctuations) = 5 one must have velocity code = 5
  • If scalar code (except pressure or fluctuations) = 6 one must have velocity code = 6
Remarks
  • Caution: to prescribe a flux (nonzero) to Rij, the viscosity to take into account is viscl even if visct exists (visct=rho cmu k2/epsilon)
  • One have the ordering array for boundary faces from the previous time step (except for the fist one, where itrifb has not been set yet).
  • The array of boundary face types itypfb has been reset before entering the subroutine.

Cell values of some variables

Cell value field ids

  • Density: icrom
  • Dynamic molecular viscosity: iviscl
  • Turbulent viscosity: ivisct
  • Specific heat: icp
  • Diffusivity(lambda): field_get_key_int(ivarfl(isca(iscal)), & kivisl, ...)

Faces identification

  • Density: field id ibrom
  • Boundary mass flux (for convecting ivar): field id iflmab using field_get_key_int(ivarfl(ivar), kbmasf, iflmab)
  • For other values: take as an approximation the value in the adjacent cell i.e. as above with iel = ifabor(ifac).

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

Function/Subroutine Documentation

◆ cs_f_user_boundary_conditions()

subroutine cs_f_user_boundary_conditions ( integer  nvar,
integer  nscal,
integer, dimension(nfabor,nvar)  icodcl,
integer, dimension(nfabor)  itrifb,
integer, dimension(nfabor)  itypfb,
integer, dimension(nfabor)  izfppp,
double precision, dimension(ncelet)  dt,
double precision, dimension(nfabor,nvar,3)  rcodcl 
)
Parameters
[in]nvartotal number of variables
[in]nscaltotal number of scalars
[out]icodclboundary condition code:
  • 1 Dirichlet
  • 2 Radiative outlet
  • 3 Neumann
  • 4 sliding and $ \vect{u} \cdot \vect{n} = 0 $
  • 5 smooth wall and $ \vect{u} \cdot \vect{n} = 0 $
  • 6 rough wall and $ \vect{u} \cdot \vect{n} = 0 $
  • 9 free inlet/outlet (input mass flux blocked to 0)
  • 13 Dirichlet for the advection operator and Neumann for the diffusion operator
[in]itrifbindirection for boundary faces ordering
[in,out]itypfbboundary face types
[out]izfpppboundary face zone number
[in]dttime step (per cell)
[in,out]rcodclboundary condition values:
  • rcodcl(1) value of the dirichlet
  • rcodcl(2) value of the exterior exchange coefficient (infinite if no exchange)
  • rcodcl(3) value flux density (negative if gain) in w/m2
    1. for the velocity $ (\mu+\mu_T) \gradt \, \vect{u} \cdot \vect{n} $
    2. for the pressure $ \Delta t \grad P \cdot \vect{n} $
    3. for a scalar $ cp \left( K + \dfrac{K_T}{\sigma_T} \right) \grad T \cdot \vect{n} $