36template<
class BasicTurbulenceModel>
42 const scalar kMin = this->kMin_.value();
56template<
class BasicTurbulenceModel>
74 const scalarField& nutw = this->nut_.boundaryField()[patchi];
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];
96 Rw[facei] = -nutw[facei]*2*
devSymm(gradUw);
103template<
class BasicTurbulenceModel>
106 const volSymmTensorField&
R
109 const label maxDiffs = 5;
117 if (maxDiffs < nDiffs)
119 Info<<
"More than " << maxDiffs <<
" times"
120 <<
" Reynolds-stress realizability checks failed."
121 <<
" Skipping further comparisons." <<
endl;
130 <<
"Reynolds stress " << r <<
" at cell " << celli
131 <<
" does not obey the constraint: Rxx >= 0"
139 <<
"Reynolds stress " << r <<
" at cell " << celli
140 <<
" does not obey the constraint: Rxx*Ryy - sqr(Rxy) >= 0"
148 <<
"Reynolds stress " << r <<
" at cell " << celli
149 <<
" does not obey the constraint: det(R) >= 0"
164template<
class BasicTurbulenceModel>
167 const word& modelName,
174 const word& propertiesName
204 this->runTime_.timeName(),
217 this->runTime_.timeName(),
229 <<
" is not in range 0 - 1" <<
nl
237template<
class BasicTurbulenceModel>
240 return BasicTurbulenceModel::read();
244template<
class BasicTurbulenceModel>
252template<
class BasicTurbulenceModel>
257 tk.ref().rename(
"k");
262template<
class BasicTurbulenceModel>
270template<
class BasicTurbulenceModel>
284 this->runTime_.timeName(),
289 this->alpha_*this->rho_*R_
290 - (this->alpha_*this->rho_*this->nu())
297template<
class BasicTurbulenceModel>
298template<
class RhoFieldType>
302 const RhoFieldType&
rho,
306 if (couplingFactor_.value() > 0.0)
312 (1.0 - couplingFactor_)*this->alpha_*
rho*this->
nut(),
333 this->alpha_*
rho*this->
nut(),
337 + fvc::div(this->alpha_*
rho*R_)
338 - fvc::div(this->alpha_*
rho*this->
nu()*dev2(
T(fvc::grad(
U))))
339 - fvm::laplacian(this->alpha_*
rho*this->nuEff(),
U)
345template<
class BasicTurbulenceModel>
356template<
class BasicTurbulenceModel>
368template<
class BasicTurbulenceModel>
375template<
class BasicTurbulenceModel>
378 BasicTurbulenceModel::correct();
#define R(A, B, C, D, E, F, K, M)
GeometricBoundaryField< symmTensor, fvPatchField, volMesh > Boundary
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word group(const word &name)
Return group (extension part of name).
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
BasicTurbulenceModel::alphaField alphaField
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
void boundNormalStress(volSymmTensorField &R) const
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual void validate()
Validate the turbulence fields after construction.
virtual bool read()=0
Re-read model coefficients if they have changed.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
tmp< fvVectorMatrix > DivDevRhoReff(const RhoFieldType &rho, volVectorField &U) const
Return the source term for the momentum equation.
virtual void correctNut()=0
Update the eddy-viscosity.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
BasicTurbulenceModel::transportModel transportModel
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.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
dimensionedScalar couplingFactor_
void correctWallShearStress(volSymmTensorField &R) const
void checkRealizabilityConditions(const volSymmTensorField &R) const
const Cmpt & yy() const noexcept
const Cmpt & xx() const noexcept
const Cmpt & xy() const noexcept
Generic dimensioned Type class.
static dimensioned< Type > getOrAddToDict(const word &name, dictionary &dict, const dimensionSet &dims=dimless, const Type &deflt=Type(Zero))
Construct dimensioned from dictionary, with default value.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
A class for managing temporary objects.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
A class for handling words, derived from Foam::string.
const polyBoundaryMesh & patches
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
dimensionedScalar det(const dimensionedSphericalTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Field< vector > vectorField
Specialisation of Field<T> for vector.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.