LESModel< BasicTurbulenceModel > Class Template Reference

Templated abstract base class for LES SGS models. More...

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

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 

Public Member Functions

 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 ()=default
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
const dimensionedScalarCe () const noexcept
 
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...
 
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...
 

Static Public Member Functions

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...
 

Protected Member Functions

virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 
 LESModel (const LESModel &)=delete
 No copy construct. More...
 
void operator= (const LESModel &)=delete
 No copy assignment. More...
 

Protected Attributes

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 Ce_
 Empirical model constant. 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...
 

Detailed Description

template<class BasicTurbulenceModel>
class Foam::LESModel< BasicTurbulenceModel >

Templated abstract base class for LES SGS models.

Source files

Definition at line 59 of file LESModel.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 109 of file LESModel.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 110 of file LESModel.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 111 of file LESModel.H.

Constructor & Destructor Documentation

◆ LESModel() [1/2]

LESModel ( const LESModel< BasicTurbulenceModel > &  )
protecteddelete

No copy construct.

◆ LESModel() [2/2]

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.

Definition at line 46 of file LESModel.C.

References Foam::dimless, Foam::dimTime, Foam::dimVelocity, Foam::New(), and Foam::sqr().

Here is the call graph for this function:

◆ ~LESModel()

virtual ~LESModel ( )
virtualdefault

Destructor.

Member Function Documentation

◆ printCoeffs()

void printCoeffs ( const word type)
protectedvirtual

Print model coefficients.

Definition at line 34 of file LESModel.C.

References Foam::endl(), and Foam::Info.

Referenced by DeardorffDiffStress< BasicTurbulenceModel >::DeardorffDiffStress().

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

◆ operator=()

void operator= ( const LESModel< BasicTurbulenceModel > &  )
protecteddelete

No copy assignment.

◆ TypeName()

TypeName ( "LES"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
LESModel< BasicTurbulenceModel >  ,
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)   
)

◆ New()

Foam::autoPtr< Foam::LESModel< BasicTurbulenceModel > > New ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport,
const word propertiesName = turbulenceModel::propertiesName 
)
static

Return a reference to the selected LES model.

Definition at line 132 of file LESModel.C.

References dict, dictionary::getCompat(), IOobject::group(), IOobject::groupName(), IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE, dictionary::subDict(), and U.

Here is the call graph for this function:

◆ read()

bool read
virtual

Read model coefficients if they have changed.

Reimplemented in DeardorffDiffStress< BasicTurbulenceModel >, linearViscousStress< LESModel< BasicTurbulenceModel > >, and ReynoldsStress< LESModel< BasicTurbulenceModel > >.

Definition at line 189 of file LESModel.C.

References Foam::type().

Here is the call graph for this function:

◆ coeffDict()

virtual const dictionary & coeffDict ( ) const
inlinevirtual

Const access to the coefficients dictionary.

Definition at line 182 of file LESModel.H.

References LESModel< BasicTurbulenceModel >::coeffDict_.

◆ Ce()

const dimensionedScalar & Ce ( ) const
inlinenoexcept

Definition at line 188 of file LESModel.H.

References LESModel< BasicTurbulenceModel >::Ce_.

◆ kMin() [1/2]

const dimensionedScalar & kMin ( ) const
inline

Return the lower allowable limit for k (default: SMALL)

Definition at line 194 of file LESModel.H.

References LESModel< BasicTurbulenceModel >::kMin_.

◆ kMin() [2/2]

dimensionedScalar & kMin ( )
inline

Allow kMin to be changed.

Definition at line 200 of file LESModel.H.

References LESModel< BasicTurbulenceModel >::kMin_.

◆ delta()

const volScalarField & delta ( ) const
inline

Access function to filter width.

Definition at line 206 of file LESModel.H.

References LESModel< BasicTurbulenceModel >::delta_.

◆ nuEff() [1/2]

virtual tmp< volScalarField > nuEff ( ) const
inlinevirtual

Return the effective viscosity.

Definition at line 213 of file LESModel.H.

References IOobject::groupName().

Here is the call graph for this function:

◆ nuEff() [2/2]

virtual tmp< scalarField > nuEff ( const label  patchi) const
inlinevirtual

Return the effective viscosity on patch.

Definition at line 226 of file LESModel.H.

References nu, and nut.

◆ epsilon()

Foam::tmp< Foam::volScalarField > epsilon
virtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 213 of file LESModel.C.

References delta, IOobject::groupName(), k, Time::New(), and Foam::pow().

Here is the call graph for this function:

◆ omega()

Foam::tmp< Foam::volScalarField > omega
virtual

Return the specific dissipation rate.

Definition at line 230 of file LESModel.C.

References Foam::dimLength, Foam::dimTime, IOobject::groupName(), Time::New(), and Foam::sqr().

Here is the call graph for this function:

◆ correct()

void correct
virtual

Solve the turbulence equations and correct the turbulence viscosity.

Reimplemented in DeardorffDiffStress< BasicTurbulenceModel >, linearViscousStress< LESModel< BasicTurbulenceModel > >, and ReynoldsStress< LESModel< BasicTurbulenceModel > >.

Definition at line 249 of file LESModel.C.

Member Data Documentation

◆ LESDict_

dictionary LESDict_
protected

LES coefficients dictionary.

Definition at line 68 of file LESModel.H.

◆ turbulence_

Switch turbulence_
protected

Turbulence on/off flag.

Definition at line 71 of file LESModel.H.

◆ printCoeffs_

Switch printCoeffs_
protected

Flag to print the model coeffs at run-time.

Definition at line 74 of file LESModel.H.

◆ coeffDict_

dictionary coeffDict_
protected

Model coefficients dictionary.

Definition at line 77 of file LESModel.H.

Referenced by LESModel< BasicTurbulenceModel >::coeffDict().

◆ Ce_

dimensionedScalar Ce_
protected

Empirical model constant.

Definition at line 80 of file LESModel.H.

Referenced by LESModel< BasicTurbulenceModel >::Ce().

◆ kMin_

dimensionedScalar kMin_
protected

Lower limit of k.

Definition at line 83 of file LESModel.H.

Referenced by LESModel< BasicTurbulenceModel >::kMin().

◆ epsilonMin_

dimensionedScalar epsilonMin_
protected

Lower limit of epsilon.

Definition at line 86 of file LESModel.H.

◆ omegaMin_

dimensionedScalar omegaMin_
protected

Lower limit for omega.

Definition at line 89 of file LESModel.H.

◆ delta_

autoPtr<Foam::LESdelta> delta_
protected

Run-time selectable delta model.

Definition at line 92 of file LESModel.H.

Referenced by LESModel< BasicTurbulenceModel >::delta().


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