Lien cubic non-linear low-Reynolds k-epsilon turbulence models for incompressible flows. More...
#include <LienCubicKE.H>


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< volScalarField > | DkEff () const |
| Return the effective diffusivity for k. | |
| tmp< volScalarField > | DepsilonEff () const |
| Return the effective diffusivity for epsilon. | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy. | |
| virtual tmp< volScalarField > | epsilon () 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< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (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< volScalarField > | nut () 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< volScalarField > | fMu () const |
| tmp< volScalarField > | f2 () const |
| tmp< volScalarField > | E (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 volScalarField & | y_ |
| 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 |
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.
Definition at line 76 of file LienCubicKE.H.
| 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 118 of file LienCubicKE.C.
References AE_, alpha, Anu_, Foam::bound(), Cbeta1_, Cbeta2_, Cbeta3_, Cbeta_, Ceps1_, Ceps2_, Cgamma1_, Cgamma2_, Cgamma4_, Cmu1_, Cmu2_, Cmu_, epsilon_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, k_, kappa_, Foam::incompressible::New(), nonlinearEddyViscosity< incompressible::RASModel >::nonlinearEddyViscosity(), phi, rho, sigmaEps_, sigmak_, timeName, U, y, and y_.

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


|
protected |
Definition at line 53 of file LienCubicKE.C.
References epsilon_, Foam::exp(), k_, nu, and Foam::sqr().
Referenced by correct(), and E().


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


|
protectedvirtual |
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 75 of file LienCubicKE.C.
References correctNonlinearStress(), and Foam::fvc::grad().

|
protectedvirtual |
Implements nonlinearEddyViscosity< incompressible::RASModel >.
Definition at line 81 of file LienCubicKE.C.
References Cbeta1_, Cbeta2_, Cbeta3_, Cbeta_, Cgamma1_, Cgamma2_, Cgamma4_, Cmu1_, Cmu2_, Foam::dev(), Foam::devSymm(), epsilon_, fMu(), Foam::innerSqr(), k_, Foam::mag(), Foam::magSqr(), nonlinearEddyViscosity< incompressible::RASModel >::nonlinearStress_, eddyViscosity< incompressible::RASModel >::nut_, Foam::pow3(), Foam::skew(), Foam::sqr(), Foam::sqrt(), Foam::symm(), and Foam::twoSymm().
Referenced by correct(), and correctNut().


| TypeName | ( | "LienCubicKE" | ) |
Runtime type information.
References alpha, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, phi, turbulenceModel::propertiesName, rho, and U.
|
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_.

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

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

|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 195 of file LienCubicKE.H.
References k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 203 of file LienCubicKE.H.
References epsilon_.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 365 of file LienCubicKE.C.
References Foam::bound(), Ceps1_, Ceps2_, eddyViscosity< BasicTurbulenceModel >::correct(), correctNonlinearStress(), Foam::fvm::ddt(), DepsilonEff(), Foam::fvm::div(), DkEff(), E(), epsilon_, f2(), Foam::fvc::grad(), k_, Foam::fvm::laplacian(), nonlinearEddyViscosity< incompressible::RASModel >::nonlinearStress_, eddyViscosity< incompressible::RASModel >::nut_, tmp< T >::ref(), solve(), Foam::fvm::Sp(), and Foam::twoSymm().

|
protected |
Definition at line 87 of file LienCubicKE.H.
Referenced by correct(), LienCubicKE(), and read().
|
protected |
Definition at line 88 of file LienCubicKE.H.
Referenced by correct(), E(), LienCubicKE(), and read().
|
protected |
Definition at line 89 of file LienCubicKE.H.
Referenced by DkEff(), LienCubicKE(), and read().
|
protected |
Definition at line 90 of file LienCubicKE.H.
Referenced by DepsilonEff(), LienCubicKE(), and read().
|
protected |
Definition at line 91 of file LienCubicKE.H.
Referenced by correctNonlinearStress(), LienCubicKE(), and read().
|
protected |
Definition at line 92 of file LienCubicKE.H.
Referenced by correctNonlinearStress(), LienCubicKE(), and read().
|
protected |
Definition at line 93 of file LienCubicKE.H.
Referenced by correctNonlinearStress(), LienCubicKE(), and read().
|
protected |
Definition at line 94 of file LienCubicKE.H.
Referenced by correctNonlinearStress(), LienCubicKE(), and read().
|
protected |
Definition at line 95 of file LienCubicKE.H.
Referenced by correctNonlinearStress(), LienCubicKE(), and read().
|
protected |
Definition at line 96 of file LienCubicKE.H.
Referenced by correctNonlinearStress(), LienCubicKE(), and read().
|
protected |
Definition at line 97 of file LienCubicKE.H.
Referenced by correctNonlinearStress(), LienCubicKE(), and read().
|
protected |
Definition at line 98 of file LienCubicKE.H.
Referenced by correctNonlinearStress(), LienCubicKE(), and read().
|
protected |
Definition at line 99 of file LienCubicKE.H.
Referenced by correctNonlinearStress(), LienCubicKE(), and read().
|
protected |
Definition at line 101 of file LienCubicKE.H.
Referenced by E(), fMu(), LienCubicKE(), and read().
|
protected |
Definition at line 102 of file LienCubicKE.H.
Referenced by E(), fMu(), LienCubicKE(), and read().
|
protected |
Definition at line 104 of file LienCubicKE.H.
Referenced by fMu(), LienCubicKE(), and read().
|
protected |
Definition at line 105 of file LienCubicKE.H.
Referenced by E(), LienCubicKE(), and read().
|
protected |
Definition at line 110 of file LienCubicKE.H.
Referenced by correct(), correctNonlinearStress(), E(), f2(), fMu(), k(), and LienCubicKE().
|
protected |
Definition at line 111 of file LienCubicKE.H.
Referenced by correct(), correctNonlinearStress(), E(), epsilon(), f2(), and LienCubicKE().
|
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().