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

Differential SGS Stress Equation Model for incompressible and compressible flows. More...

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

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from ReynoldsStress< 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 ("DeardorffDiffStress")
 
 DeardorffDiffStress (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 ~DeardorffDiffStress ()=default
 
virtual bool read ()
 
virtual void correct ()
 
- Public Member Functions inherited from ReynoldsStress< LESModel< 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 LESModel< BasicTurbulenceModel >
 TypeName ("LES")
 
 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)
 
virtual ~LESModel ()=default
 
virtual const dictionarycoeffDict () const
 
const dimensionedScalarCe () const noexcept
 
const dimensionedScalarkMin () const
 
dimensionedScalarkMin ()
 
const volScalarFielddelta () const
 
virtual tmp< volScalarFieldnuEff () const
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 
virtual tmp< volScalarFieldepsilon () const
 
virtual tmp< volScalarFieldomega () const
 

Protected Member Functions

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

Protected Attributes

dimensionedScalar Ck_
 
dimensionedScalar Cm_
 
dimensionedScalar Ce_
 
dimensionedScalar Cs_
 
- Protected Attributes inherited from ReynoldsStress< LESModel< BasicTurbulenceModel > >
dimensionedScalar couplingFactor_
 
volSymmTensorField R_
 
volScalarField nut_
 
- Protected Attributes inherited from LESModel< BasicTurbulenceModel >
dictionary LESDict_
 
Switch turbulence_
 
Switch printCoeffs_
 
dictionary coeffDict_
 
dimensionedScalar Ce_
 
dimensionedScalar kMin_
 
dimensionedScalar epsilonMin_
 
dimensionedScalar omegaMin_
 
autoPtr< Foam::LESdeltadelta_
 

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)
 

Detailed Description

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

Differential SGS Stress Equation Model for incompressible and compressible flows.

Reference:

    Deardorff, J. W. (1973).
    The use of subgrid transport equations in a three-dimensional model
    of atmospheric turbulence.
    Journal of Fluids Engineering, 95(3), 429-438.

This SGS model uses a full balance equation for the SGS stress tensor to simulate the behaviour of B.

This implementation is as described in the above paper except that the triple correlation model of Donaldson is replaced with 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.

with the default value for the coefficient Cs of 0.25 from

    Launder, B. E., Reece, G. J., & Rodi, W. (1975).
    Progress in the development of a Reynolds-stress turbulence closure.
    Journal of fluid mechanics, 68(03), 537-566.
Source files

Definition at line 81 of file DeardorffDiffStress.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 114 of file DeardorffDiffStress.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 115 of file DeardorffDiffStress.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 116 of file DeardorffDiffStress.H.

Constructor & Destructor Documentation

◆ DeardorffDiffStress()

DeardorffDiffStress ( 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 DeardorffDiffStress.C.

References Foam::type().

Here is the call graph for this function:

◆ ~DeardorffDiffStress()

virtual ~DeardorffDiffStress ( )
virtualdefault

Member Function Documentation

◆ correctNut()

void correctNut
protectedvirtual

Implements ReynoldsStress< LESModel< BasicTurbulenceModel > >.

Definition at line 36 of file DeardorffDiffStress.C.

References optionList::correct(), delta, k, options::New(), and Foam::sqrt().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "DeardorffDiffStress< BasicTurbulenceModel >"  )

◆ read()

bool read
virtual

Implements ReynoldsStress< LESModel< BasicTurbulenceModel > >.

Definition at line 121 of file DeardorffDiffStress.C.

References Foam::read().

Here is the call graph for this function:

◆ correct()

void correct
virtual

Member Data Documentation

◆ Ck_

dimensionedScalar Ck_
protected

Definition at line 100 of file DeardorffDiffStress.H.

◆ Cm_

dimensionedScalar Cm_
protected

Definition at line 101 of file DeardorffDiffStress.H.

◆ Ce_

dimensionedScalar Ce_
protected

Definition at line 102 of file DeardorffDiffStress.H.

◆ Cs_

dimensionedScalar Cs_
protected

Definition at line 103 of file DeardorffDiffStress.H.


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