Go to the documentation of this file.
33 template<
class BasicTurbulenceModel>
39 scalar kMin = this->kMin_.value();
58 template<
class BasicTurbulenceModel>
70 if (isA<wallFvPatch>(curPatch))
78 this->U_.boundaryField()[
patchi].snGrad()
82 = this->mesh_.Sf().boundaryField()[
patchi];
85 = this->mesh_.magSf().boundaryField()[
patchi];
91 = (faceAreas[facei]/magFaceAreas[facei])*snGradU[facei];
97 Rw[facei].xy() =
tauw.xy();
98 Rw[facei].xz() =
tauw.xz();
99 Rw[facei].yz() =
tauw.yz();
108 template<
class BasicTurbulenceModel>
111 const word& modelName,
118 const word& propertiesName
147 IOobject::groupName(
"R",
U.group()),
148 this->runTime_.timeName(),
160 IOobject::groupName(
"nut",
U.group()),
161 this->runTime_.timeName(),
169 if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
172 <<
"couplingFactor = " << couplingFactor_
173 <<
" is not in range 0 - 1" <<
nl
181 template<
class BasicTurbulenceModel>
188 template<
class BasicTurbulenceModel>
196 template<
class BasicTurbulenceModel>
206 template<
class BasicTurbulenceModel>
216 IOobject::groupName(
"devRhoReff", this->U_.group()),
217 this->runTime_.timeName(),
222 this->alpha_*this->rho_*R_
223 - (this->alpha_*this->rho_*this->nu())
230 template<
class BasicTurbulenceModel>
237 if (couplingFactor_.value() > 0.0)
243 this->alpha_*this->rho_*R_
250 (1.0 - couplingFactor_)*this->alpha_*this->rho_*this->
nut(),
262 fvc::div(this->alpha_*this->rho_*R_)
265 this->alpha_*this->rho_*this->
nut(),
282 template<
class BasicTurbulenceModel>
298 template<
class BasicTurbulenceModel>
305 template<
class BasicTurbulenceModel>
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from string.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
virtual bool read()=0
Re-read model coefficients if they have changed.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool read(const char *, int32_t &)
void correctWallShearStress(volSymmTensorField &R) const
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
#define R(A, B, C, D, E, F, K, M)
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
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)
Construct from components.
Pre-declare SubField and related Field type.
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...
BasicTurbulenceModel::rhoField rhoField
Generic dimensioned Type class.
virtual void validate()
Validate the turbulence fields after construction.
BasicTurbulenceModel::alphaField alphaField
errorManipArg< error, int > exit(error &err, const int errNo=1)
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
BasicTurbulenceModel::transportModel transportModel
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
void boundNormalStress(volSymmTensorField &R) const
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Generic GeometricField class.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)