Go to the documentation of this file.
42 template<
class BasicTurbulenceModel>
45 return exp(-3.4/
sqr(scalar(1) +
sqr(k_)/(this->
nu()*epsilon_)/50.0));
49 template<
class BasicTurbulenceModel>
58 template<
class BasicTurbulenceModel>
61 this->nut_ = Cmu_*fMu()*
sqr(k_)/epsilon_;
62 this->nut_.correctBoundaryConditions();
65 BasicTurbulenceModel::correctNut();
69 template<
class BasicTurbulenceModel>
77 dimVolume*this->rho_.dimensions()*k_.dimensions()
84 template<
class BasicTurbulenceModel>
92 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
101 template<
class BasicTurbulenceModel>
104 const alphaField&
alpha,
110 const word& propertiesName,
210 if (type == typeName)
219 template<
class BasicTurbulenceModel>
224 Cmu_.readIfPresent(this->coeffDict());
225 C1_.readIfPresent(this->coeffDict());
226 C2_.readIfPresent(this->coeffDict());
227 C3_.readIfPresent(this->coeffDict());
228 sigmak_.readIfPresent(this->coeffDict());
229 sigmaEps_.readIfPresent(this->coeffDict());
238 template<
class BasicTurbulenceModel>
241 if (!this->turbulence_)
247 const alphaField&
alpha = this->alpha_;
248 const rhoField&
rho = this->rho_;
254 eddyViscosity<RASModel<BasicTurbulenceModel>>
::correct();
270 tmp<fvScalarMatrix> epsEqn
284 epsEqn.ref().relax();
286 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
289 bound(epsilon_, this->epsilonMin_);
293 tmp<fvScalarMatrix> kEqn
309 bound(k_, this->kMin_);
Defines the attributes of an object for which implicit objectRegistry management is supported,...
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from Foam::string.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
void clear() const noexcept
const dimensionedScalar G
A class for managing temporary objects.
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)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar epsilonMin_
tmp< volScalarField > fMu() const
label min(const labelHashSet &set, label minValue=labelMax)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
BasicTurbulenceModel::alphaField alphaField
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 &)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
tmp< volScalarField > f2() const
BasicTurbulenceModel::rhoField rhoField
Generic dimensioned Type class.
Base-class for all transport models used by the incompressible turbulence models.
GeometricField< vector, fvPatchField, volMesh > volVectorField
virtual tmp< fvScalarMatrix > kSource() const
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)
fileName::Type type(const fileName &name, const bool followLink=true)
tmp< volScalarField > magSqrGradGrad(const GeometricField< Type, fvPatchField, volMesh > &vf)
const dimensionedScalar & D
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.
virtual void correctNut()
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)
BasicTurbulenceModel::transportModel transportModel
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Launder and Sharma low-Reynolds k-epsilon turbulence model for incompressible and compressible and co...
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)