Loading...
Searching...
No Matches
LienCubicKE Class Reference

Lien cubic non-linear low-Reynolds k-epsilon turbulence models for incompressible flows. More...

#include <LienCubicKE.H>

Inheritance diagram for LienCubicKE:
Collaboration diagram for LienCubicKE:

Public Member Functions

 TypeName ("LienCubicKE")
 Runtime type information.
 LienCubicKE (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.
virtual ~LienCubicKE ()=default
 Destructor.
virtual bool read ()
 Re-read model coefficients if they have changed.
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k.
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon.
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy.
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate.
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity.
Public Member Functions inherited from nonlinearEddyViscosity< incompressible::RASModel >
 nonlinearEddyViscosity (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.
virtual ~nonlinearEddyViscosity ()=default
 Destructor.
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor.
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation.
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.
virtual ~eddyViscosity ()=default
 Destructor.
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity.
virtual void validate ()
 Validate the turbulence fields after construction.
Public Member Functions inherited from linearViscousStress< incompressible::RASModel >
 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.
virtual ~linearViscousStress ()=default
 Destructor.

Protected Member Functions

tmp< volScalarFieldfMu () const
tmp< volScalarFieldf2 () const
tmp< volScalarFieldE (const volScalarField &f2) const
virtual void correctNut ()
virtual void correctNonlinearStress (const volTensorField &gradU)

Protected Attributes

dimensionedScalar Ceps1_
dimensionedScalar Ceps2_
dimensionedScalar sigmak_
dimensionedScalar sigmaEps_
dimensionedScalar Cmu1_
dimensionedScalar Cmu2_
dimensionedScalar Cbeta_
dimensionedScalar Cbeta1_
dimensionedScalar Cbeta2_
dimensionedScalar Cbeta3_
dimensionedScalar Cgamma1_
dimensionedScalar Cgamma2_
dimensionedScalar Cgamma4_
dimensionedScalar Cmu_
dimensionedScalar kappa_
dimensionedScalar Anu_
dimensionedScalar AE_
volScalarField k_
volScalarField epsilon_
const volScalarFieldy_
 Wall distance.
Protected Attributes inherited from nonlinearEddyViscosity< incompressible::RASModel >
volSymmTensorField nonlinearStress_
Protected Attributes inherited from eddyViscosity< incompressible::RASModel >
volScalarField nut_

Additional Inherited Members

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

Detailed Description

Lien cubic non-linear low-Reynolds k-epsilon turbulence models for incompressible flows.

This turbulence model is described in:

    Lien, F.S., Chen, W.L. & Leschziner, M.A. (1996).
    Low-Reynolds-number eddy-viscosity modeling based on non-linear
    stress-strain/vorticity relations.
    Engineering Turbulence Modelling and Experiments 3, 91-100.

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.

See also
Foam::incompressible::RASModels::ShihQuadraticKE
Source files

Definition at line 76 of file LienCubicKE.H.

Constructor & Destructor Documentation

◆ LienCubicKE()

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

◆ ~LienCubicKE()

virtual ~LienCubicKE ( )
virtualdefault

Destructor.

Member Function Documentation

◆ fMu()

tmp< volScalarField > fMu ( ) const
protected

Definition at line 43 of file LienCubicKE.C.

References Anu_, Cmu_, Foam::exp(), k_, kappa_, nu, Foam::pow(), Foam::sqrt(), and y_.

Referenced by correctNonlinearStress().

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 53 of file LienCubicKE.C.

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

Referenced by correct(), and 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 61 of file LienCubicKE.C.

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

Referenced by 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 75 of file LienCubicKE.C.

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

Here is the call graph for this function:

◆ correctNonlinearStress()

◆ TypeName()

TypeName ( "LienCubicKE" )

◆ read()

bool read ( )
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< incompressible::RASModel >.

Definition at line 336 of file LienCubicKE.C.

References AE_, Anu_, Cbeta1_, Cbeta2_, Cbeta3_, Cbeta_, Ceps1_, Ceps2_, Cgamma1_, Cgamma2_, Cgamma4_, Cmu1_, Cmu2_, Cmu_, kappa_, eddyViscosity< BasicTurbulenceModel >::read(), sigmaEps_, and sigmak_.

Here is the call graph for this function:

◆ DkEff()

tmp< volScalarField > DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 173 of file LienCubicKE.H.

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

Referenced by correct().

Here is the caller graph for this function:

◆ DepsilonEff()

tmp< volScalarField > DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 184 of file LienCubicKE.H.

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

Referenced by 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 195 of file LienCubicKE.H.

References k_.

◆ epsilon()

virtual tmp< volScalarField > epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 203 of file LienCubicKE.H.

References epsilon_.

◆ correct()

Member Data Documentation

◆ Ceps1_

dimensionedScalar Ceps1_
protected

Definition at line 87 of file LienCubicKE.H.

Referenced by correct(), LienCubicKE(), and read().

◆ Ceps2_

dimensionedScalar Ceps2_
protected

Definition at line 88 of file LienCubicKE.H.

Referenced by correct(), E(), LienCubicKE(), and read().

◆ sigmak_

dimensionedScalar sigmak_
protected

Definition at line 89 of file LienCubicKE.H.

Referenced by DkEff(), LienCubicKE(), and read().

◆ sigmaEps_

dimensionedScalar sigmaEps_
protected

Definition at line 90 of file LienCubicKE.H.

Referenced by DepsilonEff(), LienCubicKE(), and read().

◆ Cmu1_

dimensionedScalar Cmu1_
protected

Definition at line 91 of file LienCubicKE.H.

Referenced by correctNonlinearStress(), LienCubicKE(), and read().

◆ Cmu2_

dimensionedScalar Cmu2_
protected

Definition at line 92 of file LienCubicKE.H.

Referenced by correctNonlinearStress(), LienCubicKE(), and read().

◆ Cbeta_

dimensionedScalar Cbeta_
protected

Definition at line 93 of file LienCubicKE.H.

Referenced by correctNonlinearStress(), LienCubicKE(), and read().

◆ Cbeta1_

dimensionedScalar Cbeta1_
protected

Definition at line 94 of file LienCubicKE.H.

Referenced by correctNonlinearStress(), LienCubicKE(), and read().

◆ Cbeta2_

dimensionedScalar Cbeta2_
protected

Definition at line 95 of file LienCubicKE.H.

Referenced by correctNonlinearStress(), LienCubicKE(), and read().

◆ Cbeta3_

dimensionedScalar Cbeta3_
protected

Definition at line 96 of file LienCubicKE.H.

Referenced by correctNonlinearStress(), LienCubicKE(), and read().

◆ Cgamma1_

dimensionedScalar Cgamma1_
protected

Definition at line 97 of file LienCubicKE.H.

Referenced by correctNonlinearStress(), LienCubicKE(), and read().

◆ Cgamma2_

dimensionedScalar Cgamma2_
protected

Definition at line 98 of file LienCubicKE.H.

Referenced by correctNonlinearStress(), LienCubicKE(), and read().

◆ Cgamma4_

dimensionedScalar Cgamma4_
protected

Definition at line 99 of file LienCubicKE.H.

Referenced by correctNonlinearStress(), LienCubicKE(), and read().

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 101 of file LienCubicKE.H.

Referenced by E(), fMu(), LienCubicKE(), and read().

◆ kappa_

dimensionedScalar kappa_
protected

Definition at line 102 of file LienCubicKE.H.

Referenced by E(), fMu(), LienCubicKE(), and read().

◆ Anu_

dimensionedScalar Anu_
protected

Definition at line 104 of file LienCubicKE.H.

Referenced by fMu(), LienCubicKE(), and read().

◆ AE_

dimensionedScalar AE_
protected

Definition at line 105 of file LienCubicKE.H.

Referenced by E(), LienCubicKE(), and read().

◆ k_

volScalarField k_
protected

Definition at line 110 of file LienCubicKE.H.

Referenced by correct(), correctNonlinearStress(), E(), f2(), fMu(), k(), and LienCubicKE().

◆ epsilon_

volScalarField epsilon_
protected

Definition at line 111 of file LienCubicKE.H.

Referenced by correct(), correctNonlinearStress(), E(), epsilon(), f2(), and LienCubicKE().

◆ 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 119 of file LienCubicKE.H.

Referenced by E(), fMu(), and LienCubicKE().


The documentation for this class was generated from the following files:
  • src/TurbulenceModels/incompressible/turbulentTransportModels/RAS/LienCubicKE/LienCubicKE.H
  • src/TurbulenceModels/incompressible/turbulentTransportModels/RAS/LienCubicKE/LienCubicKE.C