41template<
class BasicTurbulenceModel>
45 this->
nut_.correctBoundaryConditions();
48 BasicTurbulenceModel::correctNut();
52template<
class BasicTurbulenceModel>
65template<
class BasicTurbulenceModel>
68 const alphaField&
alpha,
74 const word& propertiesName,
94 IOobject::groupName(
"k", this->alphaRhoPhi_.group()),
117 this->printCoeffs(
type);
124template<
class BasicTurbulenceModel>
129 Ck_.readIfPresent(this->coeffDict());
138template<
class BasicTurbulenceModel>
141 if (!this->turbulence_)
147 const alphaField&
alpha = this->alpha_;
148 const rhoField&
rho = this->rho_;
162 tmp<fvScalarMatrix> kEqn
179 bound(k_, this->kMin_);
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Eddy viscosity LES SGS model base class.
virtual bool read()
Read model coefficients if they have changed.
BasicTurbulenceModel::alphaField alphaField
kEqn(const kEqn &)=delete
No copy construct.
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Correct eddy-Viscosity and related properties.
virtual void correctNut()
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Generic dimensioned Type class.
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.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
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.
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Namespace for LES SGS models.
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
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)
const dimensionSet dimVolume(pow3(dimLength))