42template<
class BasicTurbulenceModel>
49template<
class BasicTurbulenceModel>
58template<
class BasicTurbulenceModel>
61 this->nut_ = Cmu_*fMu()*
sqr(k_)/epsilon_;
62 this->nut_.correctBoundaryConditions();
65 BasicTurbulenceModel::correctNut();
69template<
class BasicTurbulenceModel>
80template<
class BasicTurbulenceModel>
93template<
class BasicTurbulenceModel>
94LaunderSharmaKE<BasicTurbulenceModel>::LaunderSharmaKE
96 const alphaField&
alpha,
102 const word& propertiesName,
204 this->printCoeffs(
type);
211template<
class BasicTurbulenceModel>
216 Cmu_.readIfPresent(this->coeffDict());
217 C1_.readIfPresent(this->coeffDict());
218 C2_.readIfPresent(this->coeffDict());
219 C3_.readIfPresent(this->coeffDict());
220 sigmak_.readIfPresent(this->coeffDict());
221 sigmaEps_.readIfPresent(this->coeffDict());
230template<
class BasicTurbulenceModel>
233 if (!this->turbulence_)
239 const alphaField&
alpha = this->alpha_;
240 const rhoField&
rho = this->rho_;
246 eddyViscosity<RASModel<BasicTurbulenceModel>>
::correct();
262 tmp<fvScalarMatrix> epsEqn
276 epsEqn.ref().relax();
278 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
281 bound(epsilon_, this->epsilonMin_);
285 tmp<fvScalarMatrix> kEqn
301 bound(k_, this->kMin_);
Bound the given scalar field if it has gone unbounded.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Templated abstract base class for RAS turbulence models.
BasicTurbulenceModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< volScalarField > f2() const
dimensionedScalar sigmak_
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar sigmaEps_
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
virtual void correctNut()
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::transportModel transportModel
tmp< volScalarField > fMu() const
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
virtual tmp< volScalarField > nut() const
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
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.
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.
A class for handling words, derived from Foam::string.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< volScalarField > magSqrGradGrad(const GeometricField< Type, fvPatchField, volMesh > &vf)
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.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read(const char *buf, int32_t &val)
Same as readInt32.
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.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
const dimensionedScalar & D