FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

Base class for all compressible flow LES SGS models. More...

#include <compressibleLESModels/LESModel.H>


Detailed Description

Base class for all compressible flow LES SGS models.

This class defines the basic interface for a compressible flow SGS model, and encapsulates data of value to all possible models. In particular this includes references to all the dependent fields (rho, U, phi), the physical viscosity mu, and the LESProperties dictionary, which contains the model selection and model coefficients.

Source files

Definition at line 71 of file LESModel.H.

Inheritance diagram for LESModel:
Collaboration diagram for LESModel:

List of all members.

Public Member Functions

 TypeName ("LESModel")
 Runtime type information.
 declareRunTimeSelectionTable (autoPtr, LESModel, dictionary,(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel),(rho, U, phi, thermoPhysicalModel))
 LESModel (const word &type, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel)
 Construct from components.
virtual  ~LESModel ()
 Destructor.
const dictionary &  coeffDict () const
 Const access to the coefficients dictionary,.
const dimensionedScalar &  k0 () const
 Return the value of k0 which k is not allowed to be less than.
dimensionedScalar &  k0 ()
 Allow k0 to be changed.
const volScalarField &  delta () const
 Access function to filter width.
virtual tmp< volScalarField >  k () const =0
 Return the SGS turbulent kinetic energy.
virtual tmp< volScalarField >  epsilon () const =0
 Return the SGS turbulent dissipation.
virtual tmp< volScalarField >  muSgs () const =0
 Return the SGS turbulent viscosity.
virtual tmp< volScalarField >  muEff () const
 Return the effective viscosity.
virtual tmp< volScalarField >  alphaSgs () const =0
 Return the SGS turbulent thermal diffusivity.
virtual tmp< volScalarField >  alphaEff () const =0
 Return the SGS thermal conductivity.
virtual tmp< volSymmTensorField >  B () const =0
 Return the sub-grid stress tensor.
virtual tmp< volSymmTensorField >  devRhoBeff () const =0
 Return the deviatoric part of the effective sub-grid.
virtual tmp< fvVectorMatrix >  divDevRhoBeff (volVectorField &U) const =0
 Returns div(rho*dev(B)).
virtual tmp< volScalarField >  mut () const
 Return the turbulence viscosity.
virtual tmp< volScalarField >  alphat () const
 Return the turbulence thermal diffusivity.
virtual tmp< volSymmTensorField >  R () const
 Return the Reynolds stress tensor.
virtual tmp< volSymmTensorField >  devRhoReff () const
 Return the effective stress tensor including the laminar stress.
virtual tmp< fvVectorMatrix >  divDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation.
void  correct ()
 Correct Eddy-Viscosity and related properties.
virtual void  correct (const tmp< volTensorField > &gradU)
 Correct Eddy-Viscosity and related properties.
virtual bool  read ()=0
 Read LESProperties dictionary.

Static Public Member Functions

static autoPtr< LESModel >  New (const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel)
 Return a reference to the selected LES model.

Protected Member Functions

virtual void  printCoeffs ()
 Print model coefficients.

Protected Attributes

Switch  printCoeffs_
dictionary  coeffDict_
dimensionedScalar  k0_
autoPtr< LESdelta >  delta_

Constructor & Destructor Documentation

LESModel ( const word &   type,
const volScalarField &   rho,
const volVectorField &   U,
const surfaceScalarField &   phi,
const basicThermo &   thermoPhysicalModel  
)

Construct from components.

Definition at line 56 of file LESModel.C.

References readIfPresent().

virtual ~LESModel (  ) [inline, virtual]

Destructor.

Definition at line 155 of file LESModel.H.


Member Function Documentation

void printCoeffs (  ) [protected, virtual]

Print model coefficients.

Definition at line 44 of file LESModel.C.

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

TypeName ( "LESModel"    )

Runtime type information.

declareRunTimeSelectionTable ( autoPtr   ,
LESModel   ,
dictionary   ,
(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel)   ,
(rho, U, phi, thermoPhysicalModel)    
)
autoPtr< LESModel > New ( const volScalarField &   rho,
const volVectorField &   U,
const surfaceScalarField &   phi,
const basicThermo &   thermoPhysicalModel  
) [static]
const dictionary& coeffDict (  ) const [inline]

Const access to the coefficients dictionary,.

which provides info. about choice of models, and all related data (particularly model coefficients).

Definition at line 166 of file LESModel.H.

References LESModel::coeffDict_.

const dimensionedScalar& k0 (  ) const [inline]

Return the value of k0 which k is not allowed to be less than.

Definition at line 172 of file LESModel.H.

References LESModel::k0_.

dimensionedScalar& k0 (  ) [inline]

Allow k0 to be changed.

Definition at line 178 of file LESModel.H.

References LESModel::k0_.

const volScalarField& delta (  ) const [inline]

Access function to filter width.

Definition at line 184 of file LESModel.H.

References LESModel::delta_.

Referenced by GenSGSStress::epsilon(), and GenEddyVisc::epsilon().

virtual tmp<volScalarField> k (  ) const [pure virtual]

Return the SGS turbulent kinetic energy.

Implements turbulenceModel.

Implemented in GenEddyVisc, GenSGSStress, and SpalartAllmaras.

virtual tmp<volScalarField> epsilon (  ) const [pure virtual]

Return the SGS turbulent dissipation.

Implements turbulenceModel.

Implemented in GenEddyVisc, GenSGSStress, and SpalartAllmaras.

virtual tmp<volScalarField> muSgs (  ) const [pure virtual]
virtual tmp<volScalarField> muEff (  ) const [inline, virtual]

Return the effective viscosity.

Implements turbulenceModel.

Definition at line 200 of file LESModel.H.

References mu, and LESModel::muSgs().

virtual tmp<volScalarField> alphaSgs (  ) const [pure virtual]

Return the SGS turbulent thermal diffusivity.

Implemented in GenEddyVisc, GenSGSStress, and SpalartAllmaras.

Referenced by LESModel::alphat().

virtual tmp<volScalarField> alphaEff (  ) const [pure virtual]

Return the SGS thermal conductivity.

Implements turbulenceModel.

Implemented in GenEddyVisc, GenSGSStress, and SpalartAllmaras.

virtual tmp<volSymmTensorField> B (  ) const [pure virtual]

Return the sub-grid stress tensor.

Implemented in GenEddyVisc, GenSGSStress, and SpalartAllmaras.

Referenced by LESModel::R().

virtual tmp<volSymmTensorField> devRhoBeff (  ) const [pure virtual]

Return the deviatoric part of the effective sub-grid.

turbulence stress tensor including the laminar stress

Implemented in GenEddyVisc, GenSGSStress, and SpalartAllmaras.

Referenced by LESModel::devRhoReff().

virtual tmp<fvVectorMatrix> divDevRhoBeff ( volVectorField &   U  ) const [pure virtual]

Returns div(rho*dev(B)).

This is the additional term due to the filtering of the NSE.

Implemented in GenEddyVisc, GenSGSStress, and SpalartAllmaras.

Referenced by LESModel::divDevRhoReff().

virtual tmp<volScalarField> mut (  ) const [inline, virtual]

Return the turbulence viscosity.

Implements turbulenceModel.

Definition at line 229 of file LESModel.H.

References LESModel::muSgs().

virtual tmp<volScalarField> alphat (  ) const [inline, virtual]

Return the turbulence thermal diffusivity.

Definition at line 235 of file LESModel.H.

References LESModel::alphaSgs().

virtual tmp<volSymmTensorField> R (  ) const [inline, virtual]

Return the Reynolds stress tensor.

Implements turbulenceModel.

Definition at line 241 of file LESModel.H.

References LESModel::B().

virtual tmp<volSymmTensorField> devRhoReff (  ) const [inline, virtual]

Return the effective stress tensor including the laminar stress.

Implements turbulenceModel.

Definition at line 247 of file LESModel.H.

References LESModel::devRhoBeff().

virtual tmp<fvVectorMatrix> divDevRhoReff ( volVectorField &   U  ) const [inline, virtual]

Return the source term for the momentum equation.

Definition at line 253 of file LESModel.H.

References LESModel::divDevRhoBeff().

void correct (  ) [virtual]

Correct Eddy-Viscosity and related properties.

This calls correct(const tmp<volTensorField>& gradU) by supplying gradU calculated locally.

Implements turbulenceModel.

Definition at line 154 of file LESModel.C.

References correct(), and Foam::fvc::grad().

void correct ( const tmp< volTensorField > &   gradU  ) [virtual]

Correct Eddy-Viscosity and related properties.

Reimplemented in DeardorffDiffStress, dynOneEqEddy, GenEddyVisc, GenSGSStress, lowReOneEqEddy, oneEqEddy, Smagorinsky, and SpalartAllmaras.

Definition at line 148 of file LESModel.C.

bool read (  ) [pure virtual]

Read LESProperties dictionary.

Implements turbulenceModel.

Implemented in DeardorffDiffStress, dynOneEqEddy, GenEddyVisc, GenSGSStress, lowReOneEqEddy, oneEqEddy, Smagorinsky, and SpalartAllmaras.

Definition at line 160 of file LESModel.C.

References regIOobject::read(), readIfPresent(), and Foam::type().


Member Data Documentation

Switch printCoeffs_ [protected]

Definition at line 81 of file LESModel.H.

dictionary coeffDict_ [protected]

Definition at line 82 of file LESModel.H.

Referenced by LESModel::coeffDict().

dimensionedScalar k0_ [protected]

Definition at line 84 of file LESModel.H.

Referenced by LESModel::k0().

autoPtr<LESdelta> delta_ [protected]

Definition at line 86 of file LESModel.H.

Referenced by LESModel::delta().


The documentation for this class was generated from the following files:
  • src/turbulenceModels/compressible/LES/LESModel/LESModel.H
  • src/turbulenceModels/compressible/LES/LESModel/LESModel.C