40 template<
class BasicMomentumTransportModel>
63 gasTurbulencePtr_(
nullptr),
105 if (type == typeName)
107 this->printCoeffs(type);
114 template<
class BasicMomentumTransportModel>
119 alphaInversion_.readIfPresent(this->coeffDict());
120 Cp_.readIfPresent(this->coeffDict());
121 C3_.readIfPresent(this->coeffDict());
122 Cmub_.readIfPresent(this->coeffDict());
133 template<
class BasicMomentumTransportModel>
136 typename BasicMomentumTransportModel::transportModel
140 if (!gasTurbulencePtr_)
156 momentumTransportModel::typeName,
162 return *gasTurbulencePtr_;
166 template<
class BasicMomentumTransportModel>
170 gasTurbulence = this->gasTurbulence();
173 this->Cmu_*
sqr(this->k_)/this->epsilon_
175 *(
mag(this->U_ - gasTurbulence.
U()));
177 this->nut_.correctBoundaryConditions();
182 template<
class BasicMomentumTransportModel>
186 gasTurbulence = this->gasTurbulence();
201 +
pow(drag.
CdRe()*liquid.thermo().nu()/gas.d(), 4.0/3.0)
212 template<
class BasicMomentumTransportModel>
224 max(alphaInversion_ - alpha, scalar(0))
231 template<
class BasicMomentumTransportModel>
238 gasTurbulence = this->gasTurbulence();
240 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
244 + phaseTransferCoeff*gasTurbulence.
k()
245 -
fvm::Sp(phaseTransferCoeff, this->k_);
249 template<
class BasicMomentumTransportModel>
257 gasTurbulence = this->gasTurbulence();
259 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
262 alpha*rho*this->C3_*this->epsilon_*bubbleG()/this->k_
263 + phaseTransferCoeff*gasTurbulence.
epsilon()
264 -
fvm::Sp(phaseTransferCoeff, this->epsilon_);
268 template<
class BasicMomentumTransportModel>
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Templated abstract base class for multiphase compressible turbulence models.
const transportModel & transport() const
Access function to incompressible transport model.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual bool read()
Read model coefficients if they have changed.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Generic dimensioned Type class.
LaheyKEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Continuous-phase k-epsilon model including bubble-generated turbulence.
Info<< "Reading strained laminar flame speed field Su\"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\"<< endl;volSymmTensorField nsv(IOobject("nsv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);IOdictionary PDRProperties(IOobject("PDRProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));autoPtr< PDRDragModel > drag
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Class to represent a system of phases and model interfacial transfers between them.
A class for handling words, derived from string.
static word groupName(Name name, const word &group)
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
BasicMomentumTransportModel::alphaField alphaField
const alphaField & alpha() const
Access function to phase fraction.
tmp< volScalarField > phaseTransferCoeff() const
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual tmp< fvScalarMatrix > kSource() const
dimensionedScalar pow3(const dimensionedScalar &ds)
const volVectorField & U() const
Access function to velocity field.
const Time & time() const
Return time.
BasicMomentumTransportModel::rhoField rhoField
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
BasicMomentumTransportModel::transportModel transportModel
tmp< volScalarField > bubbleG() const
virtual void correctNut()
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual tmp< fvScalarMatrix > epsilonSource() const
A class for managing temporary objects.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volScalarField > CdRe() const =0
Drag coefficient.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
const objectRegistry & db() const
Return the local objectRegistry.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
dimensionedScalar deltaT() const
Return time step.