Go to the documentation of this file.
37 template<
class BasicTurbulenceModel>
44 BasicTurbulenceModel::correctNut();
48 template<
class BasicTurbulenceModel>
55 template<
class BasicTurbulenceModel>
65 return min(CDES*this->
delta(),
sqrt(k)/(this->betaStar_*omega));
71 template<
class BasicTurbulenceModel>
80 const word& propertiesName,
126 if (
type == typeName)
128 this->printCoeffs(
type);
135 template<
class BasicTurbulenceModel>
140 kappa_.readIfPresent(this->coeffDict());
141 CDESkom_.readIfPresent(this->coeffDict());
142 CDESkeps_.readIfPresent(this->coeffDict());
153 template<
class BasicTurbulenceModel>
156 if (!this->turbulence_)
206 + this->omegaSource()
211 omegaEqn().boundaryManipulate(omega.boundaryField());
214 bound(omega, this->omegaMin_);
239 this->correctNut(S2);
243 template<
class BasicTurbulenceModel>
264 this->mesh_.time().timeName(),
274 F1*CDESkom_ + (1 -
F1)*CDESkeps_
276 -
sqrt(
k)/(this->betaStar_*omega)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
void updateCoeffs()
Update the boundary condition coefficients.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from string.
const dimensionedScalar G
Newtonian constant of gravitation.
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool read(const char *, int32_t &)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
dimensioned< scalar > mag(const dimensioned< Type > &)
kOmegaSSTDES(const kOmegaSSTDES &)
virtual bool read()
Re-read model coefficients if they have changed.
One equation eddy-viscosity model.
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
Generic dimensioned Type class.
BasicTurbulenceModel::transportModel transportModel
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
BasicTurbulenceModel::alphaField alphaField
virtual void correctNut()
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
void clear() const
If object pointer points to valid object:
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flo...
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)