42template<
class BasicTurbulenceModel>
46 this->
nut_.correctBoundaryConditions();
49 BasicTurbulenceModel::correctNut();
55template<
class BasicTurbulenceModel>
56SSG<BasicTurbulenceModel>::SSG
58 const alphaField&
alpha,
64 const word& propertiesName,
217 this->printCoeffs(
type);
228template<
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());
254template<
class BasicTurbulenceModel>
268template<
class BasicTurbulenceModel>
282template<
class BasicTurbulenceModel>
285 if (!this->turbulence_)
291 const alphaField&
alpha = this->alpha_;
292 const rhoField&
rho = this->rho_;
307 epsilon_.boundaryFieldRef().updateCoeffs();
309 epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
323 epsEqn.ref().relax();
325 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
328 bound(epsilon_, this->epsilonMin_);
343 label celli = curPatch.
faceCells()[facei];
346 G[celli]/(0.5*
mag(
tr(P[celli])) + SMALL),
358 tmp<fvSymmTensorMatrix> REqn
366 - ((1.0/3.0)*
I)*(((2.0 - C1_)*epsilon_ - C1s_*G)*
alpha*
rho)
382 this->boundNormalStress(
R);
385 bound(k_, this->kMin_);
390 this->correctWallShearStress(
R);
#define R(A, B, C, D, E, F, K, M)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Templated abstract base class for RAS turbulence models.
BasicTurbulenceModel::alphaField alphaField
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
virtual void correctNut()
Update the eddy-viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Read model coefficients if they have changed.
Reynolds-stress turbulence model base class.
virtual tmp< volSymmTensorField > R() const
void boundNormalStress(volSymmTensorField &R) const
ReynoldsStress(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
void correctWallShearStress(volSymmTensorField &R) const
Generic dimensioned Type class.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual const labelUList & faceCells() const
Return faceCells.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
A class for managing temporary objects.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
const polyBoundaryMesh & patches
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
bool read(const char *buf, int32_t &val)
Same as readInt32.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
dimensionedTensor skew(const dimensionedTensor &dt)
#define forAll(list, i)
Loop across all elements in list.