Go to the documentation of this file.
28 #include "twoPhaseSystem.H"
29 #include "dragModel.H"
40 template<
class BasicTurbulenceModel>
49 const word& propertiesName,
65 gasTurbulencePtr_(NULL),
107 if (
type == typeName)
109 this->printCoeffs(
type);
116 template<
class BasicTurbulenceModel>
121 alphaInversion_.readIfPresent(this->coeffDict());
122 Cp_.readIfPresent(this->coeffDict());
123 C3_.readIfPresent(this->coeffDict());
124 Cmub_.readIfPresent(this->coeffDict());
135 template<
class BasicTurbulenceModel>
138 typename BasicTurbulenceModel::transportModel
142 if (!gasTurbulencePtr_)
148 refCast<const twoPhaseSystem>(liquid.fluid());
163 return *gasTurbulencePtr_;
167 template<
class BasicTurbulenceModel>
171 this->gasTurbulence();
174 this->Cmu_*
sqr(this->k_)/this->epsilon_
176 *(
mag(this->U_ - gasTurbulence.
U()));
178 this->nut_.correctBoundaryConditions();
182 template<
class BasicTurbulenceModel>
186 this->gasTurbulence();
199 +
pow(
fluid.drag(gas).CdRe()*liquid.nu()/gas.d(), 4.0/3.0)
210 template<
class BasicTurbulenceModel>
222 max(alphaInversion_ -
alpha, scalar(0))
224 *
min(gasTurbulence.
epsilon()/gasTurbulence.
k(), 1.0/
U.time().deltaT())
229 template<
class BasicTurbulenceModel>
236 this->gasTurbulence();
238 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
242 + phaseTransferCoeff*gasTurbulence.
k()
243 -
fvm::Sp(phaseTransferCoeff, this->k_);
247 template<
class BasicTurbulenceModel>
254 this->gasTurbulence();
256 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
259 alpha*
rho*this->C3_*this->epsilon_*bubbleG()/this->k_
260 + phaseTransferCoeff*gasTurbulence.
epsilon()
261 -
fvm::Sp(phaseTransferCoeff, this->epsilon_);
265 template<
class BasicTurbulenceModel>
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static word groupName(Name name, const word &group)
BasicTurbulenceModel::transportModel transportModel
A class for handling words, derived from string.
Class which solves the volume fraction equations for two phases.
tmp< volScalarField > phaseTransferCoeff() const
Templated abstract base class for multiphase compressible turbulence models.
A class for managing temporary objects.
const PhaseCompressibleTurbulenceModel< typename BasicTurbulenceModel::transportModel > & gasTurbulence() const
Return the turbulence model for the gas phase.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
static const word propertiesName
Default name of the turbulence properties dictionary.
tmp< volScalarField > bubbleG() const
const alphaField & alpha() const
Access function to phase fraction.
virtual tmp< fvScalarMatrix > kSource() const
virtual void correctNut()
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensionedScalar pow3(const dimensionedScalar &ds)
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const transportModel & transport() const
Access function to incompressible transport model.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual bool read()
Read model coefficients if they have changed.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
LaheyKEpsilon(const LaheyKEpsilon &)
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< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const volVectorField & U() const
Access function to velocity field.
BasicTurbulenceModel::alphaField alphaField