Skip to content
GitLab
Menu
Projects
Groups
Snippets
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
Development
openfoam
Commits
dcfb9280
Commit
dcfb9280
authored
Aug 02, 2013
by
Henry
Browse files
TurbulenceModels: Reorganised, cleaned, added documentation and made more consistent
parent
8751d679
Changes
17
Hide whitespace changes
Inline
Side-by-side
applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/phaseIncompressibleTurbulenceModels/kineticTheoryModels/kineticTheoryModel/kineticTheoryModel.C
View file @
dcfb9280
...
...
@@ -129,7 +129,12 @@ Foam::kineticTheoryModel::kineticTheoryModel
U
.
mesh
(),
dimensionedScalar
(
"zero"
,
dimensionSet
(
0
,
2
,
-
1
,
0
,
0
),
0
.
0
)
)
{}
{
if
(
type
==
typeName
)
{
this
->
printCoeffs
(
type
);
}
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
...
...
@@ -140,6 +145,199 @@ Foam::kineticTheoryModel::~kineticTheoryModel()
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool
Foam
::
kineticTheoryModel
::
read
()
{
if
(
eddyViscosity
<
RASModel
<
PhaseIncompressibleTurbulenceModel
<
phaseModel
>
>
>::
read
()
)
{
this
->
coeffDict
().
lookup
(
"equilibrium"
)
>>
equilibrium_
;
e_
.
readIfPresent
(
this
->
coeffDict
());
alphaMax_
.
readIfPresent
(
this
->
coeffDict
());
alphaMinFriction_
.
readIfPresent
(
this
->
coeffDict
());
viscosityModel_
->
read
();
conductivityModel_
->
read
();
radialModel_
->
read
();
granularPressureModel_
->
read
();
frictionalStressModel_
->
read
();
return
true
;
}
else
{
return
false
;
}
}
Foam
::
tmp
<
Foam
::
volScalarField
>
Foam
::
kineticTheoryModel
::
k
()
const
{
notImplemented
(
"kineticTheoryModel::k()"
);
return
nut_
;
}
Foam
::
tmp
<
Foam
::
volScalarField
>
Foam
::
kineticTheoryModel
::
epsilon
()
const
{
notImplemented
(
"kineticTheoryModel::epsilon()"
);
return
nut_
;
}
Foam
::
tmp
<
Foam
::
volSymmTensorField
>
Foam
::
kineticTheoryModel
::
R
()
const
{
return
tmp
<
volSymmTensorField
>
(
new
volSymmTensorField
(
IOobject
(
IOobject
::
groupName
(
"R"
,
this
->
U_
.
group
()),
this
->
runTime_
.
timeName
(),
this
->
mesh_
,
IOobject
::
NO_READ
,
IOobject
::
NO_WRITE
),
-
(
this
->
nut_
)
*
dev
(
twoSymm
(
fvc
::
grad
(
this
->
U_
)))
-
(
lambda_
*
fvc
::
div
(
this
->
phi_
))
*
symmTensor
::
I
)
);
}
/*
Foam::tmp<Foam::volScalarField> Foam::kineticTheoryModel::pp() const
{
// Particle pressure coefficient
// Coefficient in front of Theta (Eq. 3.22, p. 45)
volScalarField PsCoeff
(
granularPressureModel_->granularPressureCoeff
(
alpha,
gs0,
rho,
e_
)
);
// Frictional pressure
volScalarField pf
(
frictionalStressModel_->frictionalPressure
(
alpha,
alphaMinFriction_,
alphaMax_
)
);
// Return total particle pressure
return PsCoeff*Theta_ + pf;
}
*/
Foam
::
tmp
<
Foam
::
volScalarField
>
Foam
::
kineticTheoryModel
::
pPrime
()
const
{
// Local references
const
volScalarField
&
alpha
=
this
->
alpha_
;
const
volScalarField
&
rho
=
phase_
.
rho
();
return
(
Theta_
*
granularPressureModel_
->
granularPressureCoeffPrime
(
alpha
,
radialModel_
->
g0
(
alpha
,
alphaMinFriction_
,
alphaMax_
),
radialModel_
->
g0prime
(
alpha
,
alphaMinFriction_
,
alphaMax_
),
rho
,
e_
)
+
frictionalStressModel_
->
frictionalPressurePrime
(
alpha
,
alphaMinFriction_
,
alphaMax_
)
);
}
Foam
::
tmp
<
Foam
::
surfaceScalarField
>
Foam
::
kineticTheoryModel
::
pPrimef
()
const
{
// Local references
const
volScalarField
&
alpha
=
this
->
alpha_
;
const
volScalarField
&
rho
=
phase_
.
rho
();
return
fvc
::
interpolate
(
Theta_
*
granularPressureModel_
->
granularPressureCoeffPrime
(
alpha
,
radialModel_
->
g0
(
alpha
,
alphaMinFriction_
,
alphaMax_
),
radialModel_
->
g0prime
(
alpha
,
alphaMinFriction_
,
alphaMax_
),
rho
,
e_
)
+
frictionalStressModel_
->
frictionalPressurePrime
(
alpha
,
alphaMinFriction_
,
alphaMax_
)
);
}
Foam
::
tmp
<
Foam
::
volSymmTensorField
>
Foam
::
kineticTheoryModel
::
devRhoReff
()
const
{
return
tmp
<
volSymmTensorField
>
(
new
volSymmTensorField
(
IOobject
(
IOobject
::
groupName
(
"devRhoReff"
,
this
->
U_
.
group
()),
this
->
runTime_
.
timeName
(),
this
->
mesh_
,
IOobject
::
NO_READ
,
IOobject
::
NO_WRITE
),
-
(
this
->
rho_
*
this
->
nut_
)
*
dev
(
twoSymm
(
fvc
::
grad
(
this
->
U_
)))
-
((
this
->
rho_
*
lambda_
)
*
fvc
::
div
(
this
->
phi_
))
*
symmTensor
::
I
)
);
}
Foam
::
tmp
<
Foam
::
fvVectorMatrix
>
Foam
::
kineticTheoryModel
::
divDevRhoReff
(
volVectorField
&
U
)
const
{
return
(
-
fvm
::
laplacian
(
this
->
rho_
*
this
->
nut_
,
U
)
-
fvc
::
div
(
(
this
->
rho_
*
this
->
nut_
)
*
dev2
(
T
(
fvc
::
grad
(
U
)))
+
((
this
->
rho_
*
lambda_
)
*
fvc
::
div
(
this
->
phi_
))
*
dimensioned
<
symmTensor
>
(
"I"
,
dimless
,
symmTensor
::
I
)
)
);
}
void
Foam
::
kineticTheoryModel
::
correct
()
{
// Local references
...
...
@@ -343,202 +541,4 @@ void Foam::kineticTheoryModel::correct()
}
}
/*
Foam::tmp<Foam::volScalarField> Foam::kineticTheoryModel::pp() const
{
// Particle pressure coefficient
// Coefficient in front of Theta (Eq. 3.22, p. 45)
volScalarField PsCoeff
(
granularPressureModel_->granularPressureCoeff
(
alpha,
gs0,
rho,
e_
)
);
// Frictional pressure
volScalarField pf
(
frictionalStressModel_->frictionalPressure
(
alpha,
alphaMinFriction_,
alphaMax_
)
);
// Return total particle pressure
return PsCoeff*Theta_ + pf;
}
*/
Foam
::
tmp
<
Foam
::
volScalarField
>
Foam
::
kineticTheoryModel
::
pPrime
()
const
{
// Local references
const
volScalarField
&
alpha
=
this
->
alpha_
;
const
volScalarField
&
rho
=
phase_
.
rho
();
return
(
Theta_
*
granularPressureModel_
->
granularPressureCoeffPrime
(
alpha
,
radialModel_
->
g0
(
alpha
,
alphaMinFriction_
,
alphaMax_
),
radialModel_
->
g0prime
(
alpha
,
alphaMinFriction_
,
alphaMax_
),
rho
,
e_
)
+
frictionalStressModel_
->
frictionalPressurePrime
(
alpha
,
alphaMinFriction_
,
alphaMax_
)
);
}
Foam
::
tmp
<
Foam
::
surfaceScalarField
>
Foam
::
kineticTheoryModel
::
pPrimef
()
const
{
// Local references
const
volScalarField
&
alpha
=
this
->
alpha_
;
const
volScalarField
&
rho
=
phase_
.
rho
();
return
fvc
::
interpolate
(
Theta_
*
granularPressureModel_
->
granularPressureCoeffPrime
(
alpha
,
radialModel_
->
g0
(
alpha
,
alphaMinFriction_
,
alphaMax_
),
radialModel_
->
g0prime
(
alpha
,
alphaMinFriction_
,
alphaMax_
),
rho
,
e_
)
+
frictionalStressModel_
->
frictionalPressurePrime
(
alpha
,
alphaMinFriction_
,
alphaMax_
)
);
}
void
Foam
::
kineticTheoryModel
::
correctNut
()
{}
Foam
::
tmp
<
Foam
::
volScalarField
>
Foam
::
kineticTheoryModel
::
k
()
const
{
notImplemented
(
"kineticTheoryModel::k()"
);
return
nut_
;
}
Foam
::
tmp
<
Foam
::
volScalarField
>
Foam
::
kineticTheoryModel
::
epsilon
()
const
{
notImplemented
(
"kineticTheoryModel::epsilon()"
);
return
nut_
;
}
Foam
::
tmp
<
Foam
::
volSymmTensorField
>
Foam
::
kineticTheoryModel
::
R
()
const
{
return
tmp
<
volSymmTensorField
>
(
new
volSymmTensorField
(
IOobject
(
IOobject
::
groupName
(
"R"
,
this
->
U_
.
group
()),
this
->
runTime_
.
timeName
(),
this
->
mesh_
,
IOobject
::
NO_READ
,
IOobject
::
NO_WRITE
),
-
(
this
->
nut_
)
*
dev
(
twoSymm
(
fvc
::
grad
(
this
->
U_
)))
-
(
lambda_
*
fvc
::
div
(
this
->
phi_
))
*
symmTensor
::
I
)
);
}
Foam
::
tmp
<
Foam
::
volSymmTensorField
>
Foam
::
kineticTheoryModel
::
devRhoReff
()
const
{
return
tmp
<
volSymmTensorField
>
(
new
volSymmTensorField
(
IOobject
(
IOobject
::
groupName
(
"devRhoReff"
,
this
->
U_
.
group
()),
this
->
runTime_
.
timeName
(),
this
->
mesh_
,
IOobject
::
NO_READ
,
IOobject
::
NO_WRITE
),
-
(
this
->
rho_
*
this
->
nut_
)
*
dev
(
twoSymm
(
fvc
::
grad
(
this
->
U_
)))
-
((
this
->
rho_
*
lambda_
)
*
fvc
::
div
(
this
->
phi_
))
*
symmTensor
::
I
)
);
}
Foam
::
tmp
<
Foam
::
fvVectorMatrix
>
Foam
::
kineticTheoryModel
::
divDevRhoReff
(
volVectorField
&
U
)
const
{
return
(
-
fvm
::
laplacian
(
this
->
rho_
*
this
->
nut_
,
U
)
-
fvc
::
div
(
(
this
->
rho_
*
this
->
nut_
)
*
dev2
(
T
(
fvc
::
grad
(
U
)))
+
((
this
->
rho_
*
lambda_
)
*
fvc
::
div
(
this
->
phi_
))
*
dimensioned
<
symmTensor
>
(
"I"
,
dimless
,
symmTensor
::
I
)
)
);
}
bool
Foam
::
kineticTheoryModel
::
read
()
{
if
(
eddyViscosity
<
RASModel
<
PhaseIncompressibleTurbulenceModel
<
phaseModel
>
>
>::
read
()
)
{
this
->
coeffDict
().
lookup
(
"equilibrium"
)
>>
equilibrium_
;
e_
.
readIfPresent
(
this
->
coeffDict
());
alphaMax_
.
readIfPresent
(
this
->
coeffDict
());
alphaMinFriction_
.
readIfPresent
(
this
->
coeffDict
());
viscosityModel_
->
read
();
conductivityModel_
->
read
();
radialModel_
->
read
();
granularPressureModel_
->
read
();
frictionalStressModel_
->
read
();
return
true
;
}
else
{
return
false
;
}
}
// ************************************************************************* //
applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/phaseIncompressibleTurbulenceModels/kineticTheoryModels/kineticTheoryModel/kineticTheoryModel.H
View file @
dcfb9280
...
...
@@ -25,6 +25,17 @@ Class
Foam::kineticTheoryModel
Description
Kinetic theory particle phase RAS model
Reference:
\verbatim
"Derivation, implementation, and validation of computer simulation
models for gas-solid fluidized beds",
B.G.M. van Wachem,
Ph.D. Thesis, Delft University of Technology, Amsterdam, 2000.
\endverbatim
There are no default model coefficients.
SourceFiles
kineticTheoryModel.C
...
...
@@ -118,6 +129,9 @@ class kineticTheoryModel
// Private Member Functions
void
correctNut
()
{}
//- Disallow default bitwise copy construct
kineticTheoryModel
(
const
kineticTheoryModel
&
);
...
...
@@ -125,13 +139,6 @@ class kineticTheoryModel
void
operator
=
(
const
kineticTheoryModel
&
);
protected:
// Protected member functions
virtual
void
correctNut
();
public:
//- Runtime type information
...
...
@@ -160,6 +167,9 @@ public:
// Member Functions
//- Re-read model coefficients if they have changed
virtual
bool
read
();
//- Return the effective viscosity
virtual
tmp
<
volScalarField
>
nuEff
()
const
{
...
...
@@ -197,9 +207,6 @@ public:
//- Solve the kinetic theory equations and correct the viscosity
virtual
void
correct
();
//- Re-read model coefficients if they have changed
virtual
bool
read
();
};
...
...
applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/phaseIncompressibleTurbulenceModels/phaseIncompressibleTurbulenceModels.C
View file @
dcfb9280
...
...
@@ -182,7 +182,10 @@ namespace Foam
#include
"LESModel.H"
#include
"Smagorinsky.H"
#include
"SmagorinskyZhang.H"
#include
"kEqn.H"
#include
"NicenoKEqn.H"
#include
"continuousGasKEqn.H"
namespace
Foam
{
...
...
@@ -215,6 +218,21 @@ namespace Foam
);
}
namespace
LESModels
{
typedef
SmagorinskyZhang
<
incompressibleTransportTurbulenceModel
>
incompressibleSmagorinskyZhang
;
defineNamedTemplateTypeNameAndDebug
(
incompressibleSmagorinskyZhang
,
0
);
addToRunTimeSelectionTable
(
incompressibleLESModel
,
incompressibleSmagorinskyZhang
,
dictionary
);
}
namespace
LESModels
{
typedef
kEqn
<
incompressibleTransportTurbulenceModel
>
...
...
@@ -229,6 +247,36 @@ namespace Foam
dictionary
);
}
namespace
LESModels
{
typedef
NicenoKEqn
<
incompressibleTransportTurbulenceModel
>
incompressibleNicenoKEqn
;
defineNamedTemplateTypeNameAndDebug
(
incompressibleNicenoKEqn
,
0
);
addToRunTimeSelectionTable
(
incompressibleLESModel
,
incompressibleNicenoKEqn
,
dictionary
);
}
namespace
LESModels
{
typedef
continuousGasKEqn
<
incompressibleTransportTurbulenceModel
>
incompressiblecontinuousGasKEqn
;
defineNamedTemplateTypeNameAndDebug
(
incompressiblecontinuousGasKEqn
,
0
);
addToRunTimeSelectionTable
(
incompressibleLESModel
,
incompressiblecontinuousGasKEqn
,
dictionary