Go to the documentation of this file.
28 #include "twoPhaseSystem.H"
29 #include "virtualMassModel.H"
40 template<
class BasicTurbulenceModel>
49 const word& propertiesName,
65 liquidTurbulencePtr_(NULL),
72 this->runTime_.timeName(),
92 this->printCoeffs(
type);
99 template<
class BasicTurbulenceModel>
104 alphaInversion_.readIfPresent(this->coeffDict());
115 template<
class BasicTurbulenceModel>
132 exp(
min(thetal/thetag, scalar(50))),
138 nutEff_ = omega*liquidTurbulence.
nut();
142 template<
class BasicTurbulenceModel>
146 if (!liquidTurbulencePtr_)
152 refCast<const twoPhaseSystem>(gas.fluid());
155 liquidTurbulencePtr_ =
166 return *liquidTurbulencePtr_;
170 template<
class BasicTurbulenceModel>
180 (this->alpha_ - scalar(0.5))/(alphaInversion_ - 0.5),
193 + (1.0 - blend)*rhoEff()*nutEff_/this->transport().
rho()
200 template<
class BasicTurbulenceModel>
213 gas.rho() + (
fluid.virtualMass(gas).Cvm() + 3.0/20.0)*liquid.rho()
219 template<
class BasicTurbulenceModel>
231 max(alphaInversion_ -
alpha, scalar(0))
235 liquidTurbulence.
epsilon()/liquidTurbulence.
k(),
236 1.0/
U.time().deltaT()
242 template<
class BasicTurbulenceModel>
247 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
250 phaseTransferCoeff*liquidTurbulence.
k()
251 -
fvm::Sp(phaseTransferCoeff, this->k_);
255 template<
class BasicTurbulenceModel>
260 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
263 phaseTransferCoeff*liquidTurbulence.
epsilon()
264 -
fvm::Sp(phaseTransferCoeff, this->epsilon_);
268 template<
class BasicTurbulenceModel>
281 this->runTime_.timeName(),
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static word groupName(Name name, const word &group)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from string.
Class which solves the volume fraction equations for two phases.
A class for managing temporary objects.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
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)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
continuousGasKEpsilon(const continuousGasKEpsilon &)
BasicTurbulenceModel::rhoField rhoField
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()
static const sphericalTensor I(1)
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
virtual void correctNut()
Macros for easy insertion into run-time selection tables.
cellIbMask *cellIbMaskExt *faceIbMask *cellIbMask boundaryField().evaluateCoupled()
Generic dimensioned Type class.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
tmp< volScalarField > phaseTransferCoeff() const
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< fvScalarMatrix > epsilonSource() const
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
label k
Boltzmann constant.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
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.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)