Go to the documentation of this file.
37 template<
class BasicTurbulenceModel>
46 simpleFilter_(
dev(filter_(
sqr(this->U_)) - (
sqr(filter_(this->U_)))))
62 simpleFilter_(0.5*(LL && MM))
74 template<
class BasicTurbulenceModel>
83 simpleFilter_(this->nuEff()*(filter_(
magSqr(D)) -
magSqr(filter_(D))))
84 /simpleFilter_(
pow(KK, 1.5)/(2.0*this->
delta()))
92 template<
class BasicTurbulenceModel>
99 0.5*(filter_(
magSqr(this->U_)) -
magSqr(filter_(this->U_)))
107 template<
class BasicTurbulenceModel>
114 this->nut_ = Ck(D, KK)*
sqrt(k_)*this->
delta();
115 this->nut_.correctBoundaryConditions();
119 template<
class BasicTurbulenceModel>
124 0.5*(filter_(
magSqr(this->U_)) -
magSqr(filter_(this->U_)))
131 template<
class BasicTurbulenceModel>
139 dimVolume*this->rho_.dimensions()*k_.dimensions()
148 template<
class BasicTurbulenceModel>
157 const word& propertiesName,
178 this->runTime_.timeName(),
186 simpleFilter_(this->mesh_),
188 filter_(filterPtr_())
190 bound(k_, this->kMin_);
192 if (
type == typeName)
194 this->printCoeffs(
type);
201 template<
class BasicTurbulenceModel>
206 filter_.read(this->coeffDict());
217 template<
class BasicTurbulenceModel>
227 this->runTime_.timeName(),
238 template<
class BasicTurbulenceModel>
241 if (!this->turbulence_)
279 bound(k_, this->kMin_);
BasicTurbulenceModel::alphaField alphaField
dynamicKEqn(const dynamicKEqn &)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static word groupName(Name name, const word &group)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from string.
const dimensionedScalar G
Newtonian constant of gravitation.
A class for managing temporary objects.
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.
BasicTurbulenceModel::transportModel transportModel
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void correctNut()
One equation eddy-viscosity model.
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual bool read()
Read model coefficients if they have changed.
BasicTurbulenceModel::rhoField rhoField
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Eddy viscosity LES SGS model base class.
virtual void correct()
Correct Eddy-Viscosity and related properties.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void max(const dimensioned< Type > &)
volScalarField Ce() const
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static autoPtr< LESfilter > New(const fvMesh &, const dictionary &, const word &filterDictName="filter")
Return a reference to the selected LES filter.
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
dimensionedScalar sqrt(const dimensionedScalar &ds)
volScalarField Ck(const volSymmTensorField &D, const volScalarField &KK) const
Calculate Ck by filtering the velocity field U.
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
label k
Boltzmann constant.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
const dimensionSet dimVolume(pow3(dimLength))
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
void clear() const
If object pointer points to valid object:
virtual tmp< fvScalarMatrix > kSource() const
dimensioned< scalar > magSqr(const dimensioned< Type > &)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)