Introduction
This example contains 3 subroutines :
- cs_user_les_inflow_init: definition of global caracteristics of synthetic turbulence inlets
- cs_user_les_inflow_define: definition of the caracteristics of the synthetic turbulence
- cs_user_les_inflow_advanced: accurate definition of target statistics at inlet
Global caracteristics of synthetic turbulence inlets
Purpose
Generation of synthetic turbulence at LES inlets. Definition of global caracteristics of synthetic turbulence inlets: nent and isuisy might be defined.
- nent = Number of inlets
- isuisy = 1: Reading of the LES inflow module restart file
- isuisy = 0: not activated (synthetic turbulence reinitialized)
Local variables declaration
No local variable.
Initializations
Caracteristics of one specific inlet
Purpose
Generation of synthetic turbulence at LES inlets Definition of the caracteristics of the synthetic turbulence inlet 'nument' For each LES inlet, the following parameters might be defined:
- Data relatve to the method employed
- typent indicates the synthetic turbulence method:
- 0 : laminar, no turbulent fluctations
- 1 : random gaussian noise
- 2 : Batten method, based on Fourier mode decomposition
- 3 : Synthetic Eddy Method (SEM)
- nelent indicates the number of "entities" relative to the method (useful only for the Batten method and the SEM):
- for Batten : number of Fourier modes of the turbulent fluctuations
- for SEM : number of synthetic eddies building the fluctuations
- iverbo indicates the verbosity level (listing)
- = 0 no specific output
- > 0 additionnal output (only for SEM)
- Data relative to the LES inflow boundary faces
- nfbent: number of boundary faces of the LES inflow
- lfbent: list of boundary faces of the LES inflow
- Data relative to the flow
- vitent(3): reference mean velocity vector
- enrent: reference turbulent kinetic energy
- dspent: reference dissipation rate
- Note:
- dspent useful only for typent = 2 (Batten) or typent = 3 (SEM).
- Strictly positive values are required for enrent and dspent.
- Accurate specification of the statistics of the flow at LES inlet can be made via the user subroutine cs_user_les_inflow_advanced.
Local variables declaration
No local variable.
Initializations
First synthetic turbulence inlet: the Batten Method is used for boundary faces of color '1'.
if (nument.eq.1) then
typent = 2
nelent = 200
iverbo = 0
call getfbr(
'1',nfbent,lfbent)
vitent(1) = 18.d0
vitent(2) = 0.d0
vitent(3) = 0.d0
enrent = 4.d0
dspent = 4.d0
endif
Second synthetic turbulence inlet: the Synthetic Eddy Method is used for the boundary faces verifying a geometric criterion.
if (nument.eq.1) then
typent = 3
nelent = 2000
iverbo = 1
call getfbr(
'x < 0.1', nfbent, lfbent)
vitent(1) = 12.d0
vitent(2) = 0.d0
vitent(3) = 0.d0
enrent = 3.d0
dspent = 3.d0
endif
Accurate specification of target statistics at inlet
Purpose
Generation of synthetic turbulence at LES inlets. Accurate definition of mean velocity, Reynolds stresses and dissipation rate for each boundary faces of the synthetic turbulence inlet 'nument'.
Usage:
- uvwent(ndim,nfbent) : mean velocity vector
- rijent( 6,nfbent) : Reynolds stresses!
- epsent( nfbent) : dissipation rate
Local variables declaration
integer ii, ifac, iel
double precision d2s3
double precision utau, href, reyfro, yy, yplus, uplus, kplus, eplus
double precision uref2, xdh, xitur, xkent, xeent
double precision, dimension(:), pointer :: cpro_viscl
Example 1
Mean velocity, Reynolds stresses an dissipation are deduced from a wall law for the first synthetic turbulence inlet,
- no refining of the statistics of the flow is provided for the second synthetic turbulence inlet.
if (nument.eq.1) then
utau = uref/20.d0
href = 1.d0
do ii = 1, nfbent
ifac = lfbent(ii)
iel = ifabor(ifac)
reyfro = utau*href/cpro_viscl(iel)
yy = 1.d0-abs(cdgfbo(2,ifac))
uplus = log(1.d0+0.4d0*
yplus)/xkappa &
+ 7.8d0*( 1.d0 - exp(-
yplus/11.d0) &
+ (1.d0 - exp(-
yplus/20.d0))*4.5d0 &
/ (1.d0 + 4.d0*
yplus/reyfro)
eplus = (1.d0/xkappa) &
/ (
yplus**4+15.d0**4)**(0.25d0)
uvwent(1,ii) = uplus*utau
uvwent(2,ii) = 0.d0
uvwent(3,ii) = 0.d0
rijent(1,ii) = d2s3*kplus*utau**2
rijent(2,ii) = d2s3*kplus*utau**2
rijent(3,ii) = d2s3*kplus*utau**2
rijent(4,ii) = 0.d0
rijent(5,ii) = 0.d0
rijent(6,ii) = 0.d0
epsent(ii) = eplus*utau**4/cpro_viscl(iel)
enddo
endif
if (nument.eq.2) then
continue
endif
Example 2
Reynolds stresses and dissipation at the inlet are computed using the turbulence intensity and standard laws for a circular pipe for the first synthetic turbulence inlet,
- no refining of the statistics of the flow is provided for the other synthetic turbulence inlet.
if (nument.eq.1) then
do ii = 1, nfbent
ifac = lfbent(ii)
iel = ifabor(ifac)
uvwent(1,ii) = 1.1d0
uvwent(2,ii) = 1.1d0
uvwent(3,ii) = 1.1d0
uref2 = uvwent(1,ii)**2 &
+ uvwent(2,ii)**2 &
+ uvwent(3,ii)**2
uref2 = max(uref2,1.d-12)
xdh = 0.075d0
xitur = 0.02d0
xkent = epzero
xeent = epzero
call turbulence_bc_ke_turb_intensity&
( uref2, xitur, xdh, xkent, xeent )
rijent(1,ii) = d2s3*xkent
rijent(2,ii) = d2s3*xkent
rijent(3,ii) = d2s3*xkent
rijent(4,ii) = 0.d0
rijent(5,ii) = 0.d0
rijent(6,ii) = 0.d0
epsent(ii) = xeent
enddo
endif