117 template<
class BasicTurbulenceModel>
176 typedef typename BasicTurbulenceModel::alphaField
alphaField;
177 typedef typename BasicTurbulenceModel::rhoField
rhoField;
178 typedef typename BasicTurbulenceModel::transportModel
transportModel;
190 const alphaField&
alpha,
195 const transportModel& transport,
219 this->
nut_/sigmaK_ + this->
nu()
232 this->
nut_/sigmaEps_ + this->
nu()
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
dimensionedScalar sigmaEps_
volScalarField k_
Turbulence kinetic energy.
Eddy viscosity turbulence model base class.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
tmp< volScalarField > Ls() const
Return length scale, Ls.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields. ...
volScalarField v2_
Turbulence stress normal to streamlines.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
virtual ~v2f()
Destructor.
virtual tmp< volScalarField > v2() const
Return turbulence stress normal to streamlines.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar CmuKEps_
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
v2f(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
virtual void correctNut()
BasicTurbulenceModel::transportModel transportModel
volScalarField f_
Damping function.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
TypeName("v2f")
Runtime type information.
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > f() const
Return the damping function.
Lien and Kalitzin's v2-f turbulence model for incompressible and compressible flows, with a limit imposed on the turbulent viscosity given by Davidson et al.
virtual bool read()
Read RASProperties dictionary.
BasicTurbulenceModel::alphaField alphaField
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField epsilon_
Turbulence dissipation.
dimensionedScalar sigmaK_
tmp< volScalarField > Ts() const
Return time scale, Ts.