40template<
class BasicTurbulenceModel>
41const IDDESDelta& SpalartAllmarasIDDES<BasicTurbulenceModel>::setDelta()
const
46 <<
"The delta function must be set to a " << IDDESDelta::typeName
54template<
class BasicTurbulenceModel>
55tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::alpha()
const
57 return max(0.25 - this->y_/IDDESDelta_.hmax(), scalar(-5));
61template<
class BasicTurbulenceModel>
62tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::ft
67 return tanh(
pow3(
sqr(Ct_)*this->r(this->nut_, magGradU, this->y_)));
71template<
class BasicTurbulenceModel>
72tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fl
77 return tanh(
pow(
sqr(Cl_)*this->r(this->
nu(), magGradU, this->y_), 10));
81template<
class BasicTurbulenceModel>
82tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fdt
93template<
class BasicTurbulenceModel>
123 fdTilda*(1 + fe)*lRAS + (1 - fdTilda)*lLES,
132 lerp(lLES, lRAS, fdTilda),
140template<
class BasicTurbulenceModel>
141SpalartAllmarasIDDES<BasicTurbulenceModel>::SpalartAllmarasIDDES
143 const alphaField&
alpha,
149 const word& propertiesName,
211 IDDESDelta_(setDelta())
215 this->printCoeffs(
type);
222template<
class BasicTurbulenceModel>
227 Cdt1_.readIfPresent(this->coeffDict());
228 Cdt2_.readIfPresent(this->coeffDict());
229 Cl_.readIfPresent(this->coeffDict());
230 Ct_.readIfPresent(this->coeffDict());
239template<
class BasicTurbulenceModel>
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
dimensionedScalar CDES_
DES coefficient.
virtual tmp< volScalarField > psi(const volScalarField &chi, const volScalarField &fv1) const
Return the low Reynolds number correction function.
virtual bool read()
Re-read model coefficients if they have changed.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > fd() const
Return the shielding function.
BasicTurbulenceModel::transportModel transportModel
const IDDESDelta & IDDESDelta_
IDDES delta.
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Return the length scale.
virtual bool read()
Re-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,...
Generic dimensioned Type class.
A class for managing temporary objects.
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
const volScalarField & psi
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Namespace for LES SGS models.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)