Go to the documentation of this file.
31 #include "twoPhaseSystem.H"
32 #include "virtualMassModel.H"
33 #include "dragModel.H"
44 template<
class BasicTurbulenceModel>
45 continuousGasKEpsilon<BasicTurbulenceModel>::continuousGasKEpsilon
53 const word& propertiesName,
69 liquidTurbulencePtr_(
nullptr),
76 this->runTime_.timeName(),
96 this->printCoeffs(
type);
103 template<
class BasicTurbulenceModel>
108 alphaInversion_.readIfPresent(this->coeffDict());
117 template<
class BasicTurbulenceModel>
137 exp(
min(thetal/thetag, scalar(50))),
143 nutEff_ = omega*liquidTurbulence.
nut();
148 template<
class BasicTurbulenceModel>
152 if (!liquidTurbulencePtr_)
158 refCast<const twoPhaseSystem>(gas.fluid());
161 liquidTurbulencePtr_ =
172 return *liquidTurbulencePtr_;
176 template<
class BasicTurbulenceModel>
186 (this->alpha_ - scalar(0.5))/(alphaInversion_ - 0.5),
199 + (1.0 - blend)*rhoEff()*nutEff_/this->transport().rho()
206 template<
class BasicTurbulenceModel>
222 gas.rho() + (virtualMass.
Cvm() + 3.0/20.0)*liquid.rho()
228 template<
class BasicTurbulenceModel>
240 max(alphaInversion_ -
alpha, scalar(0))
244 liquidTurbulence.
epsilon()/liquidTurbulence.
k(),
245 1.0/
U.time().deltaT()
251 template<
class BasicTurbulenceModel>
256 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
259 phaseTransferCoeff*liquidTurbulence.
k()
260 -
fvm::Sp(phaseTransferCoeff, this->k_);
264 template<
class BasicTurbulenceModel>
269 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
272 phaseTransferCoeff*liquidTurbulence.
epsilon()
273 -
fvm::Sp(phaseTransferCoeff, this->epsilon_);
277 template<
class BasicTurbulenceModel>
290 this->runTime_.timeName(),
296 tk().boundaryField().types()
Defines the attributes of an object for which implicit objectRegistry management is supported,...
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from Foam::string.
Class which solves the volume fraction equations for two phases.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
A class for managing temporary objects.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
static const word propertiesName
Default name of the turbulence properties dictionary.
BasicTurbulenceModel::transportModel transportModel
dimensionedScalar exp(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
BasicTurbulenceModel::rhoField rhoField
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
BasicTurbulenceModel::alphaField alphaField
const turbulenceModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
virtual bool read()
Re-read model coefficients if they have changed.
virtual void correctNut()
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
virtual void correctNut()
Generic dimensioned Type class.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
tmp< volScalarField > phaseTransferCoeff() const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< fvScalarMatrix > epsilonSource() const
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
label k
Boltzmann constant.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
virtual tmp< fvScalarMatrix > kSource() const
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)