41#ifndef eddyViscosity_H
42#define eddyViscosity_H
55template<
class BasicTurbulenceModel>
77 typedef typename BasicTurbulenceModel::alphaField
alphaField;
78 typedef typename BasicTurbulenceModel::rhoField
rhoField;
79 typedef typename BasicTurbulenceModel::transportModel
transportModel;
87 const word& modelName,
94 const word& propertiesName
105 virtual bool read() = 0;
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
BasicTurbulenceModel::alphaField alphaField
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< scalarField > nut(const label patchi) const
Return the turbulence viscosity on patch.
virtual void validate()
Validate the turbulence fields after construction.
eddyViscosity(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 bool read()=0
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity.
virtual void correctNut()=0
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual ~eddyViscosity()=default
Destructor.
linearViscousStress(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.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField