Go to the documentation of this file.
30 template<
class BasicTurbulenceModel>
42 template<
class BasicTurbulenceModel>
52 const word& propertiesName
67 LESDict_(this->subOrEmptyDict(
"LES")),
68 turbulence_(LESDict_.lookup(
"turbulence")),
69 printCoeffs_(LESDict_.lookupOrDefault<
Switch>(
"printCoeffs",
false)),
70 coeffDict_(LESDict_.subOrEmptyDict(
type +
"Coeffs")),
98 IOobject::groupName(
"delta",
U.group()),
106 this->mesh_.deltaCoeffs();
112 template<
class BasicTurbulenceModel>
122 const word& propertiesName
133 IOobject::groupName(propertiesName,
U.group()),
136 IOobject::MUST_READ_IF_MODIFIED,
143 Info<<
"Selecting LES turbulence model " << modelType <<
endl;
145 typename dictionaryConstructorTable::iterator cstrIter =
146 dictionaryConstructorTablePtr_->find(modelType);
148 if (cstrIter == dictionaryConstructorTablePtr_->end())
151 <<
"Unknown LESModel type "
152 << modelType <<
nl <<
nl
153 <<
"Valid LESModel types:" <<
endl
154 << dictionaryConstructorTablePtr_->sortedToc()
160 cstrIter()(
alpha,
rho,
U, alphaRhoPhi,
phi, transport, propertiesName)
167 template<
class BasicTurbulenceModel>
172 LESDict_ <<= this->subDict(
"LES");
173 LESDict_.lookup(
"turbulence") >> turbulence_;
177 coeffDict_ <<= *dictPtr;
180 delta_().read(LESDict_);
182 kMin_.readIfPresent(LESDict_);
193 template<
class BasicTurbulenceModel>
BasicTurbulenceModel::rhoField rhoField
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
const dimensionSet dimVelocity
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool read(const char *, int32_t &)
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.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
BasicTurbulenceModel::transportModel transportModel
LESModel(const LESModel &)
Disallow default bitwise copy construct.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
virtual bool read()
Read model coefficients if they have changed.
virtual void printCoeffs(const word &type)
Print model coefficients.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Generic dimensioned Type class.
errorManipArg< error, int > exit(error &err, const int errNo=1)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
BasicTurbulenceModel::alphaField alphaField
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Generic GeometricField class.
const dictionary * subDictPtr(const word &) const
Find and return a sub-dictionary pointer if present.