Go to the documentation of this file.
43 template<
class BasicTurbulenceModel>
46 return nuTilda_/this->
nu();
50 template<
class BasicTurbulenceModel>
58 return chi3/(chi3 +
pow3(Cv1_));
62 template<
class BasicTurbulenceModel>
69 return scalar(1) - chi/(scalar(1) + chi*fv1);
73 template<
class BasicTurbulenceModel>
90 Omega + fv2(chi(), fv1())*nuTilda_()/
sqr(kappa_*y_),
97 template<
class BasicTurbulenceModel>
131 template<
class BasicTurbulenceModel>
134 this->nut_ = nuTilda_*this->fv1(this->chi());
135 this->nut_.correctBoundaryConditions();
138 BasicTurbulenceModel::correctNut();
144 template<
class BasicTurbulenceModel>
153 const word& propertiesName,
206 Cw1_(Cb1_/
sqr(kappa_) + (scalar(1) + Cb2_)/sigmaNut_),
249 this->runTime_.timeName(),
259 if (
type == typeName)
261 this->printCoeffs(
type);
268 template<
class BasicTurbulenceModel>
273 sigmaNut_.readIfPresent(this->coeffDict());
274 kappa_.readIfPresent(this->coeffDict());
276 Cb1_.readIfPresent(this->coeffDict());
277 Cb2_.readIfPresent(this->coeffDict());
278 Cw1_ = Cb1_/
sqr(kappa_) + (scalar(1) + Cb2_)/sigmaNut_;
280 Cw3_.readIfPresent(this->coeffDict());
281 Cv1_.readIfPresent(this->coeffDict());
282 Cs_.readIfPresent(this->coeffDict());
291 template<
class BasicTurbulenceModel>
297 (nuTilda_ + this->
nu())/sigmaNut_
302 template<
class BasicTurbulenceModel>
306 const scalar Cmu = 0.09;
313 this->runTime_.timeName(),
316 cbrt(this->fv1(this->chi()))
318 *::
sqrt(scalar(2)/Cmu)
320 this->nut_.boundaryField().types()
325 template<
class BasicTurbulenceModel>
329 const scalar Cmu = 0.09;
337 this->runTime_.timeName(),
340 pow(this->fv1(this->chi()), 0.5)
342 /(nuTilda_ + this->nut_ + nutSMALL),
343 this->nut_.boundaryField().types()
348 template<
class BasicTurbulenceModel>
352 const scalar betaStar = 0.09;
360 this->runTime_.timeName(),
363 this->
epsilon()/(betaStar*(this->
k() + k0)),
364 this->nut_.boundaryField().types()
369 template<
class BasicTurbulenceModel>
372 if (!this->turbulence_)
395 Cb1_*
alpha()*
rho()*Stilda*nuTilda_()
400 nuTildaEqn.
ref().relax();
405 nuTilda_.correctBoundaryConditions();
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from Foam::string.
virtual void correctNut()
Update nut with the latest available nuTilda.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedTensor skew(const dimensionedTensor &dt)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
tmp< volScalarField > fv1(const volScalarField &chi) const
Templated abstract base class for RAS turbulence models.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
bool read(const char *buf, int32_t &val)
Same as readInt32.
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
virtual tmp< volScalarField > k() const
Return the (estimated) turbulent kinetic energy.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
virtual tmp< volScalarField > omega() const
Return the (estimated) specific dissipation rate.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
dimensionedScalar pow6(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
bool readIfPresent(const dictionary &dict)
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)
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Bound the given scalar field if it has gone unbounded.
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void correct()
Solve the turbulence equations and correct the turbulent viscosity.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
BasicTurbulenceModel::rhoField rhoField
Generic dimensioned Type class.
const uniformDimensionedVectorField & g
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::transportModel transportModel
Spalart-Allmaras one-transport-equation linear-eddy-viscosity turbulence closure model for incompress...
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< volScalarField::Internal > Stilda() const
tmp< volScalarField::Internal > fw(const volScalarField::Internal &Stilda) const
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< volScalarField > chi() const
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
label k
Boltzmann constant.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Eddy viscosity turbulence model base class.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dimensionedScalar cbrt(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual tmp< volScalarField > epsilon() const
Return the (estimated) turbulent kinetic energy dissipation rate.
tmp< volScalarField::Internal > fv2(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const