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>
154 fdTilda*(1 + fe)*lRAS + (1 - fdTilda)*lLES,
162 template<
class BasicTurbulenceModel>
171 const word& propertiesName,
223 IDDESDelta_(setDelta())
225 if (
type == typeName)
227 this->printCoeffs(
type);
234 template<
class BasicTurbulenceModel>
239 Cdt1_.readIfPresent(this->coeffDict());
240 Cdt2_.readIfPresent(this->coeffDict());
241 Cl_.readIfPresent(this->coeffDict());
242 Ct_.readIfPresent(this->coeffDict());
A class for handling words, derived from string.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< volScalarField > alpha() const
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
BasicTurbulenceModel::rhoField rhoField
const IDDESDelta & setDelta() const
Check that the supplied delta is an IDDESDelta.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > fl(const volScalarField &magGradU) const
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)
Generic dimensioned Type class.
BasicTurbulenceModel::transportModel transportModel
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const volScalarField & psi
SpalartAllmarasIDDES(const SpalartAllmarasIDDES &)
tmp< volScalarField > fdt(const volScalarField &magGradU) const
Delay function.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
BasicTurbulenceModel::alphaField alphaField
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar neg(const dimensionedScalar &ds)
tmp< volScalarField > ft(const volScalarField &magGradU) const
tmp< volScalarField > rd(const volScalarField &nur, const volScalarField &magGradU) const
dimensionedScalar pos(const dimensionedScalar &ds)