Go to the documentation of this file.
37 template<
class BasicTurbulenceModel>
40 return nuTilda_/this->
nu();
44 template<
class BasicTurbulenceModel>
51 return chi3/(chi3 +
pow3(Cv1_));
55 template<
class BasicTurbulenceModel>
62 return 1.0 - chi/(1.0 + chi*fv1);
66 template<
class BasicTurbulenceModel>
72 return Ct3_*
exp(-Ct4_*
sqr(chi));
76 template<
class BasicTurbulenceModel>
86 template<
class BasicTurbulenceModel>
100 + fv2(chi, fv1)*nuTilda_/
sqr(kappa_*dTilda),
107 template<
class BasicTurbulenceModel>
131 tr().boundaryField() == 0.0;
137 template<
class BasicTurbulenceModel>
151 template<
class BasicTurbulenceModel>
175 if (lowReCorrection_)
188 (1 - Cb1_/(Cw1_*
sqr(kappa_)*fwStar_)*(ft2 + (1 - ft2)*fv2))
189 /
max(SMALL, (fv1*
max(scalar(1
e-10), 1 - ft2)))
198 template<
class BasicTurbulenceModel>
212 template<
class BasicTurbulenceModel>
218 this->nut_ = nuTilda_*fv1;
221 BasicTurbulenceModel::correctNut();
225 template<
class BasicTurbulenceModel>
228 correctNut(fv1(this->chi()));
234 template<
class BasicTurbulenceModel>
243 const word& propertiesName,
295 Cw1_(Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
392 this->runTime_.timeName(),
402 if (
type == typeName)
404 this->printCoeffs(
type);
411 template<
class BasicTurbulenceModel>
416 sigmaNut_.readIfPresent(this->coeffDict());
417 kappa_.readIfPresent(*
this);
419 Cb1_.readIfPresent(this->coeffDict());
420 Cb2_.readIfPresent(this->coeffDict());
421 Cw1_ = Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
423 Cw3_.readIfPresent(this->coeffDict());
424 Cv1_.readIfPresent(this->coeffDict());
425 Cs_.readIfPresent(this->coeffDict());
427 CDES_.readIfPresent(this->coeffDict());
428 ck_.readIfPresent(this->coeffDict());
430 lowReCorrection_.readIfPresent(
"lowReCorrection", this->coeffDict());
431 Ct3_.readIfPresent(this->coeffDict());
432 Ct4_.readIfPresent(this->coeffDict());
433 fwStar_.readIfPresent(this->coeffDict());
444 template<
class BasicTurbulenceModel>
455 template<
class BasicTurbulenceModel>
464 template<
class BasicTurbulenceModel>
477 this->mesh_.time().timeName(),
490 template<
class BasicTurbulenceModel>
493 if (!this->turbulence_)
513 const volScalarField Stilda(this->Stilda(chi, fv1, Omega, dTilda));
522 Cb1_*
alpha*
rho*Stilda*nuTilda_*(scalar(1) - ft2)
525 (Cw1_*fw(Stilda, dTilda) - Cb1_/
sqr(kappa_)*ft2)
531 nuTildaEqn().relax();
534 nuTilda_.correctBoundaryConditions();
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volScalarField &Omega, const volScalarField &dTilda) const
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from string.
dimensionedTensor skew(const dimensionedTensor &dt)
A class for managing temporary objects.
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
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.
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
rDeltaT dimensionedInternalField()
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &Omega, const volScalarField &dTilda) const
BasicTurbulenceModel::rhoField rhoField
const dimensionedVector & g
dimensionedScalar exp(const dimensionedScalar &ds)
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 > chi() const
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
dimensionedScalar pow6(const dimensionedScalar &ds)
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
dimensionedScalar pow3(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual bool read()
Re-read model coefficients if they have changed.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Generic dimensioned Type class.
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
BasicTurbulenceModel::transportModel transportModel
const double e
Elementary charge.
void correctBoundaryConditions()
Correct boundary field.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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)
tmp< volScalarField > psi(const volScalarField &chi, const volScalarField &fv1) const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const volScalarField & psi
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< volScalarField > ft2(const volScalarField &chi) const
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
BasicTurbulenceModel::alphaField alphaField
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)
SpalartAllmarasDES(const SpalartAllmarasDES &)
virtual void correct()
Correct nuTilda and related properties.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< volScalarField > fv1(const volScalarField &chi) const
dimensionedScalar neg(const dimensionedScalar &ds)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
tmp< volScalarField > fw(const volScalarField &Omega, const volScalarField &dTilda) const
tmp< volScalarField > Omega(const volTensorField &gradU) const
virtual void correctNut()
static Switch lookupOrAddToDict(const word &, dictionary &, const Switch &defaultValue=false)
Construct from dictionary, supplying default value so that if the.