35 #ifndef ReynoldsStress_H 36 #define ReynoldsStress_H 47 template<
class BasicMomentumTransportModel>
50 public BasicMomentumTransportModel
76 template<
class RhoFieldType>
79 const RhoFieldType&
rho,
86 typedef typename BasicMomentumTransportModel::alphaField
alphaField;
87 typedef typename BasicMomentumTransportModel::rhoField
rhoField;
88 typedef typename BasicMomentumTransportModel::transportModel
transportModel;
96 const word& modelName,
97 const alphaField&
alpha,
102 const transportModel& transport
114 virtual bool read() = 0;
void correctWallShearStress(volSymmTensorField &R) const
virtual tmp< fvVectorMatrix > divDevTau(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)
Construct from components.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual tmp< volSymmTensorField > sigma() const
Return the Reynolds stress tensor [m^2/s^2].
BasicMomentumTransportModel::rhoField rhoField
BasicMomentumTransportModel::alphaField alphaField
virtual void correctNut()=0
Update the eddy-viscosity.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
tmp< fvVectorMatrix > DivDevRhoReff(const RhoFieldType &rho, volVectorField &U) const
Return the source term for the momentum equation.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
virtual bool read()=0
Re-read model coefficients if they have changed.
BasicMomentumTransportModel::transportModel transportModel
A class for handling words, derived from string.
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
virtual ~ReynoldsStress()
Destructor.
#define R(A, B, C, D, E, F, K, M)
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Reynolds-stress turbulence model base class.