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