Go to the documentation of this file.
39 template<
class BasicTurbulenceModel>
42 this->nut_ = Cmu_*
sqr(k_)/epsilon_;
43 this->nut_.correctBoundaryConditions();
48 BasicTurbulenceModel::correctNut();
52 template<
class BasicTurbulenceModel>
60 dimVolume*this->rho_.dimensions()*k_.dimensions()
67 template<
class BasicTurbulenceModel>
75 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
84 template<
class BasicTurbulenceModel>
93 const word& propertiesName,
169 this->runTime_.timeName(),
181 this->runTime_.timeName(),
189 bound(k_, this->kMin_);
190 bound(epsilon_, this->epsilonMin_);
192 if (
type == typeName)
194 this->printCoeffs(
type);
201 template<
class BasicTurbulenceModel>
206 Cmu_.readIfPresent(this->coeffDict());
207 C1_.readIfPresent(this->coeffDict());
208 C2_.readIfPresent(this->coeffDict());
209 C3_.readIfPresent(this->coeffDict());
210 sigmak_.readIfPresent(this->coeffDict());
211 sigmaEps_.readIfPresent(this->coeffDict());
222 template<
class BasicTurbulenceModel>
225 if (!this->turbulence_)
247 epsilon_.boundaryField().updateCoeffs();
265 epsEqn().boundaryManipulate(epsilon_.boundaryField());
268 bound(epsilon_, this->epsilonMin_);
288 bound(k_, this->kMin_);
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static word groupName(Name name, const word &group)
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.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
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)
BasicTurbulenceModel::alphaField alphaField
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual bool read()
Re-read model coefficients if they have changed.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
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 > &)
virtual void correctNut()
Generic dimensioned Type class.
virtual tmp< fvScalarMatrix > epsilonSource() const
BasicTurbulenceModel::transportModel transportModel
kEpsilon(const kEpsilon &)
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)))
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Eddy viscosity turbulence model base class.
virtual tmp< fvScalarMatrix > kSource() const
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:
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)