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

The Smagorinsky SGS model. More...

Inheritance diagram for Smagorinsky< BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for Smagorinsky< BasicTurbulenceModel >:
Collaboration graph
[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 LESModel< BasicTurbulenceModel > ::alphaField alphaField
 
typedef LESModel< BasicTurbulenceModel > ::rhoField rhoField
 
typedef LESModel< BasicTurbulenceModel > ::transportModel transportModel
 
- Public Types inherited from linearViscousStress< LESModel< BasicTurbulenceModel > >
typedef LESModel< BasicTurbulenceModel > ::alphaField alphaField
 
typedef LESModel< BasicTurbulenceModel > ::rhoField rhoField
 
typedef LESModel< BasicTurbulenceModel > ::transportModel transportModel
 
- Public Types inherited from LESModel< 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 ()
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
virtual tmp< volScalarFieldk () const
 Return SGS kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return sub-grid disipation rate. 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 ()
 Destructor. 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 ()
 Destructor. 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< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
- Public Member Functions inherited from linearViscousStress< LESModel< 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 ()
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. 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...
 
- Public Member Functions inherited from LESModel< BasicTurbulenceModel >
 TypeName ("LES")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, LESModel, 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))
 
 LESModel (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 ~LESModel ()
 Destructor. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: small) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
const volScalarFielddelta () const
 Access function to filter width. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. 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...
 
- Protected Member Functions inherited from LESModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 

Protected Attributes

dimensionedScalar Ck_
 
- Protected Attributes inherited from LESeddyViscosity< BasicTurbulenceModel >
dimensionedScalar Ce_
 
- Protected Attributes inherited from eddyViscosity< LESModel< BasicTurbulenceModel > >
volScalarField nut_
 
- Protected Attributes inherited from LESModel< BasicTurbulenceModel >
dictionary LESDict_
 LES 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...
 
autoPtr< Foam::LESdeltadelta_
 Run-time selectable delta model. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from LESModel< BasicTurbulenceModel >
static autoPtr< LESModelNew (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 LES model. More...
 

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 86 of file Smagorinsky.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 116 of file Smagorinsky.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 117 of file Smagorinsky.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 118 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 83 of file Smagorinsky.C.

◆ ~Smagorinsky()

virtual ~Smagorinsky ( )
inlinevirtual

Destructor.

Definition at line 142 of file Smagorinsky.H.

References Smagorinsky< BasicTurbulenceModel >::read().

Here is the call graph for this function:

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

References Foam::constant::physicoChemical::b, Foam::constant::universal::c, 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 67 of file Smagorinsky.C.

References optionList::correct(), delta, Foam::fvc::grad(), k, options::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 126 of file Smagorinsky.C.

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

Here is the caller graph for this function:

◆ k() [2/2]

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return SGS kinetic energy.

Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.

Definition at line 152 of file Smagorinsky.H.

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

Here is the call graph for this function:

◆ epsilon()

tmp< volScalarField > epsilon ( ) const
virtual

Return sub-grid disipation rate.

Reimplemented from LESeddyViscosity< BasicTurbulenceModel >.

Definition at line 142 of file Smagorinsky.C.

References delta, Foam::fvc::grad(), IOobject::groupName(), k, and Foam::sqrt().

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

References eddyViscosity< LESModel< BasicTurbulenceModel > >::correct().

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

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

Member Data Documentation

◆ Ck_

dimensionedScalar Ck_
protected

Definition at line 101 of file Smagorinsky.H.


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