Go to the documentation of this file.
31 #include "twoPhaseSystem.H"
32 #include "dragModel.H"
43 template<
class BasicTurbulenceModel>
44 LaheyKEpsilon<BasicTurbulenceModel>::LaheyKEpsilon
52 const word& propertiesName,
68 gasTurbulencePtr_(
nullptr),
110 if (
type == typeName)
112 this->printCoeffs(
type);
119 template<
class BasicTurbulenceModel>
124 alphaInversion_.readIfPresent(this->coeffDict());
125 Cp_.readIfPresent(this->coeffDict());
126 C3_.readIfPresent(this->coeffDict());
127 Cmub_.readIfPresent(this->coeffDict());
136 template<
class BasicTurbulenceModel>
139 typename BasicTurbulenceModel::transportModel
143 if (!gasTurbulencePtr_)
149 refCast<const twoPhaseSystem>(
liquid.fluid());
164 return *gasTurbulencePtr_;
168 template<
class BasicTurbulenceModel>
172 this->gasTurbulence();
175 this->Cmu_*
sqr(this->k_)/this->epsilon_
176 + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
177 *(
mag(this->U_ - gasTurbulence.U()));
179 this->nut_.correctBoundaryConditions();
182 BasicTurbulenceModel::correctNut();
186 template<
class BasicTurbulenceModel>
190 this->gasTurbulence();
216 template<
class BasicTurbulenceModel>
228 max(alphaInversion_ -
alpha, scalar(0))
230 *
min(gasTurbulence.
epsilon()/gasTurbulence.
k(), 1.0/
U.time().deltaT())
235 template<
class BasicTurbulenceModel>
242 this->gasTurbulence();
244 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
248 + phaseTransferCoeff*gasTurbulence.k()
249 -
fvm::Sp(phaseTransferCoeff, this->k_);
253 template<
class BasicTurbulenceModel>
260 this->gasTurbulence();
262 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
265 alpha*
rho*this->C3_*this->epsilon_*bubbleG()/this->k_
266 + phaseTransferCoeff*gasTurbulence.epsilon()
267 -
fvm::Sp(phaseTransferCoeff, this->epsilon_);
271 template<
class BasicTurbulenceModel>
Info<< "Reading strained laminar flame speed field Su\n"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), 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
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
BasicTurbulenceModel::transportModel transportModel
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.
tmp< volScalarField > phaseTransferCoeff() const
Templated abstract base class for multiphase compressible turbulence models.
A class for managing temporary objects.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
static const word propertiesName
Default name of the turbulence properties dictionary.
tmp< volScalarField > bubbleG() const
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
virtual tmp< fvScalarMatrix > kSource() const
virtual void correctNut()
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensionedScalar pow3(const dimensionedScalar &ds)
Continuous-phase k-epsilon model including bubble-generated turbulence.
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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).
Generic dimensioned Type class.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Base-class for all transport models used by the incompressible turbulence models.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual bool read()
Read model coefficients if they have changed.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
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< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
BasicTurbulenceModel::alphaField alphaField