42template<
class BasicTurbulenceModel>
46 this->
nut_.correctBoundaryConditions();
49 BasicTurbulenceModel::correctNut();
55template<
class BasicTurbulenceModel>
56LRR<BasicTurbulenceModel>::LRR
58 const alphaField&
alpha,
64 const word& propertiesName,
208 this->printCoeffs(
type);
219template<
class BasicTurbulenceModel>
224 Cmu_.readIfPresent(this->coeffDict());
225 C1_.readIfPresent(this->coeffDict());
226 C2_.readIfPresent(this->coeffDict());
227 Ceps1_.readIfPresent(this->coeffDict());
228 Ceps2_.readIfPresent(this->coeffDict());
229 Cs_.readIfPresent(this->coeffDict());
230 Ceps_.readIfPresent(this->coeffDict());
232 wallReflection_.readIfPresent(
"wallReflection", this->coeffDict());
233 kappa_.readIfPresent(this->coeffDict());
234 Cref1_.readIfPresent(this->coeffDict());
235 Cref2_.readIfPresent(this->coeffDict());
244template<
class BasicTurbulenceModel>
258template<
class BasicTurbulenceModel>
272template<
class BasicTurbulenceModel>
275 if (!this->turbulence_)
281 const alphaField&
alpha = this->alpha_;
282 const rhoField&
rho = this->rho_;
297 epsilon_.boundaryFieldRef().updateCoeffs();
299 epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
313 epsEqn.ref().relax();
315 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
318 bound(epsilon_, this->epsilonMin_);
333 label celli = curPatch.
faceCells()[facei];
336 G[celli]/(0.5*
mag(
tr(P[celli])) + SMALL),
344 tmp<fvSymmTensorMatrix> REqn
352 - (2.0/3.0*(1 - C1_)*
I)*
alpha*
rho*epsilon_
365 Cref1_*
R - ((Cref2_*C2_)*(k_/epsilon_))*
dev(P)
378 this->boundNormalStress(
R);
385 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,...
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
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
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
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)
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
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
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
#define forAll(list, i)
Loop across all elements in list.