Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
SSG< BasicTurbulenceModel > Class Template Reference

Speziale, Sarkar and Gatski Reynolds-stress turbulence model for incompressible and compressible flows. More...

Inheritance diagram for SSG< BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for SSG< BasicTurbulenceModel >:
Collaboration graph
[legend]

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from ReynoldsStress< RASModel< BasicTurbulenceModel > >
typedef RASModel< BasicTurbulenceModel > ::alphaField alphaField
 
typedef RASModel< BasicTurbulenceModel > ::rhoField rhoField
 
typedef RASModel< BasicTurbulenceModel > ::transportModel transportModel
 
- Public Types inherited from RASModel< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 

Public Member Functions

 TypeName ("SSG")
 
 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)
 
virtual ~SSG ()=default
 
virtual bool read ()
 
virtual tmp< volScalarFieldk () const
 
virtual tmp< volScalarFieldepsilon () const
 
tmp< volSymmTensorFieldDREff () const
 
tmp< volSymmTensorFieldDepsilonEff () const
 
virtual void correct ()
 
- Public Member Functions inherited from ReynoldsStress< RASModel< BasicTurbulenceModel > >
Foam::tmp< Foam::fvVectorMatrixDivDevRhoReff (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)
 
virtual ~ReynoldsStress ()=default
 
virtual tmp< volScalarFieldnut () const
 
virtual tmp< scalarFieldnut (const label patchi) const
 
virtual tmp< volScalarFieldk () const
 
virtual tmp< volSymmTensorFieldR () const
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 
virtual void validate ()
 
- Public Member Functions inherited from RASModel< BasicTurbulenceModel >
 TypeName ("RAS")
 
 declareRunTimeSelectionTable (autoPtr, RASModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName))
 
 RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 
virtual ~RASModel ()=default
 
const dimensionedScalarkMin () const
 
const dimensionedScalarepsilonMin () const
 
const dimensionedScalaromegaMin () const
 
dimensionedScalarkMin ()
 
dimensionedScalarepsilonMin ()
 
dimensionedScalaromegaMin ()
 
virtual const dictionarycoeffDict () const
 
virtual tmp< volScalarFieldnuEff () const
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 
virtual tmp< volScalarFieldomega () const
 

Protected Member Functions

virtual void correctNut ()
 
- Protected Member Functions inherited from ReynoldsStress< RASModel< BasicTurbulenceModel > >
void boundNormalStress (volSymmTensorField &R) const
 
void correctWallShearStress (volSymmTensorField &R) const
 
tmp< fvVectorMatrixDivDevRhoReff (const RhoFieldType &rho, volVectorField &U) const
 
- Protected Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 
 RASModel (const RASModel &)=delete
 
void operator= (const RASModel &)=delete
 

Protected Attributes

dimensionedScalar Cmu_
 
dimensionedScalar C1_
 
dimensionedScalar C1s_
 
dimensionedScalar C2_
 
dimensionedScalar C3_
 
dimensionedScalar C3s_
 
dimensionedScalar C4_
 
dimensionedScalar C5_
 
dimensionedScalar Ceps1_
 
dimensionedScalar Ceps2_
 
dimensionedScalar Cs_
 
dimensionedScalar Ceps_
 
volScalarField k_
 
volScalarField epsilon_
 
- Protected Attributes inherited from ReynoldsStress< RASModel< BasicTurbulenceModel > >
dimensionedScalar couplingFactor_
 
volSymmTensorField R_
 
volScalarField nut_
 
- Protected Attributes inherited from RASModel< BasicTurbulenceModel >
dictionary RASDict_
 
Switch turbulence_
 
Switch printCoeffs_
 
dictionary coeffDict_
 
dimensionedScalar kMin_
 
dimensionedScalar epsilonMin_
 
dimensionedScalar omegaMin_
 

Additional Inherited Members

- Static Public Member Functions inherited from RASModel< BasicTurbulenceModel >
static autoPtr< RASModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
 

Detailed Description

template<class BasicTurbulenceModel>
class Foam::RASModels::SSG< BasicTurbulenceModel >

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;
    }
Source files

Definition at line 94 of file SSG.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 142 of file SSG.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 143 of file SSG.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 144 of file SSG.H.

Constructor & Destructor Documentation

◆ SSG()

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 
)

Definition at line 50 of file SSG.C.

References Foam::bound(), Foam::tr(), and Foam::type().

Here is the call graph for this function:

◆ ~SSG()

virtual ~SSG ( )
virtualdefault

Member Function Documentation

◆ correctNut()

void correctNut
protectedvirtual

Implements ReynoldsStress< RASModel< BasicTurbulenceModel > >.

Definition at line 36 of file SSG.C.

References optionList::correct(), options::New(), and Foam::sqr().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "SSG< BasicTurbulenceModel >"  )

◆ read()

bool read
virtual

Implements ReynoldsStress< RASModel< BasicTurbulenceModel > >.

Definition at line 222 of file SSG.C.

References Foam::read().

Here is the call graph for this function:

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Definition at line 177 of file SSG.H.

References SSG< BasicTurbulenceModel >::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Reimplemented from RASModel< BasicTurbulenceModel >.

Definition at line 183 of file SSG.H.

References SSG< BasicTurbulenceModel >::epsilon_.

◆ DREff()

Definition at line 248 of file SSG.C.

References Foam::I, and nu.

◆ DepsilonEff()

tmp< volSymmTensorField > DepsilonEff

Definition at line 262 of file SSG.C.

References Foam::I, and nu.

◆ correct()

void correct
virtual

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 113 of file SSG.H.

◆ C1_

dimensionedScalar C1_
protected

Definition at line 115 of file SSG.H.

◆ C1s_

dimensionedScalar C1s_
protected

Definition at line 116 of file SSG.H.

◆ C2_

dimensionedScalar C2_
protected

Definition at line 117 of file SSG.H.

◆ C3_

dimensionedScalar C3_
protected

Definition at line 118 of file SSG.H.

◆ C3s_

dimensionedScalar C3s_
protected

Definition at line 119 of file SSG.H.

◆ C4_

dimensionedScalar C4_
protected

Definition at line 120 of file SSG.H.

◆ C5_

dimensionedScalar C5_
protected

Definition at line 121 of file SSG.H.

◆ Ceps1_

dimensionedScalar Ceps1_
protected

Definition at line 123 of file SSG.H.

◆ Ceps2_

dimensionedScalar Ceps2_
protected

Definition at line 124 of file SSG.H.

◆ Cs_

dimensionedScalar Cs_
protected

Definition at line 125 of file SSG.H.

◆ Ceps_

dimensionedScalar Ceps_
protected

Definition at line 126 of file SSG.H.

◆ k_

volScalarField k_
protected

Definition at line 130 of file SSG.H.

Referenced by SSG< BasicTurbulenceModel >::k().

◆ epsilon_

volScalarField epsilon_
protected

Definition at line 131 of file SSG.H.

Referenced by SSG< BasicTurbulenceModel >::epsilon().


The documentation for this class was generated from the following files: