42template<
class BasicEddyViscosityModel>
72template<
class BasicEddyViscosityModel>
79 (scalar(2)/betaStar_)*
sqrt(k_)/(omega_*y_),
80 scalar(500)*(this->
mu()/this->rho_)/(
sqr(y_)*omega_)
89template<
class BasicEddyViscosityModel>
94 150*(this->
mu()/this->rho_)/(omega_*
sqr(y_)),
102template<
class BasicEddyViscosityModel>
116template<
class BasicEddyViscosityModel>
124 this->nut_.correctBoundaryConditions();
129template<
class BasicEddyViscosityModel>
136template<
class BasicEddyViscosityModel>
146template<
class BasicEddyViscosityModel>
156template<
class BasicEddyViscosityModel>
167template<
class BasicEddyViscosityModel>
182template<
class BasicEddyViscosityModel>
198template<
class BasicEddyViscosityModel>
209template<
class BasicEddyViscosityModel>
220template<
class BasicEddyViscosityModel>
238template<
class BasicEddyViscosityModel>
242 const alphaField&
alpha,
248 const word& propertiesName
251 BasicEddyViscosityModel
387 IOobject::groupName(
"k", alphaRhoPhi.group()),
399 IOobject::groupName(
"omega", alphaRhoPhi.group()),
447template<
class BasicEddyViscosityModel>
453 decayControl_.readIfPresent(
"decayControl",
dict);
458 omegaInf_.read(
dict);
460 Info<<
" Employing decay control with kInf:" << kInf_
461 <<
" and omegaInf:" << omegaInf_ <<
endl;
471template<
class BasicEddyViscosityModel>
474 if (BasicEddyViscosityModel::read())
476 alphaK1_.readIfPresent(this->coeffDict());
477 alphaK2_.readIfPresent(this->coeffDict());
478 alphaOmega1_.readIfPresent(this->coeffDict());
479 alphaOmega2_.readIfPresent(this->coeffDict());
480 gamma1_.readIfPresent(this->coeffDict());
481 gamma2_.readIfPresent(this->coeffDict());
482 beta1_.readIfPresent(this->coeffDict());
483 beta2_.readIfPresent(this->coeffDict());
484 betaStar_.readIfPresent(this->coeffDict());
485 a1_.readIfPresent(this->coeffDict());
486 b1_.readIfPresent(this->coeffDict());
487 c1_.readIfPresent(this->coeffDict());
488 F3_.readIfPresent(
"F3", this->coeffDict());
490 setDecayControl(this->coeffDict());
499template<
class BasicEddyViscosityModel>
502 if (!this->turbulence_)
508 const alphaField&
alpha = this->alpha_;
509 const rhoField&
rho = this->rho_;
515 BasicEddyViscosityModel::correct();
541 omega_.boundaryFieldRef().updateCoeffs();
546 omega_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
567 GbyNu0 = GbyNu(GbyNu0, F23(), S2());
581 alpha()*
rho()*(
F1() - scalar(1))*CDkOmega()/omega_(),
590 omegaEqn.ref().relax();
592 omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
595 bound(omega_, this->omegaMin_);
609 +
alpha()*
rho()*betaStar_*omegaInf_*kInf_
620 bound(k_, this->kMin_);
Bound the given scalar field if it has gone unbounded.
DimensionedField< scalar, volMesh > Internal
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Generic dimensioned Type class.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
void setDecayControl(const dictionary &dict)
volScalarField k_
Turbulent kinetic energy field [m^2/s^2].
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
Switch F3_
Flag to include the F3 term.
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
kOmegaSSTBase(const kOmegaSSTBase &)=delete
No copy construct.
virtual void correctNut(const volScalarField &S2)
dimensionedScalar alphaOmega2_
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volScalarField > F2() const
dimensionedScalar gamma2_
virtual tmp< volScalarField > F23() const
BasicEddyViscosityModel::transportModel transportModel
dimensionedScalar omegaInf_
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
dimensionedScalar alphaK2_
const volScalarField & y_
Wall distance.
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
virtual tmp< volScalarField::Internal > GbyNu0(const volTensorField &gradU, const volScalarField &S2) const
Return (G/nu)_0.
dimensionedScalar alphaK1_
volScalarField omega_
Specific dissipation rate field [1/s].
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
virtual tmp< fvScalarMatrix > omegaSource() const
virtual void correctNut()
virtual tmp< fvScalarMatrix > kSource() const
dimensionedScalar alphaOmega1_
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Return epsilon/k which for standard RAS is betaStar*omega.
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
dimensionedScalar gamma1_
virtual tmp< volScalarField > S2(const volTensorField &gradU) const
Return square of strain rate.
BasicEddyViscosityModel::rhoField rhoField
BasicEddyViscosityModel::alphaField alphaField
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.
dimensionedScalar betaStar_
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > F3() const
A class for managing temporary objects.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Base-class for all transport models used by the incompressible turbulence models.
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
A class for handling words, derived from Foam::string.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Ostream & endl(Ostream &os)
Add newline and flush stream.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
dimensionedScalar pow4(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)