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),
79 this->printCoeffs(
type);
86 template<
class BasicTurbulenceModel>
91 Cmub_.readIfPresent(this->coeffDict());
101 template<
class BasicTurbulenceModel>
104 typename BasicTurbulenceModel::transportModel
108 if (!gasTurbulencePtr_)
114 refCast<const twoPhaseSystem>(liquid.fluid());
129 return *gasTurbulencePtr_;
133 template<
class BasicTurbulenceModel>
137 this->gasTurbulence();
141 pow(this->betaStar_, 0.25)*this->y_*
sqrt(this->k_)/this->
nu()
148 this->a1_*this->omega_,
153 *(
mag(this->U_ - gasTurbulence.
U()));
155 this->nut_.correctBoundaryConditions();
159 template<
class BasicTurbulenceModel>
static word groupName(Name name, const word &group)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
BasicTurbulenceModel::rhoField rhoField
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from string.
Class which solves the volume fraction equations for two phases.
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows.
Templated abstract base class for multiphase compressible turbulence models.
kOmegaSSTSato(const kOmegaSSTSato &)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual bool read()
Read model coefficients if they have changed.
static const word propertiesName
Default name of the turbulence properties dictionary.
const alphaField & alpha() const
Access function to phase fraction.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
BasicTurbulenceModel::transportModel transportModel
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
const transportModel & transport() const
Access function to incompressible transport model.
virtual void correctNut()
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
BasicTurbulenceModel::alphaField alphaField
const PhaseCompressibleTurbulenceModel< typename BasicTurbulenceModel::transportModel > & gasTurbulence() const
Return the turbulence model for the gas phase.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
const volVectorField & U() const
Access function to velocity field.