Go to the documentation of this file.
42 template<
class BasicTurbulenceModel>
55 (2*
sqrt(2.0))*((S&S)&&S)
71 return 1.0/(A0_ + As*
Us*k_/epsilon_);
75 template<
class BasicTurbulenceModel>
83 this->nut_ = rCmu(gradU, S2, magS)*
sqr(k_)/epsilon_;
84 this->nut_.correctBoundaryConditions();
87 BasicTurbulenceModel::correctNut();
91 template<
class BasicTurbulenceModel>
97 correctNut(tgradU(), S2, magS);
101 template<
class BasicTurbulenceModel>
109 dimVolume*this->rho_.dimensions()*k_.dimensions()
116 template<
class BasicTurbulenceModel>
124 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
133 template<
class BasicTurbulenceModel>
136 const alphaField&
alpha,
142 const word& propertiesName,
222 if (type == typeName)
231 template<
class BasicTurbulenceModel>
236 A0_.readIfPresent(this->coeffDict());
237 C2_.readIfPresent(this->coeffDict());
238 sigmak_.readIfPresent(this->coeffDict());
239 sigmaEps_.readIfPresent(this->coeffDict());
248 template<
class BasicTurbulenceModel>
251 if (!this->turbulence_)
257 const alphaField&
alpha = this->alpha_;
258 const rhoField&
rho = this->rho_;
264 eddyViscosity<RASModel<BasicTurbulenceModel>>
::correct();
278 epsilon_.boundaryFieldRef().updateCoeffs();
292 tmp<fvScalarMatrix> epsEqn
308 epsEqn.ref().relax();
310 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
313 bound(epsilon_, this->epsilonMin_);
318 tmp<fvScalarMatrix> kEqn
335 bound(k_, this->kMin_);
337 correctNut(tgradU(), S2, magS);
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)
virtual tmp< fvScalarMatrix > epsilonSource() const
A class for handling words, derived from Foam::string.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
void clear() const noexcept
constexpr const char *const group
dimensionedTensor skew(const dimensionedTensor &dt)
const dimensionedScalar G
A class for managing temporary objects.
static constexpr const zero Zero
Templated abstract base class for RAS turbulence models.
static options & New(const fvMesh &mesh)
const dimensionedScalar alpha
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
bool read(const char *buf, int32_t &val)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
dimensionedScalar epsilonMin_
label min(const labelHashSet &set, label minValue=labelMax)
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
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 > &)
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.
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label max(const labelHashSet &set, label maxValue=labelMin)
BasicTurbulenceModel::rhoField rhoField
Generic dimensioned Type class.
virtual void correctNut()
BasicTurbulenceModel::alphaField alphaField
Base-class for all transport models used by the incompressible turbulence models.
GeometricField< vector, fvPatchField, volMesh > volVectorField
BasicTurbulenceModel::transportModel transportModel
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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)
fileName::Type type(const fileName &name, const bool followLink=true)
A special matrix type and solver, designed for finite volume solutions of scalar equations....
virtual void printCoeffs(const word &type)
Eddy viscosity turbulence model base class.
static word groupName(StringType base, const word &group)
const dimensionSet dimVolume(pow3(dimLength))
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
edgeScalarField phis(IOobject("phis", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), linearEdgeInterpolate(Us) &aMesh.Le())
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar cos(const dimensionedScalar &ds)