Go to the documentation of this file.
31 #include "twoPhaseSystem.H"
32 #include "dragModel.H"
43 template<
class BasicTurbulenceModel>
52 const word& propertiesName,
68 gasTurbulencePtr_(
nullptr),
100 if (
type == typeName)
102 this->printCoeffs(
type);
109 template<
class BasicTurbulenceModel>
114 alphaInversion_.readIfPresent(this->coeffDict());
115 Cp_.readIfPresent(this->coeffDict());
116 Cmub_.readIfPresent(this->coeffDict());
125 template<
class BasicTurbulenceModel>
126 const PhaseCompressibleTurbulenceModel
128 typename BasicTurbulenceModel::transportModel
130 NicenoKEqn<BasicTurbulenceModel>::gasTurbulence()
const
132 if (!gasTurbulencePtr_)
136 const transportModel& liquid = this->transport();
137 const twoPhaseSystem&
fluid =
138 refCast<const twoPhaseSystem>(liquid.fluid());
139 const transportModel& gas =
fluid.otherPhase(liquid);
143 .lookupObject<PhaseCompressibleTurbulenceModel<transportModel>>
153 return *gasTurbulencePtr_;
157 template<
class BasicTurbulenceModel>
161 this->gasTurbulence();
165 + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
166 *(
mag(this->U_ - gasTurbulence.U()));
168 this->nut_.correctBoundaryConditions();
171 BasicTurbulenceModel::correctNut();
175 template<
class BasicTurbulenceModel>
179 this->gasTurbulence();
183 refCast<const twoPhaseSystem>(
liquid.fluid());
199 template<
class BasicTurbulenceModel>
204 const alphaField&
alpha = this->alpha_;
205 const rhoField&
rho = this->rho_;
211 max(alphaInversion_ -
alpha, scalar(0))
215 this->Ce_*
sqrt(gasTurbulence.
k())/this->delta(),
216 1.0/
U.time().deltaT()
222 template<
class BasicTurbulenceModel>
225 const alphaField&
alpha = this->alpha_;
226 const rhoField&
rho = this->rho_;
229 this->gasTurbulence();
231 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
235 + phaseTransferCoeff*gasTurbulence.k()
236 -
fvm::Sp(phaseTransferCoeff, this->k_);
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
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)
Templated abstract base class for multiphase compressible turbulence models.
A class for managing temporary objects.
tmp< volScalarField > phaseTransferCoeff() const
static options & New(const fvMesh &mesh)
const dimensionedScalar alpha
static const word propertiesName
BasicTurbulenceModel::rhoField rhoField
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
label min(const labelHashSet &set, label minValue=labelMax)
One equation eddy-viscosity model.
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::alphaField alphaField
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
One-equation SGS model for the continuous phase in a two-phase system including bubble-generated turb...
label max(const labelHashSet &set, label maxValue=labelMin)
tmp< volScalarField > bubbleG() const
Abstract base class for turbulence models (RAS, LES and laminar).
Generic dimensioned Type class.
virtual tmp< volScalarField > k() const =0
virtual void correctNut()
Base-class for all transport models used by the incompressible turbulence models.
GeometricField< vector, fvPatchField, volMesh > volVectorField
scalar rho(scalar p, scalar T) const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
BasicTurbulenceModel::transportModel transportModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
static word groupName(StringType base, const word &group)
Generic GeometricField class.