Go to the documentation of this file.
41 template<
class BasicTurbulenceModel>
50 const word& propertiesName,
66 gasTurbulencePtr_(
nullptr),
80 this->printCoeffs(
type);
87 template<
class BasicTurbulenceModel>
92 Cmub_.readIfPresent(this->coeffDict());
101 template<
class BasicTurbulenceModel>
102 const PhaseCompressibleTurbulenceModel
104 typename BasicTurbulenceModel::transportModel
106 SmagorinskyZhang<BasicTurbulenceModel>::gasTurbulence()
const
108 if (!gasTurbulencePtr_)
112 const transportModel& liquid = this->transport();
113 const twoPhaseSystem&
fluid =
114 refCast<const twoPhaseSystem>(liquid.fluid());
115 const transportModel& gas =
fluid.otherPhase(liquid);
119 .lookupObject<PhaseCompressibleTurbulenceModel<transportModel>>
129 return *gasTurbulencePtr_;
133 template<
class BasicTurbulenceModel>
137 this->gasTurbulence();
143 + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
144 *(
mag(this->U_ - gasTurbulence.U()));
146 this->nut_.correctBoundaryConditions();
149 BasicTurbulenceModel::correctNut();
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
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)
BasicTurbulenceModel::rhoField rhoField
Templated abstract base class for multiphase compressible turbulence models.
static options & New(const fvMesh &mesh)
const dimensionedScalar alpha
static const word propertiesName
BasicTurbulenceModel::alphaField alphaField
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
The Smagorinsky SGS model including bubble-generated turbulence.
Generic dimensioned Type class.
Base-class for all transport models used by the incompressible turbulence models.
GeometricField< vector, fvPatchField, volMesh > volVectorField
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)
The Smagorinsky SGS model.
Generic GeometricField class.
BasicTurbulenceModel::transportModel transportModel
virtual void correctNut()