Go to the documentation of this file.
38 template<
class BasicTurbulenceModel>
41 this->nut_ = this->Cmu_*
sqr(k_)/epsilon_;
42 this->nut_.correctBoundaryConditions();
44 BasicTurbulenceModel::correctNut();
50 template<
class BasicTurbulenceModel>
59 const word& propertiesName,
190 this->runTime_.timeName(),
202 this->runTime_.timeName(),
210 if (
type == typeName)
212 this->printCoeffs(
type);
214 this->boundNormalStress(this->R_);
215 bound(epsilon_, this->epsilonMin_);
216 k_ = 0.5*
tr(this->R_);
223 template<
class BasicTurbulenceModel>
228 Cmu_.readIfPresent(this->coeffDict());
229 C1_.readIfPresent(this->coeffDict());
230 C1s_.readIfPresent(this->coeffDict());
231 C2_.readIfPresent(this->coeffDict());
232 C3_.readIfPresent(this->coeffDict());
233 C3s_.readIfPresent(this->coeffDict());
234 C4_.readIfPresent(this->coeffDict());
235 C5_.readIfPresent(this->coeffDict());
237 Ceps1_.readIfPresent(this->coeffDict());
238 Ceps2_.readIfPresent(this->coeffDict());
239 Cs_.readIfPresent(this->coeffDict());
240 Ceps_.readIfPresent(this->coeffDict());
251 template<
class BasicTurbulenceModel>
259 (Cs_*(this->k_/this->epsilon_))*this->R_ +
I*this->
nu()
265 template<
class BasicTurbulenceModel>
273 (Ceps_*(this->k_/this->epsilon_))*this->R_ +
I*this->
nu()
279 template<
class BasicTurbulenceModel>
282 if (!this->turbulence_)
303 epsilon_.boundaryField().updateCoeffs();
318 epsEqn().boundaryManipulate(epsilon_.boundaryField());
321 bound(epsilon_, this->epsilonMin_);
332 if (isA<wallFvPatch>(curPatch))
339 G[faceCelli]/(0.5*
mag(
tr(P[faceCelli])) + SMALL),
359 - ((1.0/3.0)*
I)*(((2.0 - C1_)*epsilon_ - C1s_*
G)*
alpha*
rho)
372 this->boundNormalStress(
R);
379 this->correctWallShearStress(
R);
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)
A class for handling words, derived from string.
dimensionedTensor skew(const dimensionedTensor &dt)
#define forAll(list, i)
Loop across all elements in list.
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.
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
#define R(A, B, C, D, E, F, K, M)
Reynolds-stress turbulence model base class.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
static const sphericalTensor I(1)
BasicTurbulenceModel::rhoField rhoField
Generic dimensioned Type class.
BasicTurbulenceModel::alphaField alphaField
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)
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void correctNut()
Update the eddy-viscosity.
virtual bool read()
Read model coefficients if they have changed.
virtual const labelUList & faceCells() const
Return faceCells.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)