LienLeschziner Class Reference

Lien and Leschziner low-Reynolds number k-epsilon turbulence model for incompressible flows. More...

Inheritance diagram for LienLeschziner:
[legend]
Collaboration diagram for LienLeschziner:
[legend]

Public Member Functions

 TypeName ("LienLeschziner")
 
 LienLeschziner (const geometricOneField &alpha, const geometricOneField &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 ~LienLeschziner ()=default
 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...
 
- Public Member Functions inherited from eddyViscosity< incompressible::RASModel >
 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 ()=default
 Destructor. More...
 
virtual bool read ()=0
 Re-read model coefficients if they have changed. 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< volScalarFieldk () const=0
 Return the turbulence kinetic energy. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
- Public Member Functions inherited from linearViscousStress< 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 ()=default
 Destructor. More...
 
virtual bool read ()=0
 Re-read model coefficients if they have changed. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field. 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...
 
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity. More...
 

Protected Member Functions

tmp< volScalarFieldfMu () const
 
tmp< volScalarFieldf2 () const
 
tmp< volScalarFieldE (const volScalarField &f2) const
 
virtual void correctNut ()
 
virtual void correctNut ()=0
 

Protected Attributes

dimensionedScalar Ceps1_
 
dimensionedScalar Ceps2_
 
dimensionedScalar sigmak_
 
dimensionedScalar sigmaEps_
 
dimensionedScalar Cmu_
 
dimensionedScalar kappa_
 
dimensionedScalar Anu_
 
dimensionedScalar Aeps_
 
dimensionedScalar AE_
 
volScalarField k_
 
volScalarField epsilon_
 
const volScalarFieldy_
 Wall distance. More...
 
- Protected Attributes inherited from eddyViscosity< incompressible::RASModel >
volScalarField nut_
 

Additional Inherited Members

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

Detailed Description

Lien and Leschziner low-Reynolds number k-epsilon turbulence model for incompressible flows.

This turbulence model is described in:

    Lien, F. S., & Leschziner, M. A. (1993).
    A pressure-velocity solution strategy for compressible flow
    and its application to shock/boundary-layer interaction
    using second-moment turbulence closure.
    Journal of fluids engineering, 115(4), 717-725.

Implemented according to the specification in: Apsley: Turbulence Models 2002

In addition to the low-Reynolds number damping functions support for wall-functions is also included to allow for low- and high-Reynolds number operation.

Source files

Definition at line 78 of file LienLeschziner.H.

Constructor & Destructor Documentation

◆ LienLeschziner()

LienLeschziner ( const geometricOneField alpha,
const geometricOneField 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 91 of file LienLeschziner.C.

References Foam::bound(), LienLeschziner::epsilon_, LienLeschziner::k_, and Foam::type().

Here is the call graph for this function:

◆ ~LienLeschziner()

virtual ~LienLeschziner ( )
virtualdefault

Destructor.

Member Function Documentation

◆ fMu()

tmp< volScalarField > fMu ( ) const
protected

Definition at line 50 of file LienLeschziner.C.

References LienLeschziner::Aeps_, LienLeschziner::Anu_, Foam::exp(), LienLeschziner::k_, nu, Foam::sqrt(), and LienLeschziner::y_.

Referenced by LienLeschziner::correctNut().

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

◆ f2()

tmp< volScalarField > f2 ( ) const
protected

Definition at line 60 of file LienLeschziner.C.

References LienLeschziner::epsilon_, Foam::exp(), LienLeschziner::k_, nu, and Foam::sqr().

Referenced by LienLeschziner::correct(), and LienLeschziner::E().

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

◆ E()

tmp< volScalarField > E ( const volScalarField f2) const
protected

Definition at line 68 of file LienLeschziner.C.

References LienLeschziner::AE_, LienLeschziner::Aeps_, LienLeschziner::Ceps2_, LienLeschziner::Cmu_, LienLeschziner::epsilon_, Foam::exp(), LienLeschziner::f2(), LienLeschziner::k_, LienLeschziner::kappa_, nu, Foam::pow(), Foam::sqr(), Foam::sqrt(), and LienLeschziner::y_.

Referenced by LienLeschziner::correct().

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

◆ correctNut()

void correctNut ( )
protectedvirtual

Implements eddyViscosity< incompressible::RASModel >.

Definition at line 82 of file LienLeschziner.C.

References LienLeschziner::Cmu_, GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), LienLeschziner::epsilon_, LienLeschziner::fMu(), LienLeschziner::k_, eddyViscosity< incompressible::RASModel >::nut_, and Foam::sqr().

Referenced by LienLeschziner::correct().

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

◆ TypeName()

TypeName ( "LienLeschziner"  )

◆ read()

bool read ( )
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< incompressible::RASModel >.

Definition at line 237 of file LienLeschziner.C.

References LienLeschziner::AE_, LienLeschziner::Aeps_, LienLeschziner::Anu_, LienLeschziner::Ceps1_, LienLeschziner::Ceps2_, LienLeschziner::Cmu_, LienLeschziner::kappa_, dimensioned< Type >::readIfPresent(), LienLeschziner::sigmaEps_, and LienLeschziner::sigmak_.

Here is the call graph for this function:

◆ DkEff()

tmp< volScalarField > DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 152 of file LienLeschziner.H.

References nu, eddyViscosity< incompressible::RASModel >::nut_, and LienLeschziner::sigmak_.

Referenced by LienLeschziner::correct().

Here is the caller graph for this function:

◆ DepsilonEff()

tmp< volScalarField > DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 161 of file LienLeschziner.H.

References nu, eddyViscosity< incompressible::RASModel >::nut_, and LienLeschziner::sigmaEps_.

Referenced by LienLeschziner::correct().

Here is the caller graph for this function:

◆ k()

virtual tmp< volScalarField > k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< incompressible::RASModel >.

Definition at line 170 of file LienLeschziner.H.

References LienLeschziner::k_.

◆ epsilon()

virtual tmp< volScalarField > epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 176 of file LienLeschziner.H.

References LienLeschziner::epsilon_.

◆ correct()

Member Data Documentation

◆ Ceps1_

dimensionedScalar Ceps1_
protected

Definition at line 89 of file LienLeschziner.H.

Referenced by LienLeschziner::correct(), and LienLeschziner::read().

◆ Ceps2_

dimensionedScalar Ceps2_
protected

◆ sigmak_

dimensionedScalar sigmak_
protected

Definition at line 91 of file LienLeschziner.H.

Referenced by LienLeschziner::DkEff(), and LienLeschziner::read().

◆ sigmaEps_

dimensionedScalar sigmaEps_
protected

Definition at line 92 of file LienLeschziner.H.

Referenced by LienLeschziner::DepsilonEff(), and LienLeschziner::read().

◆ Cmu_

dimensionedScalar Cmu_
protected

◆ kappa_

dimensionedScalar kappa_
protected

Definition at line 95 of file LienLeschziner.H.

Referenced by LienLeschziner::E(), and LienLeschziner::read().

◆ Anu_

dimensionedScalar Anu_
protected

Definition at line 97 of file LienLeschziner.H.

Referenced by LienLeschziner::fMu(), and LienLeschziner::read().

◆ Aeps_

dimensionedScalar Aeps_
protected

Definition at line 98 of file LienLeschziner.H.

Referenced by LienLeschziner::E(), LienLeschziner::fMu(), and LienLeschziner::read().

◆ AE_

dimensionedScalar AE_
protected

Definition at line 99 of file LienLeschziner.H.

Referenced by LienLeschziner::E(), and LienLeschziner::read().

◆ k_

◆ epsilon_

◆ y_

const volScalarField& y_
protected

Wall distance.

Note: different to wall distance in parent RASModel which is for near-wall cells only

Definition at line 110 of file LienLeschziner.H.

Referenced by LienLeschziner::E(), and LienLeschziner::fMu().


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