41 template<
class BasicMomentumTransportModel>
53 kEpsilon<BasicMomentumTransportModel>
64 gasTurbulencePtr_(nullptr),
106 if (
type == typeName)
108 this->printCoeffs(
type);
115 template<
class BasicMomentumTransportModel>
120 alphaInversion_.readIfPresent(this->coeffDict());
121 Cp_.readIfPresent(this->coeffDict());
122 C4_.readIfPresent(this->coeffDict());
123 Cmub_.readIfPresent(this->coeffDict());
134 template<
class BasicMomentumTransportModel>
138 if (!gasTurbulencePtr_)
143 refCast<const phaseModel>(this->properties());
155 momentumTransportModel::typeName,
161 return *gasTurbulencePtr_;
165 template<
class BasicMomentumTransportModel>
169 this->gasTurbulence();
176 this->Cmu_*
sqr(this->k_)/this->epsilon_
177 + Cmub_*gas.
d()*gasTurbulence.
alpha()
178 *(
mag(this->U_ - gasTurbulence.
U()));
180 this->nut_.correctBoundaryConditions();
185 template<
class BasicMomentumTransportModel>
189 this->gasTurbulence();
217 template<
class BasicMomentumTransportModel>
229 max(alphaInversion_ -
alpha, scalar(0))
231 *
min(gasTurbulence.
epsilon()/gasTurbulence.
k(), 1.0/
U.time().deltaT())
236 template<
class BasicMomentumTransportModel>
243 this->gasTurbulence();
245 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
249 + phaseTransferCoeff*gasTurbulence.
k()
250 -
fvm::Sp(phaseTransferCoeff, this->k_);
254 template<
class BasicMomentumTransportModel>
262 this->gasTurbulence();
264 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
267 alpha*
rho*this->C4_*this->epsilon_*bubbleG()/this->k_
268 + phaseTransferCoeff*gasTurbulence.
epsilon()
269 -
fvm::Sp(phaseTransferCoeff, this->epsilon_);
273 template<
class BasicMomentumTransportModel>
Generic GeometricField class.
static word groupName(Name name, const word &group)
Continuous-phase k-epsilon model including bubble-generated turbulence.
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
LaheyKEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
tmp< volScalarField > bubbleG() const
tmp< volScalarField > phaseTransferCoeff() const
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
virtual bool read()
Read model coefficients if they have changed.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Generic dimensioned Type class.
Class to represent a interface between phases where one phase is considered dispersed within the othe...
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Templated abstract base class for multiphase compressible turbulence models.
const alphaField & alpha() const
Access function to phase fraction.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
const word & name() const
Return the name of this phase.
Class to represent a system of phases and model interfacial transfers between them.
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
A class for managing temporary objects.
Abstract base class for all fluid physical properties.
A class for handling words, derived from string.
Info<< "Reading strained laminar flame speed field Su\n"<< endl;volScalarField Su(IOobject("Su", runTime.name(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\n"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\n"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\n"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\n"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\n"<< 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
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.