My Project
programmer's documentation
|
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) |
User subroutine which fills boundary conditions arrays (icodcl
, rcodcl
) for unknown variables.
See cs_user_boundary_conditions.f90 for examples.
Here one defines boundary conditions on a per-face basis.
Boundary faces may be selected using the getfbr subroutine.
string may contain:
These criteria may be combined using logical operators (and
,or
) and parentheses.
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 conditions may be assigned in two ways.
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:
Inletisolib:
Free outletisymet:
Symmetryiparoi:
Wall (smooth)iparug:
Rough wallThese 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:
ifac
, for the variable ivar:
rcodcl(ifac, ivar, 1)
rcodcl(ifac, iu, 1)
= fluid velocity in the x directionrcodcl(ifac, iv, 1)
= fluid velocity in the y directionrcodcl(ifac, iw, 1)
= fluid velocity in the z directionicodcl(ifac, ivar)
= 5rcodcl(ifac, ivar, 1)
= prescribed temperatureicodcl(ifac, ivar)
= 3rcodcl(ifac, ivar, 3)
= prescribed fluxrcodcl(ifac, iu, 1)
= fluid velocity in the x directionrcodcl(ifac, iv, 1)
= fluid velocity in the y directionrcodcl(ifac, iw, 1)
= fluid velocity in the z directionrcodcl(ifac, iu, 3)
rcodcl(ifac, iv, 3)
(values for iw are not used)icodcl(ifac, ivar)
= 6rcodcl(ifac, ivar, 1)
= prescribed temperatureicodcl(ifac, ivar)
= 3rcodcl(ifac, ivar, 3)
= prescribed fluxOther than (inlet, free outlet, wall, symmetry), one defines
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.icodcl(ifac, ivar)
rcodcl(ifac, ivar, 1)
rcodcl(ifac, ivar, 2)
rcodcl(ifac, ivar, 3)
The value of icodcl
is taken from the following:
The values of the 3 rcodcl
components are:
rcodcl(ifac, ivar, 1)
:icodcl(ifac, ivar)
= 1 or 13icodcl(ifac, ivar)
= 5rcodcl(ifac, ivar, 1)
is that of the resolved variable, for instance:rcodcl(ifac, ivar, 2)
: "exterior" exchange coefficient (between the prescribed value and the value at the domain boundary) rinfin = infinite by defaultrcodcl(ifac, ivar, 2)
= (viscl+visct) / drcodcl(ifac, ivar, 2)
= dt / drcodcl(ifac, ivar, 2)
= Cp*(viscls+visct/turb_schmidt) / drcodcl(ifac, ivar, 2)
= (viscls+visct/turb_schmidt) / drcodcl(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)rcodcl(ifac, ivar, 3)
= -(viscl+visct) * (grad U).nrcodcl(ifac, ivar, 3)
= -dt * (grad P).nrcodcl(ifac, ivar, 3)
= -Cp*(viscls+visct/turb_schmidt) * (grad T).nrcodcl(ifac, ivar, 3)
= -(viscls+visct/turb_schmidt) * (grad H).nrcodcl(ifac, ivar, 3)
= -(viscls+visct/turb_schmidt) * (grad F).nrcodcl(ifac, ivar, 3)
if icodcl(ifac, ivar)
= 6: Roughness for the rough wall lawrcodcl(ifac, iu, 3)
= roughdrcodcl(ifac, iv, 3)
= roughNote 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.d0rcodcl(ifac, ivar, 2)
= rinfinrcodcl(ifac, ivar, 3)
= 0.d0)Especially, one may have for example:
itypfb(ifac)
= iparoi
which prescribes default wall conditions for all variables at face ifac,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).
For compressible flows, only predefined boundary conditions may be assigned among: iparoi
, isymet
, iesicf
, isspcf
, isopcf
, iephcf
, ieqhcf
iparoi
: standard wallisymet
: standard symmetryiesicf
, isspcf
, isopcf
, iephcf
, ieqhcf
: inlet/outletFor 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 variablesisspcf:
supersonic outlet the user does not prescribe anythingisopcf:
subsonic outlet with prescribed pressure the user presribes the pressureiephcf:
mixed inlet with prescribed total pressure and enthalpy the user prescribes the total pressure and total enthalpyieqhcf:
subsonic inlet with prescribed mass and enthalpy flow to be implementedA few consistency rules between icodcl
codes for variables with non-standard boundary conditions:
itrifb
has not been set yet).itypfb
has been reset before entering the subroutine.Cell value field ids
irom
iviscl
ivisct
icp
field_get_key_int
(ivarfl(isca(iscal)), & kivisl, ...)field
id ibrom
ivar
): field id iflmab
using field_get_key_int(ivarfl(ivar), kbmasf, iflmab)
iel
= ifabor(ifac).Please refer to the boundary conditions section of the theory guide for more informations.
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 | ||
) |
[in] | nvar | total number of variables |
[in] | nscal | total number of scalars |
[out] | icodcl | boundary condition code:
|
[in] | itrifb | indirection for boundary faces ordering |
[in,out] | itypfb | boundary face types |
[out] | izfppp | boundary face zone number |
[in] | dt | time step (per cell) |
[in,out] | rcodcl | boundary condition values:
|