Speziale, Sarkar and Gatski Reynolds-stress turbulence model for incompressible and compressible flows. More...
Public Types | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Types inherited from ReynoldsStress< RASModel< BasicMomentumTransportModel > > | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Types inherited from RASModel< BasicMomentumTransportModel > | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Member Functions | |
TypeName ("SSG") | |
Runtime type information. More... | |
SSG (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName) | |
Construct from components. More... | |
SSG (const SSG &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~SSG () |
Destructor. More... | |
virtual bool | read () |
Read model coefficients if they have changed. 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 tmp< volScalarField > | omega () const |
Return the turbulence specific dissipation rate. More... | |
tmp< volSymmTensorField > | DREff () const |
Return the effective diffusivity for R. More... | |
tmp< volSymmTensorField > | DepsilonEff () const |
Return the effective diffusivity for epsilon. More... | |
virtual void | correct () |
Solve the turbulence equations and correct eddy-Viscosity and. More... | |
void | operator= (const SSG &)=delete |
Disallow default bitwise assignment. More... | |
Public Member Functions inherited from ReynoldsStress< RASModel< BasicMomentumTransportModel > > | |
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 viscosity &viscosity) | |
Construct from components. More... | |
virtual | ~ReynoldsStress () |
Destructor. 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 tmp< volSymmTensorField > | sigma () const |
Return the Reynolds stress tensor [m^2/s^2]. More... | |
virtual tmp< volSymmTensorField > | devTau () const |
Return the effective stress tensor. More... | |
virtual tmp< fvVectorMatrix > | divDevTau (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevTau (const volScalarField &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
Public Member Functions inherited from RASModel< BasicMomentumTransportModel > | |
TypeName ("RAS") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, RASModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity),(alpha, rho, U, alphaRhoPhi, phi, viscosity)) | |
RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
Construct from components. More... | |
RASModel (const RASModel &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~RASModel () |
Destructor. More... | |
virtual const dictionary & | coeffDict () const |
Const access to the coefficients dictionary. More... | |
virtual tmp< volScalarField > | nu () const |
Return the laminar viscosity. More... | |
virtual tmp< scalarField > | nu (const label patchi) const |
Return the laminar viscosity on patchi. More... | |
virtual tmp< volScalarField > | nuEff () const |
Return the effective viscosity. More... | |
virtual tmp< scalarField > | nuEff (const label patchi) const |
Return the effective viscosity on patch. More... | |
virtual void | predict () |
Predict the turbulence transport coefficients if possible. More... | |
void | operator= (const RASModel &)=delete |
Disallow default bitwise assignment. More... | |
Protected Member Functions | |
tmp< volScalarField > | boundEpsilon () |
Bound epsilon and return Cmu*sqr(k) for nut. More... | |
virtual void | correctNut () |
Correct the eddy-viscosity nut. More... | |
virtual tmp< fvScalarMatrix > | epsilonSource () const |
Source term for the epsilon equation. More... | |
Protected Member Functions inherited from ReynoldsStress< RASModel< BasicMomentumTransportModel > > | |
void | boundNormalStress (volSymmTensorField &R) const |
void | correctWallShearStress (volSymmTensorField &R) const |
virtual tmp< fvSymmTensorMatrix > | RSource () const |
Source term for the R equation. More... | |
tmp< fvVectorMatrix > | DivDevRhoReff (const RhoFieldType &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
Protected Member Functions inherited from RASModel< BasicMomentumTransportModel > | |
virtual void | printCoeffs (const word &type) |
Print model coefficients. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from RASModel< BasicMomentumTransportModel > | |
static autoPtr< RASModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) |
Return a reference to the selected RAS model. More... | |
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 generalised 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 BasicMomentumTransportModel::alphaField alphaField |
SSG | ( | const alphaField & | alpha, |
const rhoField & | rho, | ||
const volVectorField & | U, | ||
const surfaceScalarField & | alphaRhoPhi, | ||
const surfaceScalarField & | phi, | ||
const viscosity & | viscosity, | ||
const word & | type = typeName |
||
) |
Construct from components.
Definition at line 75 of file SSG.C.
References SSG< BasicMomentumTransportModel >::boundEpsilon(), ReynoldsStress< RASModel< BasicMomentumTransportModel > >::boundNormalStress(), SSG< BasicMomentumTransportModel >::k_, RASModel< BasicMomentumTransportModel >::printCoeffs(), ReynoldsStress< RASModel< BasicMomentumTransportModel > >::R_, Foam::tr(), and Foam::type().
Disallow default bitwise copy construction.
|
protected |
Bound epsilon and return Cmu*sqr(k) for nut.
Definition at line 41 of file SSG.C.
References Foam::max(), and Foam::sqr().
Referenced by SSG< BasicMomentumTransportModel >::SSG().
|
protectedvirtual |
Correct the eddy-viscosity nut.
Implements ReynoldsStress< RASModel< BasicMomentumTransportModel > >.
Definition at line 50 of file SSG.C.
References dictionary::New().
|
protectedvirtual |
Source term for the epsilon equation.
Definition at line 59 of file SSG.C.
References Foam::dimTime, and Foam::dimVolume.
TypeName | ( | "SSG< BasicMomentumTransportModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Implements ReynoldsStress< RASModel< BasicMomentumTransportModel > >.
Definition at line 246 of file SSG.C.
References Foam::read().
|
inlinevirtual |
Return the turbulence kinetic energy.
Reimplemented from ReynoldsStress< RASModel< BasicMomentumTransportModel > >.
Definition at line 174 of file SSG.H.
References SSG< BasicMomentumTransportModel >::k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 180 of file SSG.H.
References SSG< BasicMomentumTransportModel >::epsilon_.
|
inlinevirtual |
Return the turbulence specific dissipation rate.
Definition at line 186 of file SSG.H.
References SSG< BasicMomentumTransportModel >::Cmu_, SSG< BasicMomentumTransportModel >::epsilon_, SSG< BasicMomentumTransportModel >::k_, and GeometricField< Type, PatchField, GeoMesh >::New().
tmp< volSymmTensorField > DREff |
Return the effective diffusivity for R.
Definition at line 274 of file SSG.C.
References Foam::I, and GeometricField< Type, PatchField, GeoMesh >::New().
tmp< volSymmTensorField > DepsilonEff |
Return the effective diffusivity for epsilon.
Definition at line 285 of file SSG.C.
References Foam::I, and GeometricField< Type, PatchField, GeoMesh >::New().
|
virtual |
Solve the turbulence equations and correct eddy-Viscosity and.
related properties
Implements ReynoldsStress< RASModel< BasicMomentumTransportModel > >.
Definition at line 296 of file SSG.C.
References alpha(), b, fvConstraints::constrain(), Foam::MULES::correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), fvPatch::faceCells(), forAll, fvConstraints(), fvModels(), Foam::constant::universal::G, Foam::fvc::grad(), Foam::I, Foam::innerSqr(), Foam::fvm::laplacian(), Foam::mag(), Foam::min(), dictionary::New(), patches, patchi, Foam::R(), tmp< T >::ref(), rho, Foam::fvm::S(), Foam::skew(), Foam::solve(), fvModels::source(), Foam::fvm::Sp(), Foam::symm(), Foam::tr(), Foam::twoSymm(), and U.
|
delete |
Disallow default bitwise assignment.
|
protected |
Definition at line 102 of file SSG.H.
Referenced by SSG< BasicMomentumTransportModel >::omega().
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Definition at line 119 of file SSG.H.
Referenced by SSG< BasicMomentumTransportModel >::k(), SSG< BasicMomentumTransportModel >::omega(), and SSG< BasicMomentumTransportModel >::SSG().
|
protected |
Definition at line 120 of file SSG.H.
Referenced by SSG< BasicMomentumTransportModel >::epsilon(), and SSG< BasicMomentumTransportModel >::omega().