Go to the documentation of this file.
41 template<
class BasicTurbulenceModel>
44 return nuTilda_/this->
nu();
48 template<
class BasicTurbulenceModel>
55 return chi3/(chi3 +
pow3(Cv1_));
59 template<
class BasicTurbulenceModel>
66 return 1.0 - chi/(1.0 + chi*fv1);
70 template<
class BasicTurbulenceModel>
76 return Ct3_*
exp(-Ct4_*
sqr(chi));
80 template<
class BasicTurbulenceModel>
90 template<
class BasicTurbulenceModel>
104 + fv2(chi, fv1)*nuTilda_/
sqr(kappa_*dTilda),
111 template<
class BasicTurbulenceModel>
135 tr.ref().boundaryFieldRef() == 0.0;
141 template<
class BasicTurbulenceModel>
155 template<
class BasicTurbulenceModel>
179 if (lowReCorrection_)
192 (1 - Cb1_/(Cw1_*
sqr(kappa_)*fwStar_)*(ft2 + (1 - ft2)*fv2))
193 /
max(SMALL, (fv1*
max(scalar(1
e-10), 1 - ft2)))
202 template<
class BasicTurbulenceModel>
216 template<
class BasicTurbulenceModel>
222 this->nut_ = nuTilda_*fv1;
226 BasicTurbulenceModel::correctNut();
230 template<
class BasicTurbulenceModel>
233 correctNut(fv1(this->chi()));
239 template<
class BasicTurbulenceModel>
242 const alphaField&
alpha,
248 const word& propertiesName,
300 Cw1_(Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
407 if (
type == typeName)
416 template<
class BasicTurbulenceModel>
421 sigmaNut_.readIfPresent(this->coeffDict());
422 kappa_.readIfPresent(*
this);
424 Cb1_.readIfPresent(this->coeffDict());
425 Cb2_.readIfPresent(this->coeffDict());
426 Cw1_ = Cb1_/
sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
428 Cw3_.readIfPresent(this->coeffDict());
429 Cv1_.readIfPresent(this->coeffDict());
430 Cs_.readIfPresent(this->coeffDict());
432 CDES_.readIfPresent(this->coeffDict());
433 ck_.readIfPresent(this->coeffDict());
435 lowReCorrection_.readIfPresent(
"lowReCorrection", this->coeffDict());
436 Ct3_.readIfPresent(this->coeffDict());
437 Ct4_.readIfPresent(this->coeffDict());
438 fwStar_.readIfPresent(this->coeffDict());
447 template<
class BasicTurbulenceModel>
458 template<
class BasicTurbulenceModel>
471 this->mesh_.time().timeName(),
482 dTildaF = dTilda(chi, fv1,
fvc::grad(this->U_));
485 return sqr(this->
nut()/ck_/dTildaF);
489 template<
class BasicTurbulenceModel>
502 this->mesh_.time().timeName(),
515 template<
class BasicTurbulenceModel>
518 if (!this->turbulence_)
524 const alphaField&
alpha = this->alpha_;
525 const rhoField&
rho = this->rho_;
539 const volScalarField Stilda(this->Stilda(chi, fv1, Omega, dTilda));
548 Cb1_*
alpha*
rho*Stilda*nuTilda_*(scalar(1) - ft2)
551 (Cw1_*fw(Stilda, dTilda) - Cb1_/
sqr(kappa_)*ft2)
558 nuTildaEqn.ref().relax();
563 nuTilda_.correctBoundaryConditions();
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volScalarField &Omega, const volScalarField &dTilda) const
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
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.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
dimensionedTensor skew(const dimensionedTensor &dt)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
tmp< volScalarField > DnuTildaEff() const
static constexpr const zero Zero
static options & New(const fvMesh &mesh)
const dimensionedScalar alpha
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
static const wallDist & New(const fvMesh &mesh, Args &&... args)
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &Omega, const volScalarField &dTilda) const
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar exp(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
tmp< volScalarField > chi() const
virtual tmp< volScalarField > LESRegion() const
dimensionedScalar pow6(const dimensionedScalar &ds)
virtual tmp< volScalarField > k() const
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
bool readIfPresent(const dictionary &dict)
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
dimensionedScalar pow3(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual void printCoeffs(const word &type)
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
Generic dimensioned Type class.
const uniformDimensionedVectorField & g
BasicTurbulenceModel::transportModel transportModel
Base-class for all transport models used by the incompressible turbulence models.
void correctBoundaryConditions()
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
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)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
const dimensionedScalar e
tmp< volScalarField > ft2(const volScalarField &chi) const
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
BasicTurbulenceModel::alphaField alphaField
Generic GeometricField class.
const volScalarField & psi
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< volScalarField > fv1(const volScalarField &chi) const
dimensionedScalar neg(const dimensionedScalar &ds)
tmp< volScalarField > fw(const volScalarField &Omega, const volScalarField &dTilda) const
tmp< volScalarField > Omega(const volTensorField &gradU) const
virtual void correctNut()
const dimensionSet dimless
Templated abstract base class for DES models.
static Switch getOrAddToDict(const word &key, dictionary &dict, const Switch deflt=switchType::FALSE)