LESModel< BasicMomentumTransportModel > Class Template Reference

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

Inheritance diagram for LESModel< BasicMomentumTransportModel >:
Collaboration diagram for LESModel< BasicMomentumTransportModel >:

Public Types

typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::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),(alpha, rho, U, alphaRhoPhi, phi, transport))
 
 LESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
 LESModel (const LESModel &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~LESModel ()
 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 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 void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
void operator= (const LESModel &)=delete
 Disallow default bitwise assignment. 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)
 Return a reference to the selected LES model. More...
 

Protected Member Functions

virtual void printCoeffs (const word &type)
 Print model coefficients. 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 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 BasicMomentumTransportModel>
class Foam::LESModel< BasicMomentumTransportModel >

Templated abstract base class for LES SGS models.

Source files

Definition at line 56 of file LESModel.H.

Member Typedef Documentation

◆ alphaField

typedef BasicMomentumTransportModel::alphaField alphaField

Definition at line 98 of file LESModel.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 99 of file LESModel.H.

◆ transportModel

typedef BasicMomentumTransportModel::transportModel transportModel

Definition at line 100 of file LESModel.H.

Constructor & Destructor Documentation

◆ LESModel() [1/2]

LESModel ( const word type,
const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport 
)

Construct from components.

Definition at line 44 of file LESModel.C.

References LESModel< BasicMomentumTransportModel >::New().

Referenced by LESModel< BasicMomentumTransportModel >::printCoeffs().

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

◆ LESModel() [2/2]

LESModel ( const LESModel< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~LESModel()

virtual ~LESModel ( )
inlinevirtual

Destructor.

Definition at line 159 of file LESModel.H.

References LESModel< BasicMomentumTransportModel >::read().

Here is the call graph for this function:

Member Function Documentation

◆ printCoeffs()

void printCoeffs ( const word type)
protectedvirtual

Print model coefficients.

Definition at line 31 of file LESModel.C.

References Foam::endl(), Foam::Info, and LESModel< BasicMomentumTransportModel >::LESModel().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "LES"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
LESModel< BasicMomentumTransportModel >  ,
dictionary  ,
(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)  ,
(alpha, rho, U, alphaRhoPhi, phi, transport)   
)

◆ New()

Foam::autoPtr< Foam::LESModel< BasicMomentumTransportModel > > New ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport 
)
static

Return a reference to the selected LES model.

Definition at line 124 of file LESModel.C.

References alpha(), IOobject::db(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, IOobject::group(), Foam::Info, dictionary::lookupBackwardsCompatible(), Foam::nl, phi, rho, dictionary::subDict(), and U.

Referenced by LESModel< BasicMomentumTransportModel >::LESModel().

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

◆ read()

◆ coeffDict()

virtual const dictionary& coeffDict ( ) const
inlinevirtual

Const access to the coefficients dictionary.

Definition at line 172 of file LESModel.H.

References LESModel< BasicMomentumTransportModel >::coeffDict_.

◆ kMin() [1/2]

const dimensionedScalar& kMin ( ) const
inline

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

Definition at line 178 of file LESModel.H.

References LESModel< BasicMomentumTransportModel >::kMin_.

◆ kMin() [2/2]

dimensionedScalar& kMin ( )
inline

Allow kMin to be changed.

Definition at line 184 of file LESModel.H.

References LESModel< BasicMomentumTransportModel >::kMin_.

◆ delta()

const volScalarField& delta ( ) const
inline

Access function to filter width.

Definition at line 190 of file LESModel.H.

References LESModel< BasicMomentumTransportModel >::delta_.

Referenced by dynamicLagrangian< BasicMomentumTransportModel >::k().

Here is the caller graph for this function:

◆ nuEff() [1/2]

virtual tmp<volScalarField> nuEff ( ) const
inlinevirtual

Return the effective viscosity.

Definition at line 197 of file LESModel.H.

References IOobject::groupName(), GeometricField< scalar, fvPatchField, volMesh >::New(), and nut.

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 207 of file LESModel.H.

References LESModel< BasicMomentumTransportModel >::correct(), nut, and LESModel< BasicMomentumTransportModel >::operator=().

Here is the call graph for this function:

◆ correct()

◆ operator=()

void operator= ( const LESModel< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise assignment.

Referenced by LESModel< BasicMomentumTransportModel >::nuEff().

Here is the caller graph for this function:

Member Data Documentation

◆ LESDict_

dictionary LESDict_
protected

LES coefficients dictionary.

Definition at line 66 of file LESModel.H.

◆ turbulence_

Switch turbulence_
protected

Turbulence on/off flag.

Definition at line 69 of file LESModel.H.

◆ printCoeffs_

Switch printCoeffs_
protected

Flag to print the model coeffs at run-time.

Definition at line 72 of file LESModel.H.

◆ coeffDict_

dictionary coeffDict_
protected

Model coefficients dictionary.

Definition at line 75 of file LESModel.H.

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

◆ kMin_

dimensionedScalar kMin_
protected

Lower limit of k.

Definition at line 78 of file LESModel.H.

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

◆ epsilonMin_

dimensionedScalar epsilonMin_
protected

Lower limit of epsilon.

Definition at line 81 of file LESModel.H.

◆ omegaMin_

dimensionedScalar omegaMin_
protected

Lower limit for omega.

Definition at line 84 of file LESModel.H.

◆ delta_

autoPtr<Foam::LESdelta> delta_
protected

Run-time selectable delta model.

Definition at line 87 of file LESModel.H.

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


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