Loading...
Searching...
No Matches
SSG< BasicTurbulenceModel > Class Template Reference

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

#include <SSG.H>

Inheritance diagram for SSG< BasicTurbulenceModel >:
Collaboration diagram for SSG< BasicTurbulenceModel >:

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< volScalarFieldk () const
 Return the turbulence kinetic energy.
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate.
tmp< volSymmTensorFieldDREff () const
 Return the effective diffusivity for R.
tmp< volSymmTensorFieldDepsilonEff () 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::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)
 Construct from components.
virtual ~ReynoldsStress ()=default
 Destructor.
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity.
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor.
virtual tmp< fvVectorMatrixdivDevRhoReff (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

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_

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

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 148 of file SSG.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 149 of file SSG.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 150 of file SSG.H.

Constructor & Destructor Documentation

◆ SSG()

template<class BasicTurbulenceModel>
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 )

◆ ~SSG()

template<class BasicTurbulenceModel>
virtual ~SSG ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

template<class BasicTurbulenceModel>
void correctNut ( )
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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "SSG< BasicTurbulenceModel >" )

Runtime type information.

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ k()

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Reimplemented from ReynoldsStress< RASModel< BasicTurbulenceModel > >.

Definition at line 193 of file SSG.H.

◆ epsilon()

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 201 of file SSG.H.

◆ DREff()

template<class BasicTurbulenceModel>
tmp< volSymmTensorField > DREff ( ) const

Return the effective diffusivity for R.

Definition at line 248 of file SSG.C.

References Cs_, epsilon_, Foam::I, k_, nu, and ReynoldsStress< RASModel< BasicTurbulenceModel > >::R_.

Referenced by correct().

Here is the caller graph for this function:

◆ DepsilonEff()

template<class BasicTurbulenceModel>
tmp< volSymmTensorField > DepsilonEff ( ) const

Return the effective diffusivity for epsilon.

Definition at line 262 of file SSG.C.

References Ceps_, epsilon_, Foam::I, k_, nu, and ReynoldsStress< RASModel< BasicTurbulenceModel > >::R_.

Referenced by correct().

Here is the caller graph for this function:

◆ correct()

Member Data Documentation

◆ Cmu_

template<class BasicTurbulenceModel>
dimensionedScalar Cmu_
protected

Definition at line 117 of file SSG.H.

Referenced by correctNut(), read(), and SSG().

◆ C1_

template<class BasicTurbulenceModel>
dimensionedScalar C1_
protected

Definition at line 119 of file SSG.H.

Referenced by correct(), read(), and SSG().

◆ C1s_

template<class BasicTurbulenceModel>
dimensionedScalar C1s_
protected

Definition at line 120 of file SSG.H.

Referenced by correct(), read(), and SSG().

◆ C2_

template<class BasicTurbulenceModel>
dimensionedScalar C2_
protected

Definition at line 121 of file SSG.H.

Referenced by correct(), read(), and SSG().

◆ C3_

template<class BasicTurbulenceModel>
dimensionedScalar C3_
protected

Definition at line 122 of file SSG.H.

Referenced by correct(), read(), and SSG().

◆ C3s_

template<class BasicTurbulenceModel>
dimensionedScalar C3s_
protected

Definition at line 123 of file SSG.H.

Referenced by correct(), read(), and SSG().

◆ C4_

template<class BasicTurbulenceModel>
dimensionedScalar C4_
protected

Definition at line 124 of file SSG.H.

Referenced by correct(), read(), and SSG().

◆ C5_

template<class BasicTurbulenceModel>
dimensionedScalar C5_
protected

Definition at line 125 of file SSG.H.

Referenced by correct(), read(), and SSG().

◆ Ceps1_

template<class BasicTurbulenceModel>
dimensionedScalar Ceps1_
protected

Definition at line 127 of file SSG.H.

Referenced by correct(), read(), and SSG().

◆ Ceps2_

template<class BasicTurbulenceModel>
dimensionedScalar Ceps2_
protected

Definition at line 128 of file SSG.H.

Referenced by correct(), read(), and SSG().

◆ Cs_

template<class BasicTurbulenceModel>
dimensionedScalar Cs_
protected

Definition at line 129 of file SSG.H.

Referenced by DREff(), read(), and SSG().

◆ Ceps_

template<class BasicTurbulenceModel>
dimensionedScalar Ceps_
protected

Definition at line 130 of file SSG.H.

Referenced by DepsilonEff(), read(), and SSG().

◆ k_

template<class BasicTurbulenceModel>
volScalarField k_
protected

Definition at line 134 of file SSG.H.

Referenced by correct(), correctNut(), DepsilonEff(), DREff(), and SSG().

◆ epsilon_

template<class BasicTurbulenceModel>
volScalarField epsilon_
protected

Definition at line 135 of file SSG.H.

Referenced by correct(), correctNut(), DepsilonEff(), DREff(), and SSG().


The documentation for this class was generated from the following files:
  • src/TurbulenceModels/turbulenceModels/RAS/SSG/SSG.H
  • src/TurbulenceModels/turbulenceModels/RAS/SSG/SSG.C