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),
99 this->printCoeffs(
type);
106 template<
class BasicTurbulenceModel>
111 alphaInversion_.readIfPresent(this->coeffDict());
112 Cp_.readIfPresent(this->coeffDict());
113 Cmub_.readIfPresent(this->coeffDict());
124 template<
class BasicTurbulenceModel>
127 typename BasicTurbulenceModel::transportModel
131 if (!gasTurbulencePtr_)
137 refCast<const twoPhaseSystem>(liquid.fluid());
152 return *gasTurbulencePtr_;
156 template<
class BasicTurbulenceModel>
160 this->gasTurbulence();
165 *(
mag(this->U_ - gasTurbulence.
U()));
167 this->nut_.correctBoundaryConditions();
171 template<
class BasicTurbulenceModel>
175 this->gasTurbulence();
179 refCast<const twoPhaseSystem>(liquid.fluid());
186 Cp_*
sqr(magUr)*
fluid.drag(gas).K()/liquid.rho()
193 template<
class BasicTurbulenceModel>
205 max(alphaInversion_ -
alpha, scalar(0))
209 this->Ce_*
sqrt(gasTurbulence.
k())/this->delta(),
210 1.0/
U.time().deltaT()
216 template<
class BasicTurbulenceModel>
223 this->gasTurbulence();
225 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
229 + phaseTransferCoeff*gasTurbulence.
k()
230 -
fvm::Sp(phaseTransferCoeff, this->k_);
static word groupName(Name name, const word &group)
A class for handling words, derived from string.
Class which solves the volume fraction equations for two phases.
Templated abstract base class for multiphase compressible turbulence models.
A class for managing temporary objects.
tmp< volScalarField > phaseTransferCoeff() const
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
static const word propertiesName
Default name of the turbulence properties dictionary.
BasicTurbulenceModel::rhoField rhoField
const alphaField & alpha() const
Access function to phase fraction.
dimensioned< scalar > mag(const dimensioned< Type > &)
One equation eddy-viscosity model.
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::alphaField alphaField
virtual bool read()
Read model coefficients if they have changed.
tmp< volScalarField > bubbleG() const
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 correctNut()
NicenoKEqn(const NicenoKEqn &)
const PhaseCompressibleTurbulenceModel< typename BasicTurbulenceModel::transportModel > & gasTurbulence() const
Return the turbulence model for the gas phase.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
BasicTurbulenceModel::transportModel transportModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const volVectorField & U() const
Access function to velocity field.