Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
LienCubicKE Class Reference

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

Inheritance diagram for LienCubicKE:
Inheritance graph
[legend]
Collaboration diagram for LienCubicKE:
Collaboration graph
[legend]

Public Member Functions

 TypeName ("LienCubicKE")
 Runtime type information. More...
 
 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. More...
 
virtual ~LienCubicKE ()
 Destructor. More...
 
virtual bool read ()
 Read RASProperties dictionary. 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 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. More...
 
virtual ~nonlinearEddyViscosity ()
 Destructor. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. 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 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 ()
 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 void validate ()
 Validate the turbulence fields after construction. More...
 
- 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. More...
 
virtual ~linearViscousStress ()
 Destructor. More...
 

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. More...
 
- 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 77 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 
)

Construct from components.

Definition at line 123 of file LienCubicKE.C.

References Foam::bound(), LienCubicKE::epsilon_, and LienCubicKE::k_.

Referenced by LienCubicKE::correctNonlinearStress().

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

◆ ~LienCubicKE()

virtual ~LienCubicKE ( )
inlinevirtual

Destructor.

Definition at line 152 of file LienCubicKE.H.

References LienCubicKE::read().

Here is the call graph for this function:

Member Function Documentation

◆ fMu()

tmp< volScalarField > fMu ( ) const
protected

Definition at line 47 of file LienCubicKE.C.

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

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

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

Referenced by LienCubicKE::correct().

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

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

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

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

Here is the call graph for this function:

◆ correctNonlinearStress()

void correctNonlinearStress ( const volTensorField gradU)
protectedvirtual

◆ TypeName()

TypeName ( "LienCubicKE"  )

Runtime type information.

◆ read()

bool read ( )
virtual

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 162 of file LienCubicKE.H.

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

Referenced by LienCubicKE::correct().

Here is the caller graph for this function:

◆ DepsilonEff()

tmp<volScalarField> DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 171 of file LienCubicKE.H.

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

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

References LienCubicKE::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 186 of file LienCubicKE.H.

References LienCubicKE::correct(), and LienCubicKE::epsilon_.

Here is the call graph for this function:

◆ correct()

void correct ( )
virtual

Member Data Documentation

◆ Ceps1_

dimensionedScalar Ceps1_
protected

Definition at line 88 of file LienCubicKE.H.

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

◆ Ceps2_

dimensionedScalar Ceps2_
protected

Definition at line 89 of file LienCubicKE.H.

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

◆ sigmak_

dimensionedScalar sigmak_
protected

Definition at line 90 of file LienCubicKE.H.

Referenced by LienCubicKE::read().

◆ sigmaEps_

dimensionedScalar sigmaEps_
protected

Definition at line 91 of file LienCubicKE.H.

Referenced by LienCubicKE::read().

◆ Cmu1_

dimensionedScalar Cmu1_
protected

Definition at line 92 of file LienCubicKE.H.

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

◆ Cmu2_

dimensionedScalar Cmu2_
protected

Definition at line 93 of file LienCubicKE.H.

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

◆ Cbeta_

dimensionedScalar Cbeta_
protected

Definition at line 94 of file LienCubicKE.H.

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

◆ Cbeta1_

dimensionedScalar Cbeta1_
protected

Definition at line 95 of file LienCubicKE.H.

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

◆ Cbeta2_

dimensionedScalar Cbeta2_
protected

Definition at line 96 of file LienCubicKE.H.

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

◆ Cbeta3_

dimensionedScalar Cbeta3_
protected

Definition at line 97 of file LienCubicKE.H.

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

◆ Cgamma1_

dimensionedScalar Cgamma1_
protected

Definition at line 98 of file LienCubicKE.H.

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

◆ Cgamma2_

dimensionedScalar Cgamma2_
protected

Definition at line 99 of file LienCubicKE.H.

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

◆ Cgamma4_

dimensionedScalar Cgamma4_
protected

Definition at line 100 of file LienCubicKE.H.

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

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 102 of file LienCubicKE.H.

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

◆ kappa_

dimensionedScalar kappa_
protected

Definition at line 103 of file LienCubicKE.H.

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

◆ Anu_

dimensionedScalar Anu_
protected

Definition at line 105 of file LienCubicKE.H.

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

◆ AE_

dimensionedScalar AE_
protected

Definition at line 106 of file LienCubicKE.H.

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

◆ k_

volScalarField k_
protected

◆ epsilon_

volScalarField epsilon_
protected

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

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


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