Go to the documentation of this file.
37 template<
class BasicTurbulenceModel>
40 if (!isA<IDDESDelta>(this->delta_()))
43 <<
"The delta function must be set to a " << IDDESDelta::typeName
47 return refCast<const IDDESDelta>(this->delta_());
51 template<
class BasicTurbulenceModel>
54 return max(0.25 - this->y_/IDDESDelta_.hmax(), scalar(-5));
58 template<
class BasicTurbulenceModel>
64 return tanh(
pow3(
sqr(Ct_)*rd(this->nut_, magGradU)));
68 template<
class BasicTurbulenceModel>
78 template<
class BasicTurbulenceModel>
96 *
sqr(this->kappa_*this->y_)
101 tr().boundaryField() == 0.0;
109 template<
class BasicTurbulenceModel>
115 return 1 -
tanh(
pow(Cdt1_*rd(this->nut_, magGradU), Cdt2_));
119 template<
class BasicTurbulenceModel>
153 fdTilda*(1 + fe)*lRAS + (1 - fdTilda)*lLES,
161 template<
class BasicTurbulenceModel>
170 const word& propertiesName,
222 IDDESDelta_(setDelta())
224 if (
type == typeName)
226 this->printCoeffs(
type);
233 template<
class BasicTurbulenceModel>
238 Cdt1_.readIfPresent(this->coeffDict());
239 Cdt2_.readIfPresent(this->coeffDict());
240 Cl_.readIfPresent(this->coeffDict());
241 Ct_.readIfPresent(this->coeffDict());
tmp< volScalarField > fl(const volScalarField &magGradU) const
A class for handling words, derived from string.
k-omega-SST DES turbulence model for incompressible and compressible flows
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
tmp< volScalarField > fdt(const volScalarField &magGradU) const
Delay function.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
dimensionedScalar exp(const dimensionedScalar &ds)
tmp< volScalarField > rd(const volScalarField &nur, const volScalarField &magGradU) const
virtual bool read()
Re-read model coefficients if they have changed.
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< volScalarField > ft(const volScalarField &magGradU) const
Generic dimensioned Type class.
errorManipArg< error, int > exit(error &err, const int errNo=1)
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
BasicTurbulenceModel::alphaField alphaField
tmp< volScalarField > alpha() const
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
const IDDESDelta & setDelta() const
Check that the supplied delta is an IDDESDelta.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
kOmegaSSTIDDES(const kOmegaSSTIDDES &)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar pos(const dimensionedScalar &ds)