Skip to content
GitLab
Menu
Projects
Groups
Snippets
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
Development
openfoam
Commits
2fe4ba82
Commit
2fe4ba82
authored
Oct 29, 2008
by
henry
Browse files
Added SAS LES model.
parent
a7090d8f
Changes
2
Hide whitespace changes
Inline
Side-by-side
src/turbulenceModels/LES/incompressible/kOmegaSSTSAS/kOmegaSSTSAS.C
0 → 100644
View file @
2fe4ba82
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2007 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "kOmegaSSTSAS.H"
#include "addToRunTimeSelectionTable.H"
#include "wallDist.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace
Foam
{
namespace
incompressible
{
namespace
LESModels
{
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug
(
kOmegaSSTSAS
,
0
);
addToRunTimeSelectionTable
(
LESModel
,
kOmegaSSTSAS
,
dictionary
);
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
tmp
<
volScalarField
>
kOmegaSSTSAS
::
F1
(
const
volScalarField
&
CDkOmega
)
const
{
volScalarField
CDkOmegaPlus
=
max
(
CDkOmega
,
dimensionedScalar
(
"1.0e-10"
,
dimless
/
sqr
(
dimTime
),
1.0e-10
)
);
volScalarField
arg1
=
min
(
min
(
max
(
(
scalar
(
1
)
/
betaStar_
)
*
sqrt
(
k_
)
/
(
omega_
*
y_
),
scalar
(
500
)
*
nu
()
/
(
sqr
(
y_
)
*
omega_
)
),
(
4
*
alphaOmega2_
)
*
k_
/
(
CDkOmegaPlus
*
sqr
(
y_
))
),
scalar
(
10
)
);
return
tanh
(
pow4
(
arg1
));
}
tmp
<
volScalarField
>
kOmegaSSTSAS
::
F2
()
const
{
volScalarField
arg2
=
min
(
max
(
(
scalar
(
2
)
/
betaStar_
)
*
sqrt
(
k_
)
/
(
omega_
*
y_
),
scalar
(
500
)
*
nu
()
/
(
sqr
(
y_
)
*
omega_
)
),
scalar
(
100
)
);
return
tanh
(
sqr
(
arg2
));
}
tmp
<
volScalarField
>
kOmegaSSTSAS
::
Lvk2
(
const
volScalarField
&
S2
)
const
{
return
kappa_
*
sqrt
(
S2
)
/
(
mag
(
fvc
::
laplacian
(
U
()))
+
dimensionedScalar
(
"ROOTVSMALL"
,
dimensionSet
(
0
,
-
1
,
-
1
,
0
,
0
,
0
,
0
),
ROOTVSMALL
)
);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
kOmegaSSTSAS
::
kOmegaSSTSAS
(
const
volVectorField
&
U
,
const
surfaceScalarField
&
phi
,
transportModel
&
transport
,
const
word
&
modelName
)
:
LESModel
(
modelName
,
U
,
phi
,
transport
),
alphaK1_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"alphaK1"
,
coeffDict
(),
0
.
85034
)
),
alphaK2_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"alphaK2"
,
coeffDict
(),
1
.
0
)
),
alphaOmega1_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"alphaOmega1"
,
coeffDict
(),
0
.
5
)
),
alphaOmega2_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"alphaOmega2"
,
coeffDict
(),
0
.
85616
)
),
gamma1_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"gamma1"
,
coeffDict
(),
0
.
5532
)
),
gamma2_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"gamma2"
,
coeffDict
(),
0
.
4403
)
),
beta1_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"beta1"
,
coeffDict
(),
0
.
075
)
),
beta2_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"beta2"
,
coeffDict
(),
0
.
0828
)
),
betaStar_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"betaStar"
,
coeffDict
(),
0
.
09
)
),
a1_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"a1"
,
coeffDict
(),
0
.
31
)
),
c1_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"c1"
,
coeffDict
(),
10
.
0
)
),
alphaPhi_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"alphaPhi"
,
coeffDict
(),
0
.
666667
)
),
zetaTilda2_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"zetaTilda2"
,
coeffDict
(),
1
.
755
)
),
FSAS_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"FSAS"
,
coeffDict
(),
1
.
25
)
),
omega0_
(
"omega0"
,
dimless
/
dimTime
,
SMALL
),
omegaSmall_
(
"omegaSmall"
,
dimless
/
dimTime
,
SMALL
),
y_
(
mesh_
),
Cmu_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"Cmu"
,
coeffDict
(),
0
.
09
)
),
kappa_
(
dimensioned
<
scalar
>::
lookupOrAddToDict
(
"kappa"
,
*
this
,
0
.
4187
)
),
k_
(
IOobject
(
"k"
,
runTime_
.
timeName
(),
mesh_
,
IOobject
::
MUST_READ
,
IOobject
::
AUTO_WRITE
),
mesh_
),
omega_
(
IOobject
(
"omega"
,
runTime_
.
timeName
(),
mesh_
,
IOobject
::
MUST_READ
,
IOobject
::
AUTO_WRITE
),
mesh_
),
nuSgs_
(
IOobject
(
"nuSgs"
,
runTime_
.
timeName
(),
mesh_
,
IOobject
::
MUST_READ
,
IOobject
::
AUTO_WRITE
),
mesh_
)
{
printCoeffs
();
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void
kOmegaSSTSAS
::
correct
(
const
tmp
<
volTensorField
>&
gradU
)
{
LESModel
::
correct
(
gradU
);
if
(
mesh_
.
changing
())
{
y_
.
correct
();
}
volScalarField
S2
=
magSqr
(
symm
(
gradU
()));
gradU
.
clear
();
volVectorField
gradK
=
fvc
::
grad
(
k_
);
volVectorField
gradOmega
=
fvc
::
grad
(
omega_
);
volScalarField
L
=
sqrt
(
k_
)
/
(
pow
(
Cmu_
,
0
.
25
)
*
(
omega_
+
omegaSmall_
));
volScalarField
CDkOmega
=
(
2
.
0
*
alphaOmega2_
)
*
(
gradK
&
gradOmega
)
/
(
omega_
+
omegaSmall_
);
volScalarField
F1
=
this
->
F1
(
CDkOmega
);
volScalarField
G
=
nuSgs_
*
2
.
0
*
S2
;
// Turbulent kinetic energy equation
solve
(
fvm
::
ddt
(
k_
)
+
fvm
::
div
(
phi
(),
k_
)
-
fvm
::
Sp
(
fvc
::
div
(
phi
()),
k_
)
-
fvm
::
laplacian
(
DkEff
(
F1
),
k_
)
==
min
(
G
,
c1_
*
betaStar_
*
k_
*
omega_
)
-
fvm
::
Sp
(
betaStar_
*
omega_
,
k_
)
);
bound
(
k_
,
k0
());
volScalarField
grad_omega_k
=
max
(
magSqr
(
gradOmega
)
/
sqr
(
omega_
+
omegaSmall_
),
magSqr
(
gradK
)
/
sqr
(
k_
+
k0
())
);
// Turbulent frequency equation
solve
(
fvm
::
ddt
(
omega_
)
+
fvm
::
div
(
phi
(),
omega_
)
-
fvm
::
Sp
(
fvc
::
div
(
phi
()),
omega_
)
-
fvm
::
laplacian
(
DomegaEff
(
F1
),
omega_
)
==
gamma
(
F1
)
*
2
.
0
*
S2
-
fvm
::
Sp
(
beta
(
F1
)
*
omega_
,
omega_
)
-
fvm
::
SuSp
// cross diffusion term
(
(
F1
-
scalar
(
1
))
*
CDkOmega
/
omega_
,
omega_
)
+
FSAS_
*
max
(
dimensionedScalar
(
"zero"
,
dimensionSet
(
0
,
0
,
-
2
,
0
,
0
),
0
.
),
zetaTilda2_
*
kappa_
*
S2
*
(
L
/
Lvk2
(
S2
))
-
2
.
0
/
alphaPhi_
*
k_
*
grad_omega_k
)
);
bound
(
omega_
,
omega0_
);
// Re-calculate viscosity
nuSgs_
==
a1_
*
k_
/
max
(
a1_
*
omega_
,
F2
()
*
sqrt
(
S2
));
nuSgs_
.
correctBoundaryConditions
();
}
tmp
<
volScalarField
>
kOmegaSSTSAS
::
epsilon
()
const
{
return
2
.
0
*
nuEff
()
*
magSqr
(
symm
(
fvc
::
grad
(
U
())));
}
tmp
<
volSymmTensorField
>
kOmegaSSTSAS
::
B
()
const
{
return
((
2
.
0
/
3
.
0
)
*
I
)
*
k
()
-
nuSgs
()
*
twoSymm
(
fvc
::
grad
(
U
()));
}
tmp
<
volSymmTensorField
>
kOmegaSSTSAS
::
devBeff
()
const
{
return
-
nuEff
()
*
dev
(
twoSymm
(
fvc
::
grad
(
U
())));
}
tmp
<
fvVectorMatrix
>
kOmegaSSTSAS
::
divDevBeff
(
volVectorField
&
U
)
const
{
return
(
-
fvm
::
laplacian
(
nuEff
(),
U
)
-
fvc
::
div
(
nuEff
()
*
dev
(
fvc
::
grad
(
U
)().
T
()))
);
}
bool
kOmegaSSTSAS
::
read
()
{
if
(
LESModel
::
read
())
{
alphaK1_
.
readIfPresent
(
coeffDict
());
alphaK2_
.
readIfPresent
(
coeffDict
());
alphaOmega1_
.
readIfPresent
(
coeffDict
());
alphaOmega2_
.
readIfPresent
(
coeffDict
());
gamma1_
.
readIfPresent
(
coeffDict
());
gamma2_
.
readIfPresent
(
coeffDict
());
beta1_
.
readIfPresent
(
coeffDict
());
beta2_
.
readIfPresent
(
coeffDict
());
betaStar_
.
readIfPresent
(
coeffDict
());
a1_
.
readIfPresent
(
coeffDict
());
c1_
.
readIfPresent
(
coeffDict
());
alphaPhi_
.
readIfPresent
(
coeffDict
());
zetaTilda2_
.
readIfPresent
(
coeffDict
());
FSAS_
.
readIfPresent
(
coeffDict
());
return
true
;
}
else
{
return
false
;
}
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
}
// End namespace LESModels
}
// End namespace incompressible
}
// End namespace Foam
// ************************************************************************* //
src/turbulenceModels/LES/incompressible/kOmegaSSTSAS/kOmegaSSTSAS.H
0 → 100644
View file @
2fe4ba82
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2007 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Class
Foam::LESmodels::kOmegaSSTSAS
Description
kOmegaSSTSAS LES turbulence model for incompressible flows
Note: does not have an explicit dependency on spatial discretisation
i.e. LESdelta not used.
SourceFiles
kOmegaSSTSAS.C
\*---------------------------------------------------------------------------*/
#ifndef kOmegaSSTSAS_H
#define kOmegaSSTSAS_H
#include "LESModel.H"
#include "volFields.H"
#include "wallDist.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace
Foam
{
namespace
incompressible
{
namespace
LESModels
{
/*---------------------------------------------------------------------------*\
Class kOmegaSSTSAS Declaration
\*---------------------------------------------------------------------------*/
class
kOmegaSSTSAS
:
public
LESModel
{
// Private member functions
// Disallow default bitwise copy construct and assignment
kOmegaSSTSAS
(
const
kOmegaSSTSAS
&
);
kOmegaSSTSAS
&
operator
=
(
const
kOmegaSSTSAS
&
);
protected:
// Protected data
// Model constants
dimensionedScalar
alphaK1_
;
dimensionedScalar
alphaK2_
;
dimensionedScalar
alphaOmega1_
;
dimensionedScalar
alphaOmega2_
;
dimensionedScalar
gamma1_
;
dimensionedScalar
gamma2_
;
dimensionedScalar
beta1_
;
dimensionedScalar
beta2_
;
dimensionedScalar
betaStar_
;
dimensionedScalar
a1_
;
dimensionedScalar
c1_
;
dimensionedScalar
alphaPhi_
;
dimensionedScalar
zetaTilda2_
;
dimensionedScalar
FSAS_
;
dimensionedScalar
omega0_
;
dimensionedScalar
omegaSmall_
;
wallDist
y_
;
dimensionedScalar
Cmu_
;
dimensionedScalar
kappa_
;
// Fields
volScalarField
k_
;
volScalarField
omega_
;
volScalarField
nuSgs_
;
// Protected member functions
tmp
<
volScalarField
>
Lvk2
(
const
volScalarField
&
S2
)
const
;
tmp
<
volScalarField
>
F1
(
const
volScalarField
&
CDkOmega
)
const
;
tmp
<
volScalarField
>
F2
()
const
;
tmp
<
volScalarField
>
blend
(
const
volScalarField
&
F1
,
const
dimensionedScalar
&
psi1
,
const
dimensionedScalar
&
psi2
)
const
{
return
F1
*
(
psi1
-
psi2
)
+
psi2
;
}
tmp
<
volScalarField
>
alphaK
(
const
volScalarField
&
F1
)
const
{
return
blend
(
F1
,
alphaK1_
,
alphaK2_
);
}
tmp
<
volScalarField
>
alphaOmega
(
const
volScalarField
&
F1
)
const
{
return
blend
(
F1
,
alphaOmega1_
,
alphaOmega2_
);
}
tmp
<
volScalarField
>
beta
(
const
volScalarField
&
F1
)
const
{
return
blend
(
F1
,
beta1_
,
beta2_
);
}
tmp
<
volScalarField
>
gamma
(
const
volScalarField
&
F1