Go to the documentation of this file.
38 template<
class BasicTurbulenceModel>
60 template<
class BasicTurbulenceModel>
69 pow(k_, 1.5)/epsilon_,
85 template<
class BasicTurbulenceModel>
88 this->nut_ =
min(CmuKEps_*
sqr(k_)/epsilon_, this->Cmu_*v2_*Ts());
89 this->nut_.correctBoundaryConditions();
91 BasicTurbulenceModel::correctNut();
97 template<
class BasicTurbulenceModel>
106 const word& propertiesName,
219 this->runTime_.timeName(),
231 this->runTime_.timeName(),
243 this->runTime_.timeName(),
255 this->runTime_.timeName(),
265 bound(k_, this->kMin_);
266 bound(epsilon_, this->epsilonMin_);
270 if (
type == typeName)
272 this->printCoeffs(
type);
279 template<
class BasicTurbulenceModel>
284 Cmu_.readIfPresent(this->coeffDict());
285 CmuKEps_.readIfPresent(this->coeffDict());
286 C1_.readIfPresent(this->coeffDict());
287 C2_.readIfPresent(this->coeffDict());
288 CL_.readIfPresent(this->coeffDict());
289 Ceta_.readIfPresent(this->coeffDict());
290 Ceps2_.readIfPresent(this->coeffDict());
291 Ceps3_.readIfPresent(this->coeffDict());
292 sigmaK_.readIfPresent(this->coeffDict());
293 sigmaEps_.readIfPresent(this->coeffDict());
304 template<
class BasicTurbulenceModel>
307 if (!this->turbulence_)
335 1.0/Ts*((C1_ -
N)*v2_ - 2.0/3.0*k_*(C1_ - 1.0))
341 1.4*(1.0 + 0.05*
min(
sqrt(k_/v2_), scalar(100.0)))
345 epsilon_.boundaryField().updateCoeffs();
361 epsEqn().boundaryManipulate(epsilon_.boundaryField());
364 bound(epsilon_, this->epsilonMin_);
381 bound(k_, this->kMin_);
390 - 1.0/L2/k_*(v2fAlpha - C2_*
G)
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)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
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.
Templated abstract base class for RAS turbulence models.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool read(const char *, int32_t &)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
virtual bool read()
Re-read model coefficients if they have changed.
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< volScalarField > Ts() const
Return time scale, Ts.
dimensionedScalar pow025(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Bound the given scalar field if it has gone unbounded.
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
v2f(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
BasicTurbulenceModel::rhoField rhoField
Generic dimensioned Type class.
BasicTurbulenceModel::alphaField alphaField
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
BasicTurbulenceModel::transportModel transportModel
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void correctNut()
Eddy viscosity turbulence model base class.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
tmp< volScalarField > Ls() const
Return length scale, Ls.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)