Go to the documentation of this file.
38 template<
class BasicTurbulenceModel>
51 (2*
sqrt(2.0))*((S&S)&&S)
67 return 1.0/(A0_ + As*Us*k_/epsilon_);
71 template<
class BasicTurbulenceModel>
79 this->nut_ = rCmu(gradU, S2, magS)*
sqr(k_)/epsilon_;
80 this->nut_.correctBoundaryConditions();
82 BasicTurbulenceModel::correctNut();
86 template<
class BasicTurbulenceModel>
92 correctNut(tgradU(), S2, magS);
96 template<
class BasicTurbulenceModel>
104 dimVolume*this->rho_.dimensions()*k_.dimensions()
111 template<
class BasicTurbulenceModel>
119 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
128 template<
class BasicTurbulenceModel>
137 const word& propertiesName,
204 this->runTime_.timeName(),
216 this->runTime_.timeName(),
224 bound(k_, this->kMin_);
225 bound(epsilon_, this->epsilonMin_);
227 if (
type == typeName)
229 this->printCoeffs(
type);
236 template<
class BasicTurbulenceModel>
241 Cmu_.readIfPresent(this->coeffDict());
242 A0_.readIfPresent(this->coeffDict());
243 C2_.readIfPresent(this->coeffDict());
244 sigmak_.readIfPresent(this->coeffDict());
245 sigmaEps_.readIfPresent(this->coeffDict());
256 template<
class BasicTurbulenceModel>
259 if (!this->turbulence_)
285 epsilon_.boundaryField().updateCoeffs();
316 epsEqn().boundaryManipulate(epsilon_.boundaryField());
319 bound(epsilon_, this->epsilonMin_);
338 bound(k_, this->kMin_);
340 correctNut(tgradU(), S2, magS);
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
virtual tmp< fvScalarMatrix > epsilonSource() const
A class for handling words, derived from string.
dimensionedTensor skew(const dimensionedTensor &dt)
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.
virtual bool read()
Re-read model coefficients if they have changed.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Dimension set for the base types.
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
virtual tmp< fvScalarMatrix > kSource() const
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.
BasicTurbulenceModel::rhoField rhoField
Generic dimensioned Type class.
virtual void correctNut()
BasicTurbulenceModel::alphaField alphaField
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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)
dimensionedScalar acos(const dimensionedScalar &ds)
realizableKE(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.
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)
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
void clear() const
If object pointer points to valid object:
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
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)
dimensionedScalar cos(const dimensionedScalar &ds)