Go to the documentation of this file.
38 template<
class BasicTurbulenceModel>
41 this->nut_ = Cmu_*
sqr(k_)/epsilon_;
42 this->nut_.correctBoundaryConditions();
44 BasicTurbulenceModel::correctNut();
48 template<
class BasicTurbulenceModel>
56 dimVolume*this->rho_.dimensions()*k_.dimensions()
63 template<
class BasicTurbulenceModel>
71 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
80 template<
class BasicTurbulenceModel>
89 const word& propertiesName,
183 this->runTime_.timeName(),
195 this->runTime_.timeName(),
203 bound(k_, this->kMin_);
204 bound(epsilon_, this->epsilonMin_);
206 if (
type == typeName)
208 this->printCoeffs(
type);
215 template<
class BasicTurbulenceModel>
220 Cmu_.readIfPresent(this->coeffDict());
221 C1_.readIfPresent(this->coeffDict());
222 C2_.readIfPresent(this->coeffDict());
223 C3_.readIfPresent(this->coeffDict());
224 sigmak_.readIfPresent(this->coeffDict());
225 sigmaEps_.readIfPresent(this->coeffDict());
226 eta0_.readIfPresent(this->coeffDict());
227 beta_.readIfPresent(this->coeffDict());
238 template<
class BasicTurbulenceModel>
241 if (!this->turbulence_)
268 ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
272 epsilon_.boundaryField().updateCoeffs();
289 epsEqn().boundaryManipulate(epsilon_.boundaryField());
292 bound(epsilon_, this->epsilonMin_);
311 bound(k_, this->kMin_);
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.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool read(const char *, int32_t &)
virtual tmp< fvScalarMatrix > kSource() const
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 tmp< fvScalarMatrix > epsilonSource() const
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
#define R(A, B, C, D, E, F, K, M)
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
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 > &)
BasicTurbulenceModel::rhoField rhoField
Generic dimensioned Type class.
BasicTurbulenceModel::alphaField alphaField
virtual void correctNut()
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Re-read model coefficients if they have changed.
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)
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Eddy viscosity turbulence model base class.
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.
RNGkEpsilon(const RNGkEpsilon &)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)