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:
Henry
2015-01-25 12:37:37 +00:00
parent 2fb88d3e4e
commit f50d23f859
2 changed files with 36 additions and 23 deletions

View File

@ -24,9 +24,9 @@ License
\*---------------------------------------------------------------------------*/
#include "buoyantKEpsilon.H"
#include "addToRunTimeSelectionTable.H"
#include "uniformDimensionedFields.H"
#include "fvcGrad.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -68,7 +68,7 @@ buoyantKEpsilon<BasicTurbulenceModel>::buoyantKEpsilon
(
"Cg",
this->coeffDict_,
0.85
1.0
)
)
{
@ -107,8 +107,8 @@ buoyantKEpsilon<BasicTurbulenceModel>::Gcoef() const
lookupObject<uniformDimensionedVectorField>("g");
return
this->alpha_*(this->Cmu_/this->sigmak_)*this->k_
*(g & fvc::grad(this->rho_))/(this->epsilon_ + this->epsilonMin_);
(Cg_*this->Cmu_)*this->alpha_*this->k_*(g & fvc::grad(this->rho_))
/(this->epsilon_ + this->epsilonMin_);
}
@ -124,7 +124,19 @@ template<class BasicTurbulenceModel>
tmp<fvScalarMatrix>
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_);
}

View File

@ -28,33 +28,34 @@ Group
grpRASTurbulence
Description
k-epsilon model for the gas-phase in a two-phase system
supporting phase-inversion.
Additional buoyancy generation/dissipation term applied to the
k and epsilon equations of the standard k-epsilon model.
In the limit that the gas-phase fraction approaches zero a contribution from
the other phase is blended into the k and epsilon equations up to the
phase-fraction of alphaInversion at which point phase-inversion is
considered to have occurred and the model reverts to the pure single-phase
form.
Reference:
\verbatim
Henkes, R.A.W.M., Van Der Vlugt, F.F. & Hoogendoorn, C.J. (1991).
Natural Convection Flow in a Square Cavity Calculated with
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
on which a more physical model may be built.
This implementation is based on the density rather than temperature gradient
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
buoyantKEpsilonCoeffs
{
Cmu 0.09;
C1 1.44;
C2 1.92;
C3 -0.33;
Cg 0.85;
sigmak 1.0;
sigmaEps 1.3;
alphaInversion 0.7;
Cg 1.0;
}
\endverbatim
SeeAlso
Foam::RASModels::kEpsilon
SourceFiles
buoyantKEpsilon.C