Go to the documentation of this file.
41 template<
class BasicTurbulenceModel>
50 simpleFilter_(
dev(filter_(
sqr(this->U_)) - (
sqr(filter_(this->U_)))))
66 simpleFilter_(0.5*(LL && MM))
78 template<
class BasicTurbulenceModel>
87 simpleFilter_(this->nuEff()*(filter_(
magSqr(
D)) -
magSqr(filter_(
D))))
88 /simpleFilter_(
pow(KK, 1.5)/(2.0*this->
delta()))
96 template<
class BasicTurbulenceModel>
103 0.5*(filter_(
magSqr(this->U_)) -
magSqr(filter_(this->U_)))
111 template<
class BasicTurbulenceModel>
118 this->nut_ = Ck(
D, KK)*
sqrt(k_)*this->
delta();
119 this->nut_.correctBoundaryConditions();
122 BasicTurbulenceModel::correctNut();
126 template<
class BasicTurbulenceModel>
131 0.5*(filter_(
magSqr(this->U_)) -
magSqr(filter_(this->U_)))
138 template<
class BasicTurbulenceModel>
146 dimVolume*this->rho_.dimensions()*k_.dimensions()
155 template<
class BasicTurbulenceModel>
164 const word& propertiesName,
185 this->runTime_.timeName(),
193 simpleFilter_(this->mesh_),
195 filter_(filterPtr_())
197 bound(k_, this->kMin_);
199 if (
type == typeName)
201 this->printCoeffs(
type);
208 template<
class BasicTurbulenceModel>
213 filter_.read(this->coeffDict());
222 template<
class BasicTurbulenceModel>
232 this->runTime_.timeName(),
243 template<
class BasicTurbulenceModel>
253 this->runTime_.timeName(),
261 template<
class BasicTurbulenceModel>
264 if (!this->turbulence_)
306 bound(k_, this->kMin_);
BasicTurbulenceModel::alphaField alphaField
Defines the attributes of an object for which implicit objectRegistry management is supported,...
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 Foam::string.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
void clear() const noexcept
const dimensionedScalar G
Newtonian constant of gravitation.
static constexpr const zero Zero
Global zero (0)
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
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)
virtual void correctNut()
One equation eddy-viscosity model.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
virtual bool read()
Read model coefficients if they have changed.
Dynamic one equation eddy-viscosity model.
BasicTurbulenceModel::rhoField rhoField
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
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)
Find the max value in labelHashSet, optionally limited by second argument.
Eddy viscosity LES SGS model base class.
virtual void correct()
Correct Eddy-Viscosity and related properties.
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.
virtual tmp< volScalarField > omega() const
Return sub-grid specific dissipation rate.
dimensionedScalar sqrt(const dimensionedScalar &ds)
volScalarField Ck(const volSymmTensorField &D, const volScalarField &KK) const
Calculate Ck by filtering the velocity field U.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
virtual tmp< volScalarField > epsilon() const
Return sub-grid dissipation rate.
label k
Boltzmann constant.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
const dimensionedScalar & D
A special matrix type and solver, designed for finite volume solutions of scalar equations....
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const dimensionSet dimVolume(pow3(dimLength))
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual tmp< fvScalarMatrix > kSource() const
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)