Go to the documentation of this file.
30 #include "twoPhaseSystem.H"
41 template<
class BasicTurbulenceModel>
42 continuousGasKEqn<BasicTurbulenceModel>::continuousGasKEqn
50 const word& propertiesName,
66 liquidTurbulencePtr_(
nullptr),
80 this->printCoeffs(
type);
87 template<
class BasicTurbulenceModel>
92 alphaInversion_.readIfPresent(this->coeffDict());
101 template<
class BasicTurbulenceModel>
105 if (!liquidTurbulencePtr_)
111 refCast<const twoPhaseSystem>(gas.fluid());
114 liquidTurbulencePtr_ =
125 return *liquidTurbulencePtr_;
129 template<
class BasicTurbulenceModel>
141 max(alphaInversion_ -
alpha, scalar(0))
145 this->Ce_*
sqrt(liquidTurbulence.
k())/this->
delta(),
146 1.0/U.time().deltaT()
152 template<
class BasicTurbulenceModel>
157 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
160 phaseTransferCoeff*liquidTurbulence.
k()
161 -
fvm::Sp(phaseTransferCoeff, this->k_);
BasicTurbulenceModel::alphaField alphaField
A class for handling words, derived from Foam::string.
Class which solves the volume fraction equations for two phases.
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< fvScalarMatrix > kSource() const
static const word propertiesName
Default name of the turbulence properties dictionary.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
virtual bool read()
Read model coefficients if they have changed.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
One equation eddy-viscosity model.
tmp< volScalarField > phaseTransferCoeff() const
const turbulenceModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Abstract base class for turbulence models (RAS, LES and laminar).
Generic dimensioned Type class.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
BasicTurbulenceModel::rhoField rhoField
BasicTurbulenceModel::transportModel transportModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.