50 #ifndef ReynoldsStress_H 51 #define ReynoldsStress_H 62 template<
class BasicTurbulenceModel>
65 public BasicTurbulenceModel
92 template<
class RhoFieldType>
95 const RhoFieldType&
rho,
102 typedef typename BasicTurbulenceModel::alphaField
alphaField;
103 typedef typename BasicTurbulenceModel::rhoField
rhoField;
104 typedef typename BasicTurbulenceModel::transportModel
transportModel;
112 const word& modelName,
119 const word& propertiesName
130 virtual bool read() = 0;
void correctWallShearStress(volSymmTensorField &R) const
BasicTurbulenceModel::alphaField alphaField
virtual void correctNut()=0
Update the eddy-viscosity.
tmp< fvVectorMatrix > DivDevRhoReff(const RhoFieldType &rho, volVectorField &U) const
Return the source term for the momentum equation.
virtual ~ReynoldsStress()=default
Destructor.
BasicTurbulenceModel::rhoField rhoField
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
BasicTurbulenceModel::transportModel transportModel
BasicTurbulenceModel::transportModel transportModel
void checkRealizabilityConditions(const volSymmTensorField &R) const
virtual bool read()=0
Re-read model coefficients if they have changed.
A class for handling words, derived from Foam::string.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
BasicTurbulenceModel::rhoField rhoField
virtual void validate()
Validate the turbulence fields after construction.
dimensionedScalar couplingFactor_
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity.
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)
Construct from components.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
BasicTurbulenceModel::alphaField alphaField
Reynolds-stress turbulence model base class.