Go to the documentation of this file.
39 template<
class BasicTurbulenceModel>
42 return nuTilda_/this->
nu();
46 template<
class BasicTurbulenceModel>
53 return chi3/(chi3 +
pow3(Cv1_));
57 template<
class BasicTurbulenceModel>
64 return 1.0 - chi/(1.0 + chi*fv1);
68 template<
class BasicTurbulenceModel>
82 + fv2(chi, fv1)*nuTilda_/
sqr(kappa_*y_),
89 template<
class BasicTurbulenceModel>
119 template<
class BasicTurbulenceModel>
125 this->nut_ = nuTilda_*fv1;
128 BasicTurbulenceModel::correctNut();
132 template<
class BasicTurbulenceModel>
135 correctNut(fv1(this->chi()));
141 template<
class BasicTurbulenceModel>
150 const word& propertiesName,
203 Cw1_(Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
246 this->runTime_.timeName(),
256 if (
type == typeName)
258 this->printCoeffs(
type);
265 template<
class BasicTurbulenceModel>
270 sigmaNut_.readIfPresent(this->coeffDict());
271 kappa_.readIfPresent(this->coeffDict());
273 Cb1_.readIfPresent(this->coeffDict());
274 Cb2_.readIfPresent(this->coeffDict());
275 Cw1_ = Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
277 Cw3_.readIfPresent(this->coeffDict());
278 Cv1_.readIfPresent(this->coeffDict());
279 Cs_.readIfPresent(this->coeffDict());
290 template<
class BasicTurbulenceModel>
300 template<
class BasicTurbulenceModel>
310 this->runTime_.timeName(),
320 template<
class BasicTurbulenceModel>
324 <<
"Turbulence kinetic energy dissipation rate not defined for "
325 <<
"Spalart-Allmaras model. Returning zero field"
335 this->runTime_.timeName(),
345 template<
class BasicTurbulenceModel>
348 if (!this->turbulence_)
376 nuTildaEqn().relax();
379 nuTilda_.correctBoundaryConditions();
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from string.
virtual void correctNut()
dimensionedTensor skew(const dimensionedTensor &dt)
A class for managing temporary objects.
tmp< volScalarField > fv1(const volScalarField &chi) const
Templated abstract base class for RAS turbulence models.
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 tmp< volScalarField > k() const
Return the turbulence kinetic energy.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
const dimensionedVector & g
Ostream & endl(Ostream &os)
Add newline and flush stream.
Dimension set for the base types.
dimensioned< scalar > mag(const dimensioned< Type > &)
static const wallDist & New(const fvMesh &mesh)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1) const
SpalartAllmaras(const SpalartAllmaras &)
tmp< volScalarField > fw(const volScalarField &Stilda) const
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual bool read()
Re-read model coefficients if they have changed.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Bound the given scalar field if it has gone unbounded.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
BasicTurbulenceModel::rhoField rhoField
Generic dimensioned Type class.
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
BasicTurbulenceModel::alphaField alphaField
void correctBoundaryConditions()
Correct boundary field.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< volScalarField > chi() const
Eddy viscosity turbulence model base class.
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)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define WarningInFunction
Report a warning using Foam::Warning.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
dimensioned< scalar > magSqr(const dimensioned< Type > &)