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 RASDict_(this->subOrEmptyDict(
"RAS")),
70 turbulence_(RASDict_.getOrDefault<Switch>(
"turbulence", true)),
71 printCoeffs_(RASDict_.getOrDefault<Switch>(
"printCoeffs", false)),
72 coeffDict_(RASDict_.optionalSubDict(type +
"Coeffs")),
75 dimensioned<scalar>::getOrAddToDict
85 dimensioned<scalar>::getOrAddToDict
89 kMin_.dimensions()/dimTime,
95 dimensioned<scalar>::getOrAddToDict
106 this->mesh_.deltaCoeffs();
112template<
class BasicTurbulenceModel>
116 const alphaField&
alpha,
122 const word& propertiesName
143 dict.getCompat<
word>(
"model", {{
"RASModel", -2006}})
146 Info<<
"Selecting RAS turbulence model " << modelType <<
endl;
148 auto* ctorPtr = dictionaryConstructorTable(modelType);
157 *dictionaryConstructorTablePtr_
161 return autoPtr<RASModel>
163 ctorPtr(
alpha,
rho,
U, alphaRhoPhi,
phi, transport, propertiesName)
170template<
class BasicTurbulenceModel>
173 if (BasicTurbulenceModel::read())
175 RASDict_ <<= this->subDict(
"RAS");
176 RASDict_.readEntry(
"turbulence", turbulence_);
178 coeffDict_ <<= RASDict_.optionalSubDict(
type() +
"Coeffs");
180 kMin_.readIfPresent(RASDict_);
181 epsilonMin_.readIfPresent(RASDict_);
182 omegaMin_.readIfPresent(RASDict_);
190template<
class BasicTurbulenceModel>
194 const scalar Cmu = 0.09;
200 Cmu*this->k()*this->omega()
205template<
class BasicTurbulenceModel>
209 const scalar betaStar = 0.09;
216 this->epsilon()/(betaStar*(this->k() + k0))
221template<
class BasicTurbulenceModel>
224 BasicTurbulenceModel::correct();
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
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.
virtual void printCoeffs(const word &type)
Print model coefficients.
BasicTurbulenceModel::rhoField rhoField
RASModel(const RASModel &)=delete
No copy construct.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dictionary RASDict_
RAS coefficients dictionary.
dimensionedScalar omegaMin_
Lower limit for omega.
dimensionedScalar epsilonMin_
Lower limit of epsilon.
static autoPtr< RASModel > 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 RAS model.
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.
Switch printCoeffs_
Flag to print the model coeffs at run-time.
virtual bool read()
Read model coefficients if they have changed.
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,...
Generic dimensioned Type class.
A class for managing temporary objects.
Base-class for all transport models used by the incompressible turbulence models.
static const word propertiesName
Default name of the turbulence properties dictionary.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
const surfaceScalarField & alphaRhoPhi_
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.
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)