LaunderSharmaKE< BasicTurbulenceModel > Class Template Reference

Launder and Sharma low-Reynolds k-epsilon turbulence model for incompressible and compressible and combusting flows including rapid distortion theory (RDT) based compression term. More...

Inheritance diagram for LaunderSharmaKE< BasicTurbulenceModel >:
Collaboration diagram for LaunderSharmaKE< BasicTurbulenceModel >:

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
typedef RASModel< BasicTurbulenceModel > ::alphaField alphaField
 
typedef RASModel< BasicTurbulenceModel > ::rhoField rhoField
 
typedef RASModel< BasicTurbulenceModel > ::transportModel transportModel
 
- Public Types inherited from linearViscousStress< RASModel< BasicTurbulenceModel > >
typedef RASModel< BasicTurbulenceModel > ::alphaField alphaField
 
typedef RASModel< BasicTurbulenceModel > ::rhoField rhoField
 
typedef RASModel< BasicTurbulenceModel > ::transportModel transportModel
 
- Public Types inherited from RASModel< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 

Public Member Functions

 TypeName ("LaunderSharmaKE")
 Runtime type information. More...
 
 LaunderSharmaKE (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...
 
 LaunderSharmaKE (const LaunderSharmaKE &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~LaunderSharmaKE ()
 Destructor. More...
 
virtual bool read ()
 Re-read model coefficients if they have changed. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
void operator= (const LaunderSharmaKE &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from eddyViscosity< RASModel< 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< RASModel< 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 RASModel< BasicTurbulenceModel >
 TypeName ("RAS")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, RASModel, 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))
 
 RASModel (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...
 
 RASModel (const RASModel &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~RASModel ()
 Destructor. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: small) More...
 
const dimensionedScalarepsilonMin () const
 Return the lower allowable limit for epsilon (default: small) More...
 
const dimensionedScalaromegaMin () const
 Return the lower allowable limit for omega (default: small) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
dimensionedScalarepsilonMin ()
 Allow epsilonMin to be changed. More...
 
dimensionedScalaromegaMin ()
 Allow omegaMin to be changed. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 
void operator= (const RASModel &)=delete
 Disallow default bitwise assignment. More...
 

Protected Member Functions

tmp< volScalarFieldfMu () const
 
tmp< volScalarFieldf2 () const
 
virtual void correctNut ()
 
virtual tmp< fvScalarMatrixkSource () const
 
virtual tmp< fvScalarMatrixepsilonSource () const
 
- Protected Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 

Protected Attributes

dimensionedScalar Cmu_
 
dimensionedScalar C1_
 
dimensionedScalar C2_
 
dimensionedScalar C3_
 
dimensionedScalar sigmak_
 
dimensionedScalar sigmaEps_
 
volScalarField k_
 
volScalarField epsilon_
 
- Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_
 
- Protected Attributes inherited from RASModel< BasicTurbulenceModel >
dictionary RASDict_
 RAS 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...
 

Additional Inherited Members

- Static Public Member Functions inherited from RASModel< BasicTurbulenceModel >
static autoPtr< RASModelNew (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 RAS model. More...
 

Detailed Description

template<class BasicTurbulenceModel>
class Foam::RASModels::LaunderSharmaKE< BasicTurbulenceModel >

Launder and Sharma low-Reynolds k-epsilon turbulence model for incompressible and compressible and combusting flows including rapid distortion theory (RDT) based compression term.

References:

    Launder, B. E., & Sharma, B. I. (1974).
    Application of the energy-dissipation model of turbulence to the
    calculation of flow near a spinning disc.
    Letters in heat and mass transfer, 1(2), 131-137.

    For the RDT-based compression term:
    El Tahry, S. H. (1983).
    k-epsilon equation for compressible reciprocating engine flows.
    Journal of Energy, 7(4), 345-353.

The default model coefficients are

    LaunderSharmaKECoeffs
    {
        Cmu         0.09;
        C1          1.44;
        C2          1.92;
        C3          -0.33;
        alphah      1.0;    // only for compressible
        alphahk     1.0;    // only for compressible
        alphaEps    0.76923;
    }
Source files

Definition at line 81 of file LaunderSharmaKE.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 117 of file LaunderSharmaKE.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 118 of file LaunderSharmaKE.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 119 of file LaunderSharmaKE.H.

Constructor & Destructor Documentation

◆ LaunderSharmaKE() [1/2]

LaunderSharmaKE ( 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 100 of file LaunderSharmaKE.C.

References Foam::bound().

Referenced by LaunderSharmaKE< BasicTurbulenceModel >::epsilonSource().

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

◆ LaunderSharmaKE() [2/2]

LaunderSharmaKE ( const LaunderSharmaKE< BasicTurbulenceModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~LaunderSharmaKE()

virtual ~LaunderSharmaKE ( )
inlinevirtual

Destructor.

Definition at line 146 of file LaunderSharmaKE.H.

References LaunderSharmaKE< BasicTurbulenceModel >::read().

Here is the call graph for this function:

Member Function Documentation

◆ fMu()

tmp< volScalarField > fMu ( ) const
protected

Definition at line 40 of file LaunderSharmaKE.C.

References Foam::exp(), nu, and Foam::sqr().

Here is the call graph for this function:

◆ f2()

tmp< volScalarField > f2 ( ) const
protected

Definition at line 47 of file LaunderSharmaKE.C.

References Foam::exp(), Foam::min(), nu, and Foam::sqr().

Here is the call graph for this function:

◆ correctNut()

void correctNut ( )
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 56 of file LaunderSharmaKE.C.

References optionList::correct(), options::New(), and Foam::sqr().

Here is the call graph for this function:

◆ kSource()

tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Definition at line 67 of file LaunderSharmaKE.C.

References Foam::dimTime, and Foam::dimVolume.

◆ epsilonSource()

tmp< fvScalarMatrix > epsilonSource ( ) const
protectedvirtual

Definition at line 82 of file LaunderSharmaKE.C.

References Foam::dimTime, Foam::dimVolume, and LaunderSharmaKE< BasicTurbulenceModel >::LaunderSharmaKE().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "LaunderSharmaKE< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 217 of file LaunderSharmaKE.C.

References Foam::read().

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

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

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 156 of file LaunderSharmaKE.H.

References GeometricField< scalar, fvPatchField, volMesh >::New(), nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.

Here is the call graph for this function:

◆ DepsilonEff()

tmp<volScalarField> DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 166 of file LaunderSharmaKE.H.

References GeometricField< scalar, fvPatchField, volMesh >::New(), nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.

Here is the call graph for this function:

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 176 of file LaunderSharmaKE.H.

References LaunderSharmaKE< BasicTurbulenceModel >::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 182 of file LaunderSharmaKE.H.

References LaunderSharmaKE< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::epsilon_, and LaunderSharmaKE< BasicTurbulenceModel >::operator=().

Here is the call graph for this function:

◆ correct()

◆ operator=()

void operator= ( const LaunderSharmaKE< BasicTurbulenceModel > &  )
delete

Disallow default bitwise assignment.

Referenced by LaunderSharmaKE< BasicTurbulenceModel >::epsilon().

Here is the caller graph for this function:

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 91 of file LaunderSharmaKE.H.

◆ C1_

dimensionedScalar C1_
protected

Definition at line 92 of file LaunderSharmaKE.H.

◆ C2_

dimensionedScalar C2_
protected

Definition at line 93 of file LaunderSharmaKE.H.

◆ C3_

dimensionedScalar C3_
protected

Definition at line 94 of file LaunderSharmaKE.H.

◆ sigmak_

dimensionedScalar sigmak_
protected

Definition at line 95 of file LaunderSharmaKE.H.

◆ sigmaEps_

dimensionedScalar sigmaEps_
protected

Definition at line 96 of file LaunderSharmaKE.H.

◆ k_

volScalarField k_
protected

Definition at line 101 of file LaunderSharmaKE.H.

Referenced by LaunderSharmaKE< BasicTurbulenceModel >::k().

◆ epsilon_

volScalarField epsilon_
protected

Definition at line 102 of file LaunderSharmaKE.H.

Referenced by LaunderSharmaKE< BasicTurbulenceModel >::epsilon().


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