Go to the documentation of this file.
42 template<
class BasicTurbulenceModel>
45 this->nut_ = this->Cmu_*
sqr(k_)/epsilon_;
46 this->nut_.correctBoundaryConditions();
49 BasicTurbulenceModel::correctNut();
55 template<
class BasicTurbulenceModel>
64 const word& propertiesName,
195 this->runTime_.timeName(),
207 this->runTime_.timeName(),
215 if (
type == typeName)
217 this->printCoeffs(
type);
219 this->boundNormalStress(this->R_);
220 bound(epsilon_, this->epsilonMin_);
221 k_ = 0.5*
tr(this->R_);
228 template<
class BasicTurbulenceModel>
233 Cmu_.readIfPresent(this->coeffDict());
234 C1_.readIfPresent(this->coeffDict());
235 C1s_.readIfPresent(this->coeffDict());
236 C2_.readIfPresent(this->coeffDict());
237 C3_.readIfPresent(this->coeffDict());
238 C3s_.readIfPresent(this->coeffDict());
239 C4_.readIfPresent(this->coeffDict());
240 C5_.readIfPresent(this->coeffDict());
242 Ceps1_.readIfPresent(this->coeffDict());
243 Ceps2_.readIfPresent(this->coeffDict());
244 Cs_.readIfPresent(this->coeffDict());
245 Ceps_.readIfPresent(this->coeffDict());
254 template<
class BasicTurbulenceModel>
262 (Cs_*(this->k_/this->epsilon_))*this->R_ +
I*this->
nu()
268 template<
class BasicTurbulenceModel>
276 (Ceps_*(this->k_/this->epsilon_))*this->R_ +
I*this->
nu()
282 template<
class BasicTurbulenceModel>
285 if (!this->turbulence_)
307 epsilon_.boundaryFieldRef().updateCoeffs();
321 epsEqn.
ref().relax();
323 epsEqn.
ref().boundaryManipulate(epsilon_.boundaryFieldRef());
326 bound(epsilon_, this->epsilonMin_);
337 if (isA<wallFvPatch>(curPatch))
341 label celli = curPatch.
faceCells()[facei];
344 G[celli]/(0.5*
mag(
tr(P[celli])) + SMALL),
364 - ((1.0/3.0)*
I)*(((2.0 - C1_)*epsilon_ - C1s_*
G)*
alpha*
rho)
380 this->boundNormalStress(
R);
383 bound(k_, this->kMin_);
388 this->correctWallShearStress(
R);
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)
A class for handling words, derived from Foam::string.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
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.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
bool read(const char *buf, int32_t &val)
Same as readInt32.
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
#define forAll(list, i)
Loop across all elements in list.
#define R(A, B, C, D, E, F, K, M)
Reynolds-stress turbulence model base class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
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.
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
BasicTurbulenceModel::rhoField rhoField
Generic dimensioned Type class.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::transportModel transportModel
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.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
const polyBoundaryMesh & patches
virtual bool read()
Read model coefficients if they have changed.
virtual const labelUList & faceCells() const
Return faceCells.
Speziale, Sarkar and Gatski Reynolds-stress turbulence model for incompressible and compressible flow...
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)