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

Dynamic SGS model with Lagrangian averaging. More...

Inheritance diagram for dynamicLagrangian< BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for dynamicLagrangian< 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 ("dynamicLagrangian")
 Runtime type information. More...
 
 dynamicLagrangian (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 ~dynamicLagrangian ()
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
tmp< volScalarFieldk (const tmp< volTensorField > &gradU) const
 Return SGS kinetic energy. More...
 
virtual tmp< volScalarFieldk () const
 Return SGS kinetic energy. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. 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...
 
virtual tmp< volScalarFieldepsilon () const
 Return sub-grid disipation rate. 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

void correctNut (const tmp< volTensorField > &gradU)
 Update sub-grid eddy-viscosity. More...
 
virtual void correctNut ()
 
- Protected Member Functions inherited from LESModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 

Protected Attributes

volScalarField flm_
 
volScalarField fmm_
 
dimensionedScalar theta_
 
simpleFilter simpleFilter_
 
autoPtr< LESfilterfilterPtr_
 
LESfilterfilter_
 
dimensionedScalar flm0_
 
dimensionedScalar fmm0_
 
- 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::dynamicLagrangian< BasicTurbulenceModel >

Dynamic SGS model with Lagrangian averaging.

Reference:

    Meneveau, C., Lund, T. S., & Cabot, W. H. (1996).
    A Lagrangian dynamic subgrid-scale model of turbulence.
    Journal of Fluid Mechanics, 319, 353-385.
Source files

Definition at line 59 of file dynamicLagrangian.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 97 of file dynamicLagrangian.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 98 of file dynamicLagrangian.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 99 of file dynamicLagrangian.H.

Constructor & Destructor Documentation

◆ dynamicLagrangian()

dynamicLagrangian ( 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 63 of file dynamicLagrangian.C.

◆ ~dynamicLagrangian()

virtual ~dynamicLagrangian ( )
inlinevirtual

Destructor.

Definition at line 122 of file dynamicLagrangian.H.

References dynamicLagrangian< BasicTurbulenceModel >::read().

Here is the call graph for this function:

Member Function Documentation

◆ correctNut() [1/2]

void correctNut ( const tmp< volTensorField > &  gradU)
protected

Update sub-grid eddy-viscosity.

Definition at line 40 of file dynamicLagrangian.C.

References optionList::correct(), delta, Foam::dev(), Foam::mag(), options::New(), Foam::sqr(), and Foam::symm().

Here is the call graph for this function:

◆ correctNut() [2/2]

void correctNut ( )
protectedvirtual

Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.

Definition at line 53 of file dynamicLagrangian.C.

References Foam::fvc::grad().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "dynamicLagrangian< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Read model coefficients if they have changed.

Reimplemented from LESeddyViscosity< BasicTurbulenceModel >.

Definition at line 137 of file dynamicLagrangian.C.

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

Here is the caller graph for this function:

◆ k() [1/2]

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

Return SGS kinetic energy.

Definition at line 132 of file dynamicLagrangian.H.

References LESeddyViscosity< BasicTurbulenceModel >::Ce_, LESModel< BasicTurbulenceModel >::delta(), Foam::dev(), Foam::magSqr(), Foam::pow(), Foam::sqr(), and Foam::symm().

Here is the call graph for this function:

◆ k() [2/2]

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return SGS kinetic energy.

Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.

Definition at line 141 of file dynamicLagrangian.H.

References Foam::fvc::grad().

Here is the call graph for this function:

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 147 of file dynamicLagrangian.H.

References dynamicLagrangian< BasicTurbulenceModel >::correct(), nu, and eddyViscosity< LESModel< BasicTurbulenceModel > >::nut_.

Here is the call graph for this function:

◆ correct()

void correct ( )
virtual

Member Data Documentation

◆ flm_

volScalarField flm_
protected

Definition at line 74 of file dynamicLagrangian.H.

◆ fmm_

volScalarField fmm_
protected

Definition at line 75 of file dynamicLagrangian.H.

◆ theta_

dimensionedScalar theta_
protected

Definition at line 77 of file dynamicLagrangian.H.

◆ simpleFilter_

simpleFilter simpleFilter_
protected

Definition at line 79 of file dynamicLagrangian.H.

◆ filterPtr_

autoPtr<LESfilter> filterPtr_
protected

Definition at line 80 of file dynamicLagrangian.H.

◆ filter_

LESfilter& filter_
protected

Definition at line 81 of file dynamicLagrangian.H.

◆ flm0_

dimensionedScalar flm0_
protected

Definition at line 83 of file dynamicLagrangian.H.

◆ fmm0_

dimensionedScalar fmm0_
protected

Definition at line 84 of file dynamicLagrangian.H.


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