Shih's quadratic algebraic Reynolds stress k-epsilon turbulence model for incompressible flows. More...
Public Member Functions | |
TypeName ("ShihQuadraticKE") | |
Runtime type information. More... | |
ShihQuadraticKE (const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName) | |
Construct from components. More... | |
virtual | ~ShihQuadraticKE () |
Destructor. More... | |
virtual bool | read () |
Read RASProperties dictionary. More... | |
tmp< volScalarField > | DkEff () const |
Return the effective diffusivity for k. More... | |
tmp< volScalarField > | DepsilonEff () const |
Return the effective diffusivity for epsilon. More... | |
virtual tmp< volScalarField > | k () const |
Return the turbulence kinetic energy. More... | |
virtual tmp< volScalarField > | epsilon () const |
Return the turbulence kinetic energy dissipation rate. More... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
Public Member Functions inherited from nonlinearEddyViscosity< incompressible::RASModel > | |
nonlinearEddyViscosity (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. More... | |
virtual | ~nonlinearEddyViscosity () |
Destructor. More... | |
virtual tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
Public Member Functions inherited from eddyViscosity< incompressible::RASModel > | |
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. More... | |
virtual | ~eddyViscosity () |
Destructor. More... | |
volScalarField & | evNut () |
Return non-const access to the turbulence viscosity. More... | |
virtual tmp< volScalarField > | nut () const |
Return the turbulence viscosity. More... | |
virtual tmp< scalarField > | nut (const label patchi) const |
Return the turbulence viscosity on patch. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
Public Member Functions inherited from linearViscousStress< incompressible::RASModel > | |
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. More... | |
virtual | ~linearViscousStress () |
Destructor. More... | |
Protected Member Functions | |
virtual void | correctNut () |
virtual void | correctNonlinearStress (const volTensorField &gradU) |
Protected Attributes | |
dimensionedScalar | Ceps1_ |
dimensionedScalar | Ceps2_ |
dimensionedScalar | sigmak_ |
dimensionedScalar | sigmaEps_ |
dimensionedScalar | Cmu1_ |
dimensionedScalar | Cmu2_ |
dimensionedScalar | Cbeta_ |
dimensionedScalar | Cbeta1_ |
dimensionedScalar | Cbeta2_ |
dimensionedScalar | Cbeta3_ |
volScalarField | k_ |
volScalarField | epsilon_ |
Protected Attributes inherited from nonlinearEddyViscosity< incompressible::RASModel > | |
volSymmTensorField | nonlinearStress_ |
Protected Attributes inherited from eddyViscosity< incompressible::RASModel > | |
volScalarField | nut_ |
Additional Inherited Members | |
Public Types inherited from nonlinearEddyViscosity< incompressible::RASModel > | |
typedef incompressible::RASModel::alphaField | alphaField |
typedef incompressible::RASModel::rhoField | rhoField |
typedef incompressible::RASModel::transportModel | transportModel |
Public Types inherited from eddyViscosity< incompressible::RASModel > | |
typedef incompressible::RASModel::alphaField | alphaField |
typedef incompressible::RASModel::rhoField | rhoField |
typedef incompressible::RASModel::transportModel | transportModel |
Public Types inherited from linearViscousStress< incompressible::RASModel > | |
typedef incompressible::RASModel::alphaField | alphaField |
typedef incompressible::RASModel::rhoField | rhoField |
typedef incompressible::RASModel::transportModel | transportModel |
Shih's quadratic algebraic Reynolds stress k-epsilon turbulence model for incompressible flows.
This turbulence model is described in:
Shih, T. H., Zhu, J., & Lumley, J. L. (1993). A realizable Reynolds stress algebraic equation model. NASA technical memorandum 105993.
Implemented according to the specification in: Apsley: Turbulence Models 2002
Definition at line 69 of file ShihQuadraticKE.H.
ShihQuadraticKE | ( | const geometricOneField & | alpha, |
const geometricOneField & | rho, | ||
const volVectorField & | U, | ||
const surfaceScalarField & | alphaRhoPhi, | ||
const surfaceScalarField & | phi, | ||
const transportModel & | transport, | ||
const word & | propertiesName = turbulenceModel::propertiesName , |
||
const word & | type = typeName |
||
) |
Construct from components.
Definition at line 80 of file ShihQuadraticKE.C.
References Foam::bound(), ShihQuadraticKE::epsilon_, and ShihQuadraticKE::k_.
Referenced by ShihQuadraticKE::correctNonlinearStress().
|
inlinevirtual |
Destructor.
Definition at line 127 of file ShihQuadraticKE.H.
References ShihQuadraticKE::read().
|
protectedvirtual |
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 48 of file ShihQuadraticKE.C.
References ShihQuadraticKE::correctNonlinearStress(), and Foam::fvc::grad().
|
protectedvirtual |
Implements nonlinearEddyViscosity< incompressible::RASModel >.
Definition at line 54 of file ShihQuadraticKE.C.
References ShihQuadraticKE::Cbeta1_, ShihQuadraticKE::Cbeta2_, ShihQuadraticKE::Cbeta3_, ShihQuadraticKE::Cbeta_, ShihQuadraticKE::Cmu1_, ShihQuadraticKE::Cmu2_, GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dev(), ShihQuadraticKE::epsilon_, Foam::innerSqr(), ShihQuadraticKE::k_, Foam::mag(), nonlinearEddyViscosity< incompressible::RASModel >::nonlinearStress_, eddyViscosity< incompressible::RASModel >::nut_, Foam::pow3(), ShihQuadraticKE::ShihQuadraticKE(), Foam::skew(), Foam::sqr(), Foam::sqrt(), Foam::symm(), and Foam::twoSymm().
Referenced by ShihQuadraticKE::correct(), and ShihQuadraticKE::correctNut().
TypeName | ( | "ShihQuadraticKE" | ) |
Runtime type information.
|
virtual |
Read RASProperties dictionary.
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 232 of file ShihQuadraticKE.C.
References ShihQuadraticKE::Cbeta1_, ShihQuadraticKE::Cbeta2_, ShihQuadraticKE::Cbeta3_, ShihQuadraticKE::Cbeta_, ShihQuadraticKE::Ceps1_, ShihQuadraticKE::Ceps2_, ShihQuadraticKE::Cmu1_, ShihQuadraticKE::Cmu2_, dimensioned< Type >::readIfPresent(), ShihQuadraticKE::sigmaEps_, and ShihQuadraticKE::sigmak_.
Referenced by ShihQuadraticKE::~ShihQuadraticKE().
|
inline |
Return the effective diffusivity for k.
Definition at line 137 of file ShihQuadraticKE.H.
References nu, and eddyViscosity< incompressible::RASModel >::nut_.
Referenced by ShihQuadraticKE::correct().
|
inline |
Return the effective diffusivity for epsilon.
Definition at line 146 of file ShihQuadraticKE.H.
References nu, and eddyViscosity< incompressible::RASModel >::nut_.
Referenced by ShihQuadraticKE::correct().
|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 155 of file ShihQuadraticKE.H.
References ShihQuadraticKE::k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 161 of file ShihQuadraticKE.H.
References ShihQuadraticKE::correct(), and ShihQuadraticKE::epsilon_.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 256 of file ShihQuadraticKE.C.
References Foam::bound(), GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), ShihQuadraticKE::Ceps1_, ShihQuadraticKE::Ceps2_, eddyViscosity< BasicTurbulenceModel >::correct(), ShihQuadraticKE::correctNonlinearStress(), Foam::fvm::ddt(), ShihQuadraticKE::DepsilonEff(), Foam::fvm::div(), ShihQuadraticKE::DkEff(), ShihQuadraticKE::epsilon_, Foam::constant::universal::G, Foam::fvc::grad(), ShihQuadraticKE::k_, Foam::fvm::laplacian(), nonlinearEddyViscosity< incompressible::RASModel >::nonlinearStress_, eddyViscosity< incompressible::RASModel >::nut_, tmp< T >::ref(), Foam::solve(), Foam::fvm::Sp(), and Foam::twoSymm().
Referenced by ShihQuadraticKE::epsilon().
|
protected |
Definition at line 80 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::correct(), and ShihQuadraticKE::read().
|
protected |
Definition at line 81 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::correct(), and ShihQuadraticKE::read().
|
protected |
Definition at line 82 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::read().
|
protected |
Definition at line 83 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::read().
|
protected |
Definition at line 84 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::correctNonlinearStress(), and ShihQuadraticKE::read().
|
protected |
Definition at line 85 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::correctNonlinearStress(), and ShihQuadraticKE::read().
|
protected |
Definition at line 86 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::correctNonlinearStress(), and ShihQuadraticKE::read().
|
protected |
Definition at line 87 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::correctNonlinearStress(), and ShihQuadraticKE::read().
|
protected |
Definition at line 88 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::correctNonlinearStress(), and ShihQuadraticKE::read().
|
protected |
Definition at line 89 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::correctNonlinearStress(), and ShihQuadraticKE::read().
|
protected |
Definition at line 94 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::correct(), ShihQuadraticKE::correctNonlinearStress(), ShihQuadraticKE::k(), and ShihQuadraticKE::ShihQuadraticKE().
|
protected |
Definition at line 95 of file ShihQuadraticKE.H.
Referenced by ShihQuadraticKE::correct(), ShihQuadraticKE::correctNonlinearStress(), ShihQuadraticKE::epsilon(), and ShihQuadraticKE::ShihQuadraticKE().