Smagorinsky< BasicTurbulenceModel > Class Template Reference

The Smagorinsky SGS model. More...

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

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from LESeddyViscosity< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from eddyViscosity< LESModel< BasicTurbulenceModel > >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from linearViscousStress< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 

Public Member Functions

 TypeName ("Smagorinsky")
 Runtime type information. More...
 
 Smagorinsky (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 ~Smagorinsky ()=default
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
virtual tmp< volScalarFieldk () const
 Return SGS kinetic energy. More...
 
virtual void correct ()
 Correct Eddy-Viscosity and related properties. More...
 
- Public Member Functions inherited from LESeddyViscosity< BasicTurbulenceModel >
 LESeddyViscosity (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
 Construct from components. More...
 
virtual ~LESeddyViscosity ()=default
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
- Public Member Functions inherited from eddyViscosity< LESModel< BasicTurbulenceModel > >
 eddyViscosity (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 ~eddyViscosity ()=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=0
 Return the turbulence kinetic energy. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. 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 linearViscousStress< BasicTurbulenceModel >
 linearViscousStress (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 ~linearViscousStress ()=default
 Destructor. More...
 
virtual bool read ()=0
 Re-read model coefficients if they have changed. 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 correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity. More...
 

Protected Member Functions

tmp< volScalarFieldk (const tmp< volTensorField > &gradU) const
 Return SGS kinetic energy. More...
 
virtual void correctNut ()
 Update the SGS eddy viscosity. More...
 
virtual void correctNut ()=0
 

Protected Attributes

dimensionedScalar Ck_
 
- Protected Attributes inherited from eddyViscosity< LESModel< BasicTurbulenceModel > >
volScalarField nut_
 

Detailed Description

template<class BasicTurbulenceModel>
class Foam::LESModels::Smagorinsky< BasicTurbulenceModel >

The Smagorinsky SGS model.

Reference:

    Smagorinsky, J. (1963).
    General circulation experiments with the primitive equations: I.
    The basic experiment*.
    Monthly weather review, 91(3), 99-164.

The form of the Smagorinsky model implemented is obtained from the k-equation model assuming local equilibrium which provides estimates of both k and epsilon separate from the sub-grid scale viscosity:

    B = 2/3*k*I - 2*nuSgs*dev(D)

where

    D = symm(grad(U));
    k from D:B + Ce*k^3/2/delta = 0
    nuSgs = Ck*sqrt(k)*delta

The default model coefficients are

    SmagorinskyCoeffs
    {
        Ck                  0.094;
        Ce                  1.048;
    }
See also
Foam::LESModels::kEqn
Source files

Definition at line 92 of file Smagorinsky.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 124 of file Smagorinsky.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 125 of file Smagorinsky.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 126 of file Smagorinsky.H.

Constructor & Destructor Documentation

◆ Smagorinsky()

Smagorinsky ( 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 85 of file Smagorinsky.C.

References Foam::type().

Here is the call graph for this function:

◆ ~Smagorinsky()

virtual ~Smagorinsky ( )
virtualdefault

Destructor.

Member Function Documentation

◆ k() [1/2]

tmp< volScalarField > k ( const tmp< volTensorField > &  gradU) const
protected

Return SGS kinetic energy.

calculated from the given velocity gradient

Definition at line 42 of file Smagorinsky.C.

References b, D, delta, Foam::dev(), IOobject::groupName(), Foam::sqr(), Foam::sqrt(), Foam::symm(), and Foam::tr().

Here is the call graph for this function:

◆ correctNut()

void correctNut
protectedvirtual

Update the SGS eddy viscosity.

Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.

Reimplemented in SmagorinskyZhang< BasicTurbulenceModel >.

Definition at line 70 of file Smagorinsky.C.

References delta, Foam::fvc::grad(), k, Time::New(), and Foam::sqrt().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "Smagorinsky< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Read model coefficients if they have changed.

Reimplemented from LESeddyViscosity< BasicTurbulenceModel >.

Reimplemented in SmagorinskyZhang< BasicTurbulenceModel >.

Definition at line 129 of file Smagorinsky.C.

◆ k() [2/2]

virtual tmp< volScalarField > k ( ) const
inlinevirtual

Return SGS kinetic energy.

Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.

Definition at line 159 of file Smagorinsky.H.

References Foam::fvc::grad(), and Smagorinsky< BasicTurbulenceModel >::k().

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

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

◆ correct()

void correct
virtual

Correct Eddy-Viscosity and related properties.

Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.

Definition at line 143 of file Smagorinsky.C.

References kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct().

Here is the call graph for this function:

Member Data Documentation

◆ Ck_

dimensionedScalar Ck_
protected

Definition at line 109 of file Smagorinsky.H.


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