41template<
class BasicTurbulenceModel>
47 return dev(
symm(gradU & gradU));
51template<
class BasicTurbulenceModel>
66 this->runTime_.timeName(),
76 +
pow(magSqrSd, 5.0/4.0)
91template<
class BasicTurbulenceModel>
95 this->nut_.correctBoundaryConditions();
98 BasicTurbulenceModel::correctNut();
104template<
class BasicTurbulenceModel>
113 const word& propertiesName,
149 if (
type == typeName)
151 this->printCoeffs(
type);
158template<
class BasicTurbulenceModel>
163 Ck_.readIfPresent(this->coeffDict());
164 Cw_.readIfPresent(this->coeffDict());
173template<
class BasicTurbulenceModel>
176 if (!this->turbulence_)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Eddy viscosity LES SGS model base class.
The Wall-adapting local eddy-viscosity (WALE) SGS model.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Correct Eddy-Viscosity and related properties.
virtual void correctNut()
Update the SGS eddy-viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
virtual bool read()
Read model coefficients if they have changed.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
vector Sd() const
Face area normal for side d.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)