37 template<
class BasicTurbulenceModel>
38 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::alpha()
const 42 0.25 - this->y_/static_cast<const volScalarField&>(IDDESDelta_.hmax()),
48 template<
class BasicTurbulenceModel>
49 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::ft
54 return tanh(
pow3(
sqr(ct_)*rd(this->nut_, magGradU)));
58 template<
class BasicTurbulenceModel>
59 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fl
68 template<
class BasicTurbulenceModel>
69 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::rd
83 )*
sqr(this->kappa_*this->y_)
90 template<
class BasicTurbulenceModel>
91 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fd
96 return 1 -
tanh(
pow3(8*rd(this->nuEff(), magGradU)));
102 template<
class BasicTurbulenceModel>
115 2*(
pos0(alpha)*
pow(expTerm, -11.09) +
neg(alpha)*
pow(expTerm, -9.0));
132 - this->Cb1_*this->fv2(chi, fv1)
133 /(this->Cw1_*
sqr(this->kappa_)*fwStar_)
142 fHyb*(1 + fRestore*Psi)*this->y_
143 + (1 - fHyb)*this->CDES_*Psi*this->
delta()
150 template<
class BasicTurbulenceModel>
159 const word& propertiesName,
200 IDDESDelta_(refCast<IDDESDelta>(this->delta_()))
206 template<
class BasicTurbulenceModel>
211 fwStar_.readIfPresent(this->coeffDict());
212 cl_.readIfPresent(this->coeffDict());
213 ct_.readIfPresent(this->coeffDict());
dimensionedScalar tanh(const dimensionedScalar &ds)
BasicTurbulenceModel::transportModel transportModel
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
dimensionedScalar sqrt(const dimensionedScalar &ds)
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Generic dimensioned Type class.
virtual bool read()
Read model coefficients if they have changed.
dimensionedScalar neg(const dimensionedScalar &ds)
BasicTurbulenceModel::alphaField alphaField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar exp(const dimensionedScalar &ds)
SpalartAllmarasIDDES(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
A class for handling words, derived from string.
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
BasicTurbulenceModel::rhoField rhoField
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.