Speziale, Sarkar and Gatski Reynolds-stress turbulence model for incompressible and compressible flows. More...
#include <SSG.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from ReynoldsStress< RASModel< BasicTurbulenceModel > > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("SSG") | |
| Runtime type information. | |
| SSG (const alphaField &alpha, const rhoField &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. | |
| virtual | ~SSG ()=default |
| Destructor. | |
| virtual bool | read () |
| Read model coefficients if they have changed. | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy. | |
| virtual tmp< volScalarField > | epsilon () const |
| Return the turbulence kinetic energy dissipation rate. | |
| tmp< volSymmTensorField > | DREff () const |
| Return the effective diffusivity for R. | |
| tmp< volSymmTensorField > | DepsilonEff () const |
| Return the effective diffusivity for epsilon. | |
| virtual void | correct () |
| Solve the turbulence equations and correct eddy-Viscosity and. | |
| Public Member Functions inherited from ReynoldsStress< RASModel< BasicTurbulenceModel > > | |
| Foam::tmp< Foam::fvVectorMatrix > | DivDevRhoReff (const RhoFieldType &rho, volVectorField &U) 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 | ~ReynoldsStress ()=default |
| Destructor. | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity. | |
| virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual void | validate () |
| Validate the turbulence fields after construction. | |
Protected Member Functions | |
| virtual void | correctNut () |
| Update the eddy-viscosity. | |
| Protected Member Functions inherited from ReynoldsStress< RASModel< BasicTurbulenceModel > > | |
| void | boundNormalStress (volSymmTensorField &R) const |
| void | correctWallShearStress (volSymmTensorField &R) const |
| void | checkRealizabilityConditions (const volSymmTensorField &R) const |
Speziale, Sarkar and Gatski Reynolds-stress turbulence model for incompressible and compressible flows.
Reference:
Speziale, C. G., Sarkar, S., & Gatski, T. B. (1991).
Modelling the pressure-strain correlation of turbulence:
an invariant dynamical systems approach.
Journal of Fluid Mechanics, 227, 245-272.
Including the generalized gradient diffusion model of Daly and Harlow:
Daly, B. J., & Harlow, F. H. (1970).
Transport equations in turbulence.
Physics of Fluids (1958-1988), 13(11), 2634-2649.
The default model coefficients are:
SSGCoeffs
{
Cmu 0.09;
C1 3.4;
C1s 1.8;
C2 4.2;
C3 0.8;
C3s 1.3;
C4 1.25;
C5 0.4;
Ceps1 1.44;
Ceps2 1.92;
Cs 0.25;
Ceps 0.15;
couplingFactor 0.0;
}
| typedef BasicTurbulenceModel::alphaField alphaField |
| typedef BasicTurbulenceModel::rhoField rhoField |
| typedef BasicTurbulenceModel::transportModel transportModel |
| SSG | ( | const alphaField & | alpha, |
| const rhoField & | 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 49 of file SSG.C.
References alpha, Foam::bound(), ReynoldsStress< RASModel< BasicTurbulenceModel > >::boundNormalStress(), C1_, C1s_, C2_, C3_, C3s_, C4_, C5_, Ceps1_, Ceps2_, Ceps_, Cmu_, Cs_, epsilon_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, k_, phi, ReynoldsStress< RASModel< BasicTurbulenceModel > >::R_, ReynoldsStress< RASModel< BasicTurbulenceModel > >::ReynoldsStress(), rho, timeName, Foam::tr(), and U.

|
virtualdefault |
Destructor.
|
protectedvirtual |
Update the eddy-viscosity.
Implements ReynoldsStress< RASModel< BasicTurbulenceModel > >.
Definition at line 36 of file SSG.C.
References Cmu_, optionList::correct(), epsilon_, k_, options::New(), ReynoldsStress< RASModel< BasicTurbulenceModel > >::nut_, and Foam::sqr().
Referenced by correct().


| TypeName | ( | "SSG< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Implements ReynoldsStress< RASModel< BasicTurbulenceModel > >.
Definition at line 222 of file SSG.C.
References C1_, C1s_, C2_, C3_, C3s_, C4_, C5_, Ceps1_, Ceps2_, Ceps_, Cmu_, Cs_, read(), and ReynoldsStress< RASModel< BasicTurbulenceModel > >::ReynoldsStress().
Referenced by read().


|
inlinevirtual |
Return the turbulence kinetic energy.
Reimplemented from ReynoldsStress< RASModel< BasicTurbulenceModel > >.
|
inlinevirtual |
| tmp< volSymmTensorField > DREff | ( | ) | const |
| tmp< volSymmTensorField > DepsilonEff | ( | ) | const |
|
virtual |
Solve the turbulence equations and correct eddy-Viscosity and.
related properties
Implements ReynoldsStress< RASModel< BasicTurbulenceModel > >.
Definition at line 276 of file SSG.C.
References alpha, b, Foam::bound(), ReynoldsStress< RASModel< BasicTurbulenceModel > >::boundNormalStress(), C1_, C1s_, C2_, C3_, C3s_, C4_, C5_, Ceps1_, Ceps2_, correct(), correctNut(), ReynoldsStress< RASModel< BasicTurbulenceModel > >::correctWallShearStress(), Foam::fvm::ddt(), DepsilonEff(), Foam::dev(), Foam::devTwoSymm(), Foam::fvm::div(), DREff(), epsilon_, fvPatch::faceCells(), forAll, fvOptions, Foam::fvc::grad(), Foam::I, Foam::innerSqr(), Foam::isA(), k_, Foam::fvm::laplacian(), Foam::mag(), Foam::min(), options::New(), patches, ReynoldsStress< RASModel< BasicTurbulenceModel > >::R(), ReynoldsStress< RASModel< BasicTurbulenceModel > >::R_, tmp< T >::ref(), ReynoldsStress< RASModel< BasicTurbulenceModel > >::ReynoldsStress(), rho, Foam::skew(), solve(), Foam::fvm::Sp(), Foam::symm(), Foam::tr(), Foam::twoSymm(), and U.

|
protected |
Definition at line 117 of file SSG.H.
Referenced by correctNut(), read(), and SSG().
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Definition at line 130 of file SSG.H.
Referenced by DepsilonEff(), read(), and SSG().
|
protected |
Definition at line 134 of file SSG.H.
Referenced by correct(), correctNut(), DepsilonEff(), DREff(), and SSG().
|
protected |
Definition at line 135 of file SSG.H.
Referenced by correct(), correctNut(), DepsilonEff(), DREff(), and SSG().