mirror of
https://github.com/OpenFOAM/OpenFOAM-6.git
synced 2025-12-08 06:57:46 +00:00
buoyantKEpsilon: Changed the additional buoyancy generation/dissipation
term to the more commonly used form of Henkes, R.A.W.M., Van Der Vlugt, F.F. & Hoogendoorn, C.J. (1991).
This commit is contained in:
@ -24,9 +24,9 @@ License
|
|||||||
\*---------------------------------------------------------------------------*/
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
#include "buoyantKEpsilon.H"
|
#include "buoyantKEpsilon.H"
|
||||||
#include "addToRunTimeSelectionTable.H"
|
|
||||||
#include "uniformDimensionedFields.H"
|
#include "uniformDimensionedFields.H"
|
||||||
#include "fvcGrad.H"
|
#include "fvcGrad.H"
|
||||||
|
#include "addToRunTimeSelectionTable.H"
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
@ -68,7 +68,7 @@ buoyantKEpsilon<BasicTurbulenceModel>::buoyantKEpsilon
|
|||||||
(
|
(
|
||||||
"Cg",
|
"Cg",
|
||||||
this->coeffDict_,
|
this->coeffDict_,
|
||||||
0.85
|
1.0
|
||||||
)
|
)
|
||||||
)
|
)
|
||||||
{
|
{
|
||||||
@ -107,8 +107,8 @@ buoyantKEpsilon<BasicTurbulenceModel>::Gcoef() const
|
|||||||
lookupObject<uniformDimensionedVectorField>("g");
|
lookupObject<uniformDimensionedVectorField>("g");
|
||||||
|
|
||||||
return
|
return
|
||||||
this->alpha_*(this->Cmu_/this->sigmak_)*this->k_
|
(Cg_*this->Cmu_)*this->alpha_*this->k_*(g & fvc::grad(this->rho_))
|
||||||
*(g & fvc::grad(this->rho_))/(this->epsilon_ + this->epsilonMin_);
|
/(this->epsilon_ + this->epsilonMin_);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
@ -124,7 +124,19 @@ template<class BasicTurbulenceModel>
|
|||||||
tmp<fvScalarMatrix>
|
tmp<fvScalarMatrix>
|
||||||
buoyantKEpsilon<BasicTurbulenceModel>::epsilonSource() const
|
buoyantKEpsilon<BasicTurbulenceModel>::epsilonSource() const
|
||||||
{
|
{
|
||||||
return -fvm::SuSp(this->C1_*(1.0 - this->Cg_)*Gcoef(), this->epsilon_);
|
const uniformDimensionedVectorField& g =
|
||||||
|
this->mesh_.objectRegistry::template
|
||||||
|
lookupObject<uniformDimensionedVectorField>("g");
|
||||||
|
|
||||||
|
vector gHat(g.value()/mag(g.value()));
|
||||||
|
|
||||||
|
volScalarField v(gHat & this->U_);
|
||||||
|
volScalarField u
|
||||||
|
(
|
||||||
|
mag(this->U_ - gHat*v) + dimensionedScalar("SMALL", dimVelocity, SMALL)
|
||||||
|
);
|
||||||
|
|
||||||
|
return -fvm::SuSp(this->C1_*tanh(mag(v)/u)*Gcoef(), this->epsilon_);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
|||||||
@ -28,33 +28,34 @@ Group
|
|||||||
grpRASTurbulence
|
grpRASTurbulence
|
||||||
|
|
||||||
Description
|
Description
|
||||||
k-epsilon model for the gas-phase in a two-phase system
|
Additional buoyancy generation/dissipation term applied to the
|
||||||
supporting phase-inversion.
|
k and epsilon equations of the standard k-epsilon model.
|
||||||
|
|
||||||
In the limit that the gas-phase fraction approaches zero a contribution from
|
Reference:
|
||||||
the other phase is blended into the k and epsilon equations up to the
|
\verbatim
|
||||||
phase-fraction of alphaInversion at which point phase-inversion is
|
Henkes, R.A.W.M., Van Der Vlugt, F.F. & Hoogendoorn, C.J. (1991).
|
||||||
considered to have occurred and the model reverts to the pure single-phase
|
Natural Convection Flow in a Square Cavity Calculated with
|
||||||
form.
|
Low-Reynolds-Number Turbulence Models.
|
||||||
|
Int. J. Heat Mass Transfer, 34, 1543-1557.
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
This model is unpublished and is provided as a stable numerical framework
|
This implementation is based on the density rather than temperature gradient
|
||||||
on which a more physical model may be built.
|
extending the applicability to systems in which the density gradient may be
|
||||||
|
generated by variation of composition rather than temperature. Further, the
|
||||||
|
1/Prt coefficient is replaced by Cg to provide more general control of
|
||||||
|
model.
|
||||||
|
|
||||||
The default model coefficients correspond to the following:
|
The default model coefficients are
|
||||||
\verbatim
|
\verbatim
|
||||||
buoyantKEpsilonCoeffs
|
buoyantKEpsilonCoeffs
|
||||||
{
|
{
|
||||||
Cmu 0.09;
|
Cg 1.0;
|
||||||
C1 1.44;
|
|
||||||
C2 1.92;
|
|
||||||
C3 -0.33;
|
|
||||||
Cg 0.85;
|
|
||||||
sigmak 1.0;
|
|
||||||
sigmaEps 1.3;
|
|
||||||
alphaInversion 0.7;
|
|
||||||
}
|
}
|
||||||
\endverbatim
|
\endverbatim
|
||||||
|
|
||||||
|
SeeAlso
|
||||||
|
Foam::RASModels::kEpsilon
|
||||||
|
|
||||||
SourceFiles
|
SourceFiles
|
||||||
buoyantKEpsilon.C
|
buoyantKEpsilon.C
|
||||||
|
|
||||||
|
|||||||
Reference in New Issue
Block a user