33template<
class BasicTurbulenceModel>
45template<
class BasicTurbulenceModel>
49 const alphaField&
alpha,
51 const volVectorField&
U,
52 const surfaceScalarField& alphaRhoPhi,
53 const surfaceScalarField&
phi,
54 const transportModel& transport,
55 const word& propertiesName
69 LESDict_(this->subOrEmptyDict(
"LES")),
70 turbulence_(LESDict_.getOrDefault<Switch>(
"turbulence", true)),
71 printCoeffs_(LESDict_.getOrDefault<Switch>(
"printCoeffs", false)),
72 coeffDict_(LESDict_.optionalSubDict(type +
"Coeffs")),
75 dimensioned<scalar>::getOrAddToDict
84 dimensioned<scalar>::getOrAddToDict
94 dimensioned<scalar>::getOrAddToDict
98 kMin_.dimensions()/dimTime,
104 dimensioned<scalar>::getOrAddToDict
116 IOobject::groupName(
"delta", alphaRhoPhi.group()),
124 this->mesh_.deltaCoeffs();
130template<
class BasicTurbulenceModel>
134 const alphaField&
alpha,
140 const word& propertiesName
161 dict.getCompat<
word>(
"model", {{
"LESModel", -2006}})
164 Info<<
"Selecting LES turbulence model " << modelType <<
endl;
166 auto* ctorPtr = dictionaryConstructorTable(modelType);
175 *dictionaryConstructorTablePtr_
181 ctorPtr(
alpha,
rho,
U, alphaRhoPhi,
phi, transport, propertiesName)
188template<
class BasicTurbulenceModel>
191 if (BasicTurbulenceModel::read())
193 LESDict_ <<= this->subDict(
"LES");
194 LESDict_.readEntry(
"turbulence", turbulence_);
196 coeffDict_ <<= LESDict_.optionalSubDict(
type() +
"Coeffs");
198 delta_().read(LESDict_);
200 Ce_.readIfPresent(LESDict_);
202 kMin_.readIfPresent(LESDict_);
211template<
class BasicTurbulenceModel>
220 this->mesh_.time().timeName(),
228template<
class BasicTurbulenceModel>
232 const scalar betaStar = 0.09;
240 this->mesh_.time().timeName(),
243 this->epsilon()/(betaStar*(this->k() + k0))
248template<
class BasicTurbulenceModel>
252 BasicTurbulenceModel::correct();
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word group(const word &name)
Return group (extension part of name).
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
BasicTurbulenceModel::alphaField alphaField
dictionary coeffDict_
Model coefficients dictionary.
static autoPtr< LESModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected LES model.
virtual void printCoeffs(const word &type)
Print model coefficients.
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar Ce_
Empirical model constant.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar omegaMin_
Lower limit for omega.
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
dimensionedScalar epsilonMin_
Lower limit of epsilon.
dictionary LESDict_
LES coefficients dictionary.
Switch turbulence_
Turbulence on/off flag.
dimensionedScalar kMin_
Lower limit of k.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
LESModel(const LESModel &)=delete
No copy construct.
Switch printCoeffs_
Flag to print the model coeffs at run-time.
virtual bool read()
Read model coefficients if they have changed.
Abstract base class for LES deltas.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Generic dimensioned Type class.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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)
Ostream & endl(Ostream &os)
Add newline and flush stream.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)