Go to the documentation of this file.
40 template<
class BasicTurbulenceModel>
41 const IDDESDelta& SpalartAllmarasIDDES<BasicTurbulenceModel>::setDelta()
const
43 if (!isA<IDDESDelta>(this->delta_()))
46 <<
"The delta function must be set to a " << IDDESDelta::typeName
50 return refCast<const IDDESDelta>(this->delta_());
54 template<
class BasicTurbulenceModel>
55 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::alpha()
const
57 return max(0.25 - this->y_/IDDESDelta_.hmax(), scalar(-5));
61 template<
class BasicTurbulenceModel>
62 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::ft
67 return tanh(
pow3(
sqr(Ct_)*rd(this->nut_, magGradU)));
71 template<
class BasicTurbulenceModel>
72 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fl
81 template<
class BasicTurbulenceModel>
82 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::rd
88 tmp<volScalarField>
tr
99 *
sqr(this->kappa_*this->y_)
104 tr.ref().boundaryFieldRef() == 0.0;
110 template<
class BasicTurbulenceModel>
111 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fdt
116 return 1 -
tanh(
pow(Cdt1_*rd(this->nut_, magGradU), Cdt2_));
122 template<
class BasicTurbulenceModel>
157 fdTilda*(1 + fe)*lRAS + (1 - fdTilda)*lLES,
165 template<
class BasicTurbulenceModel>
168 const alphaField&
alpha,
174 const word& propertiesName,
226 IDDESDelta_(setDelta())
228 if (
type == typeName)
237 template<
class BasicTurbulenceModel>
242 Cdt1_.readIfPresent(this->coeffDict());
243 Cdt2_.readIfPresent(this->coeffDict());
244 Cl_.readIfPresent(this->coeffDict());
245 Ct_.readIfPresent(this->coeffDict());
SpalartAllmaras IDDES turbulence model for incompressible and compressible flows.
A class for handling words, derived from Foam::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
const dimensionedScalar alpha
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
IDDESDelta used by the IDDES (improved low Re Spalart-Allmaras DES model) The min and max delta are c...
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual void printCoeffs(const word &type)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
Generic dimensioned Type class.
BasicTurbulenceModel::transportModel transportModel
errorManipArg< error, int > exit(error &err, const int errNo=1)
Base-class for all transport models used by the incompressible turbulence models.
#define FatalErrorInFunction
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
BasicTurbulenceModel::alphaField alphaField
Generic GeometricField class.
const volScalarField & psi
dimensionedScalar neg(const dimensionedScalar &ds)