44 #ifndef SpalartAllmarasDES_H 45 #define SpalartAllmarasDES_H 60 template<
class BasicTurbulenceModel>
145 typedef typename BasicTurbulenceModel::alphaField
alphaField;
146 typedef typename BasicTurbulenceModel::rhoField
rhoField;
147 typedef typename BasicTurbulenceModel::transportModel
transportModel;
159 const alphaField&
alpha,
164 const transportModel& transport,
Eddy viscosity LES SGS model base class.
TypeName("SpalartAllmarasDES")
Runtime type information.
SpalartAllmarasDES(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.
void operator=(const SpalartAllmarasDES &)=delete
Disallow default bitwise assignment.
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volScalarField &Omega, const volScalarField &dTilda) const
const volScalarField & y_
Wall distance.
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &Omega, const volScalarField &dTilda) const
virtual bool read()
Read model coefficients if they have changed.
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
virtual ~SpalartAllmarasDES()
Destructor.
tmp< volScalarField > fv1(const volScalarField &chi) const
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
tmp< volScalarField > S(const volTensorField &gradU) const
BasicTurbulenceModel::rhoField rhoField
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
tmp< volScalarField > nuTilda() const
BasicTurbulenceModel::alphaField alphaField
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.
dimensionedScalar sigmaNut_
tmp< volScalarField > chi() const
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
virtual void correct()
Correct nuTilda and related properties.
virtual void correctNut()
tmp< volScalarField > fw(const volScalarField &Omega, const volScalarField &dTilda) const
tmp< volScalarField > Omega(const volTensorField &gradU) const
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.