OpenFOAM logo
The Open Source CFD Toolbox
  Source Guide OpenCFD Solutions Contact OpenFOAM

LESModel Class Reference

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

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 dictionarycoeffDict () const
 Const access to the coefficients dictionary,.
const dimensionedScalark0 () const
 Return the value of k0 which k is not allowed to be less than.
dimensionedScalark0 ()
 Allow k0 to be changed.
const volScalarFielddelta () const
 Access function to filter width.
virtual tmp< volScalarFieldk () const =0
 Return the SGS turbulent kinetic energy.
virtual tmp< volScalarFieldepsilon () const =0
 Return the SGS turbulent dissipation.
virtual tmp< volScalarFieldmuSgs () const =0
 Return the SGS turbulent viscosity.
virtual tmp< volScalarFieldmuEff () const
 Return the effective viscosity.
virtual tmp< volScalarFieldalphaSgs () const =0
 Return the SGS turbulent thermal diffusivity.
virtual tmp< volScalarFieldalphaEff () const =0
 Return the SGS thermal conductivity.
virtual tmp< volSymmTensorFieldB () const =0
 Return the sub-grid stress tensor.
virtual tmp< volSymmTensorFielddevRhoBeff () const =0
 Return the deviatoric part of the effective sub-grid.
virtual tmp< fvVectorMatrixdivDevRhoBeff (volVectorField &U) const =0
 Returns div(rho*dev(B)).
virtual tmp< volScalarFieldmut () const
 Return the turbulence viscosity.
virtual tmp< volScalarFieldalphat () const
 Return the turbulence thermal diffusivity.
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor including the laminar stress.
virtual tmp< fvVectorMatrixdivDevRhoReff (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< LESModelNew (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< LESdeltadelta_

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


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 48 of file LESModel.C.

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

Here is the call graph for this function:

virtual ~LESModel (  )  [inline, virtual]

Destructor.

Definition at line 151 of file LESModel.H.


Member Function Documentation

void printCoeffs (  )  [protected, virtual]

Print model coefficients.

Definition at line 36 of file LESModel.C.

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]

Return a reference to the selected LES model.

Definition at line 88 of file LESModel.C.

References turbulenceModel::mesh_.

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

References LESModel::coeffDict().

Referenced by LESModel::coeffDict().

Here is the call graph for this function:

Here is the caller graph for this function:

const dimensionedScalar& k0 (  )  const [inline]

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

Definition at line 170 of file LESModel.H.

References LESModel::coeffDict_.

Referenced by oneEqEddy::read(), and lowReOneEqEddy::read().

Here is the caller graph for this function:

dimensionedScalar& k0 (  )  [inline]

Allow k0 to be changed.

Definition at line 176 of file LESModel.H.

References LESModel::k0_.

const volScalarField& delta (  )  const [inline]

Access function to filter width.

Definition at line 182 of file LESModel.H.

References LESModel::k0_.

Referenced by GenSGSStress::muSgs(), GenEddyVisc::muSgs(), dynOneEqEddy::read(), and lowReOneEqEddy::read().

Here is the caller graph for this function:

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]

Return the SGS turbulent viscosity.

Implemented in GenEddyVisc, GenSGSStress, and SpalartAllmaras.

Referenced by LESModel::muEff(), and LESModel::mut().

Here is the caller graph for this function:

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

Return the effective viscosity.

Implements turbulenceModel.

Definition at line 198 of file LESModel.H.

References mu, and LESModel::muSgs().

Here is the call graph for this function:

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

Return the SGS turbulent thermal diffusivity.

Implemented in GenEddyVisc, GenSGSStress, and SpalartAllmaras.

Referenced by LESModel::alphat().

Here is the caller graph for this function:

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

Here is the caller graph for this function:

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

Here is the caller graph for this function:

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

Here is the caller graph for this function:

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

Return the turbulence viscosity.

Implements turbulenceModel.

Definition at line 231 of file LESModel.H.

References LESModel::muSgs(), and LESModel::mut().

Referenced by LESModel::mut().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Return the turbulence thermal diffusivity.

Definition at line 237 of file LESModel.H.

References LESModel::alphaSgs(), and LESModel::alphat().

Referenced by LESModel::alphat().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Return the Reynolds stress tensor.

Implements turbulenceModel.

Definition at line 243 of file LESModel.H.

References LESModel::B(), and LESModel::R().

Referenced by LESModel::R().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Return the effective stress tensor including the laminar stress.

Implements turbulenceModel.

Definition at line 249 of file LESModel.H.

References LESModel::devRhoBeff(), and LESModel::devRhoReff().

Referenced by LESModel::devRhoReff().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Return the source term for the momentum equation.

Definition at line 255 of file LESModel.H.

References LESModel::divDevRhoBeff(), and LESModel::divDevRhoReff().

Referenced by LESModel::divDevRhoReff().

Here is the call graph for this function:

Here is the caller graph for this function:

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 146 of file LESModel.C.

References LESModel::delta_.

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 140 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 152 of file LESModel.C.


Member Data Documentation

Switch printCoeffs_ [protected]

Definition at line 77 of file LESModel.H.

dictionary coeffDict_ [protected]

Definition at line 78 of file LESModel.H.

Referenced by LESModel::k0(), and LESModel::LESModel().

dimensionedScalar k0_ [protected]

Definition at line 80 of file LESModel.H.

Referenced by LESModel::delta(), and LESModel::k0().

autoPtr<LESdelta> delta_ [protected]

Definition at line 82 of file LESModel.H.

Referenced by LESModel::correct().


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
Copyright © 2000-2009 OpenCFD Ltd