SSG< BasicTurbulenceModel > Class Template Reference

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

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

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 Types inherited from RASModel< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 

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 transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
 Construct from components. More...
 
virtual ~SSG ()=default
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 
tmp< volSymmTensorFieldDREff () const
 Return the effective diffusivity for R. More...
 
tmp< volSymmTensorFieldDepsilonEff () const
 Return the effective diffusivity for epsilon. More...
 
virtual void correct ()
 Solve the turbulence equations and correct eddy-Viscosity and. More...
 
- 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. More...
 
virtual ~ReynoldsStress ()=default
 Destructor. More...
 
virtual bool read ()=0
 Re-read model coefficients if they have changed. More...
 
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity. More...
 
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (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...
 
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
- Public Member Functions inherited from RASModel< BasicTurbulenceModel >
 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 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)
 Construct from components. More...
 
virtual ~RASModel ()=default
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: SMALL) More...
 
const dimensionedScalarepsilonMin () const
 Return the lower allowable limit for epsilon (default: SMALL) More...
 
const dimensionedScalaromegaMin () const
 Return the lower allowable limit for omega (default: SMALL) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
dimensionedScalarepsilonMin ()
 Allow epsilonMin to be changed. More...
 
dimensionedScalaromegaMin ()
 Allow omegaMin to be changed. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 
virtual tmp< volScalarFieldomega () const
 Return the specific dissipation rate. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 

Protected Member Functions

virtual void correctNut ()
 Update the eddy-viscosity. More...
 
- Protected Member Functions inherited from ReynoldsStress< RASModel< BasicTurbulenceModel > >
void boundNormalStress (volSymmTensorField &R) const
 
void correctWallShearStress (volSymmTensorField &R) const
 
void checkRealizabilityConditions (const volSymmTensorField &R) const
 
virtual void correctNut ()=0
 Update the eddy-viscosity. More...
 
tmp< fvVectorMatrixDivDevRhoReff (const RhoFieldType &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
- Protected Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 
 RASModel (const RASModel &)=delete
 No copy construct. More...
 
void operator= (const RASModel &)=delete
 No copy assignment. More...
 

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_
 RAS coefficients dictionary. More...
 
Switch turbulence_
 Turbulence on/off flag. More...
 
Switch printCoeffs_
 Flag to print the model coeffs at run-time. More...
 
dictionary coeffDict_
 Model coefficients dictionary. More...
 
dimensionedScalar kMin_
 Lower limit of k. More...
 
dimensionedScalar epsilonMin_
 Lower limit of epsilon. More...
 
dimensionedScalar omegaMin_
 Lower limit for omega. More...
 

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)
 Return a reference to the selected RAS model. More...
 

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 98 of file SSG.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 146 of file SSG.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 147 of file SSG.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 148 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 
)

Construct from components.

Definition at line 56 of file SSG.C.

References Foam::bound(), ReynoldsStress< RASModel< BasicTurbulenceModel > >::boundNormalStress(), SSG< BasicTurbulenceModel >::epsilon_, RASModel< BasicTurbulenceModel >::epsilonMin_, SSG< BasicTurbulenceModel >::k_, RASModel< BasicTurbulenceModel >::printCoeffs(), ReynoldsStress< RASModel< BasicTurbulenceModel > >::R_, Foam::tr(), and Foam::type().

Here is the call graph for this function:

◆ ~SSG()

virtual ~SSG ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

void correctNut
protectedvirtual

Update the eddy-viscosity.

Implements ReynoldsStress< RASModel< BasicTurbulenceModel > >.

Definition at line 43 of file SSG.C.

References Time::New(), and Foam::sqr().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "SSG< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Read model coefficients if they have changed.

Implements ReynoldsStress< RASModel< BasicTurbulenceModel > >.

Definition at line 229 of file SSG.C.

References Foam::read().

Here is the call graph for this function:

◆ k()

virtual tmp< volScalarField > k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Reimplemented from ReynoldsStress< RASModel< BasicTurbulenceModel > >.

Definition at line 181 of file SSG.H.

References SSG< BasicTurbulenceModel >::k_.

◆ epsilon()

virtual tmp< volScalarField > epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Reimplemented from RASModel< BasicTurbulenceModel >.

Definition at line 187 of file SSG.H.

References SSG< BasicTurbulenceModel >::epsilon_.

◆ DREff()

Return the effective diffusivity for R.

Definition at line 255 of file SSG.C.

References Foam::I, and nu.

◆ DepsilonEff()

tmp< volSymmTensorField > DepsilonEff

Return the effective diffusivity for epsilon.

Definition at line 269 of file SSG.C.

References Foam::I, and nu.

◆ correct()

void correct
virtual

Solve the turbulence equations and correct eddy-Viscosity and.

related properties

Implements ReynoldsStress< RASModel< BasicTurbulenceModel > >.

Definition at line 283 of file SSG.C.

References alpha, b, Foam::bound(), correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), fvPatch::faceCells(), forAll, fvOptions, Foam::fvc::grad(), Foam::I, Foam::innerSqr(), Foam::fvm::laplacian(), Foam::mag(), Foam::min(), Time::New(), patches, R, tmp< T >::ref(), rho, Foam::skew(), solve(), Foam::fvm::Sp(), Foam::symm(), Foam::tr(), Foam::twoSymm(), and U.

Here is the call graph for this function:

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 117 of file SSG.H.

◆ C1_

dimensionedScalar C1_
protected

Definition at line 119 of file SSG.H.

◆ C1s_

dimensionedScalar C1s_
protected

Definition at line 120 of file SSG.H.

◆ C2_

dimensionedScalar C2_
protected

Definition at line 121 of file SSG.H.

◆ C3_

dimensionedScalar C3_
protected

Definition at line 122 of file SSG.H.

◆ C3s_

dimensionedScalar C3s_
protected

Definition at line 123 of file SSG.H.

◆ C4_

dimensionedScalar C4_
protected

Definition at line 124 of file SSG.H.

◆ C5_

dimensionedScalar C5_
protected

Definition at line 125 of file SSG.H.

◆ Ceps1_

dimensionedScalar Ceps1_
protected

Definition at line 127 of file SSG.H.

◆ Ceps2_

dimensionedScalar Ceps2_
protected

Definition at line 128 of file SSG.H.

◆ Cs_

dimensionedScalar Cs_
protected

Definition at line 129 of file SSG.H.

◆ Ceps_

dimensionedScalar Ceps_
protected

Definition at line 130 of file SSG.H.

◆ k_

volScalarField k_
protected

Definition at line 134 of file SSG.H.

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

◆ epsilon_

volScalarField epsilon_
protected

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