40template<
class BasicTurbulenceModel>
41const IDDESDelta& kOmegaSSTIDDES<BasicTurbulenceModel>::setDelta()
const
43 if (!isA<IDDESDelta>(this->delta_()))
50 return refCast<const IDDESDelta>(this->delta_());
54template<
class BasicTurbulenceModel>
55tmp<volScalarField> kOmegaSSTIDDES<BasicTurbulenceModel>::alpha()
const
57 return max(0.25 - this->y_/IDDESDelta_.hmax(), scalar(-5));
61template<
class BasicTurbulenceModel>
62tmp<volScalarField> kOmegaSSTIDDES<BasicTurbulenceModel>::ft
67 return tanh(
pow3(
sqr(Ct_)*rd(this->nut_, magGradU)));
71template<
class BasicTurbulenceModel>
72tmp<volScalarField> kOmegaSSTIDDES<BasicTurbulenceModel>::fl
81template<
class BasicTurbulenceModel>
82tmp<volScalarField> kOmegaSSTIDDES<BasicTurbulenceModel>::rd
88 tmp<volScalarField>
tr
99 *
sqr(this->kappa_*this->y_)
104 tr.ref().boundaryFieldRef() == 0.0;
112template<
class BasicTurbulenceModel>
113tmp<volScalarField> kOmegaSSTIDDES<BasicTurbulenceModel>::fdt
118 return 1 -
tanh(
pow(Cdt1_*rd(this->nut_, magGradU), Cdt2_));
122template<
class BasicTurbulenceModel>
156 fdTilda*(1 + fe)*lRAS + (1 - fdTilda)*lLES,
164template<
class BasicTurbulenceModel>
173 const word& propertiesName,
225 IDDESDelta_(setDelta())
227 if (
type == typeName)
229 this->printCoeffs(
type);
236template<
class BasicTurbulenceModel>
241 Cdt1_.readIfPresent(this->coeffDict());
242 Cdt2_.readIfPresent(this->coeffDict());
243 Cl_.readIfPresent(this->coeffDict());
244 Ct_.readIfPresent(this->coeffDict());
IDDESDelta used by the IDDES (improved low Re Spalart-Allmaras DES model) The min and max delta are c...
k-omega-SST DES turbulence model for incompressible and compressible flows
k-omega-SST IDDES turbulence model for incompressible and compressible flows
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Re-read model coefficients if they have changed.
Generic dimensioned Type class.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
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)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
dimensionedScalar neg(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
static const char *const typeName
The type name used in ensight case files.