Lien cubic non-linear low-Reynolds k-epsilon turbulence models for incompressible flows. More...
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< volScalarField > | DkEff () const |
Return the effective diffusivity for k. More... | |
tmp< volScalarField > | DepsilonEff () const |
Return the effective diffusivity for epsilon. More... | |
virtual tmp< volScalarField > | k () const |
Return the turbulence kinetic energy. More... | |
virtual tmp< volScalarField > | epsilon () 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< volSymmTensorField > | R () const |
Return the Reynolds stress tensor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (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... | |
volScalarField & | evNut () |
Return non-const access to the turbulence viscosity. More... | |
virtual tmp< volScalarField > | nut () const |
Return the turbulence viscosity. More... | |
virtual tmp< scalarField > | nut (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< 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. 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 |
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 77 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 123 of file LienCubicKE.C.
References Foam::bound(), LienCubicKE::epsilon_, and LienCubicKE::k_.
Referenced by LienCubicKE::correctNonlinearStress().
|
inlinevirtual |
Destructor.
Definition at line 152 of file LienCubicKE.H.
References LienCubicKE::read().
|
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().
|
protected |
Definition at line 57 of file LienCubicKE.C.
References LienCubicKE::epsilon_, Foam::exp(), LienCubicKE::k_, nu, and Foam::sqr().
Referenced by LienCubicKE::correct().
|
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().
|
protectedvirtual |
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 79 of file LienCubicKE.C.
References LienCubicKE::correctNonlinearStress(), and Foam::fvc::grad().
|
protectedvirtual |
Implements nonlinearEddyViscosity< incompressible::RASModel >.
Definition at line 85 of file LienCubicKE.C.
References LienCubicKE::Cbeta1_, LienCubicKE::Cbeta2_, LienCubicKE::Cbeta3_, LienCubicKE::Cbeta_, LienCubicKE::Cgamma1_, LienCubicKE::Cgamma2_, LienCubicKE::Cgamma4_, LienCubicKE::Cmu1_, LienCubicKE::Cmu2_, GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dev(), LienCubicKE::epsilon_, LienCubicKE::fMu(), Foam::innerSqr(), LienCubicKE::k_, LienCubicKE::LienCubicKE(), 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 LienCubicKE::correct(), and LienCubicKE::correctNut().
TypeName | ( | "LienCubicKE" | ) |
Runtime type information.
|
virtual |
Read RASProperties dictionary.
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 340 of file LienCubicKE.C.
References LienCubicKE::AE_, LienCubicKE::Anu_, LienCubicKE::Cbeta1_, LienCubicKE::Cbeta2_, LienCubicKE::Cbeta3_, LienCubicKE::Cbeta_, LienCubicKE::Ceps1_, LienCubicKE::Ceps2_, LienCubicKE::Cgamma1_, LienCubicKE::Cgamma2_, LienCubicKE::Cgamma4_, LienCubicKE::Cmu1_, LienCubicKE::Cmu2_, LienCubicKE::Cmu_, LienCubicKE::kappa_, dimensioned< Type >::readIfPresent(), LienCubicKE::sigmaEps_, and LienCubicKE::sigmak_.
Referenced by LienCubicKE::~LienCubicKE().
|
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().
|
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().
|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 180 of file LienCubicKE.H.
References LienCubicKE::k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 186 of file LienCubicKE.H.
References LienCubicKE::correct(), and LienCubicKE::epsilon_.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 371 of file LienCubicKE.C.
References Foam::bound(), GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), LienCubicKE::Ceps1_, LienCubicKE::Ceps2_, eddyViscosity< BasicTurbulenceModel >::correct(), LienCubicKE::correctNonlinearStress(), Foam::fvm::ddt(), LienCubicKE::DepsilonEff(), Foam::fvm::div(), LienCubicKE::DkEff(), LienCubicKE::E(), LienCubicKE::epsilon_, LienCubicKE::f2(), Foam::constant::universal::G, Foam::fvc::grad(), LienCubicKE::k_, Foam::fvm::laplacian(), nonlinearEddyViscosity< incompressible::RASModel >::nonlinearStress_, eddyViscosity< incompressible::RASModel >::nut_, tmp< T >::ref(), Foam::solve(), Foam::fvm::Sp(), and Foam::twoSymm().
Referenced by LienCubicKE::epsilon().
|
protected |
Definition at line 88 of file LienCubicKE.H.
Referenced by LienCubicKE::correct(), and LienCubicKE::read().
|
protected |
Definition at line 89 of file LienCubicKE.H.
Referenced by LienCubicKE::correct(), LienCubicKE::E(), and LienCubicKE::read().
|
protected |
Definition at line 90 of file LienCubicKE.H.
Referenced by LienCubicKE::read().
|
protected |
Definition at line 91 of file LienCubicKE.H.
Referenced by LienCubicKE::read().
|
protected |
Definition at line 92 of file LienCubicKE.H.
Referenced by LienCubicKE::correctNonlinearStress(), and LienCubicKE::read().
|
protected |
Definition at line 93 of file LienCubicKE.H.
Referenced by LienCubicKE::correctNonlinearStress(), and LienCubicKE::read().
|
protected |
Definition at line 94 of file LienCubicKE.H.
Referenced by LienCubicKE::correctNonlinearStress(), and LienCubicKE::read().
|
protected |
Definition at line 95 of file LienCubicKE.H.
Referenced by LienCubicKE::correctNonlinearStress(), and LienCubicKE::read().
|
protected |
Definition at line 96 of file LienCubicKE.H.
Referenced by LienCubicKE::correctNonlinearStress(), and LienCubicKE::read().
|
protected |
Definition at line 97 of file LienCubicKE.H.
Referenced by LienCubicKE::correctNonlinearStress(), and LienCubicKE::read().
|
protected |
Definition at line 98 of file LienCubicKE.H.
Referenced by LienCubicKE::correctNonlinearStress(), and LienCubicKE::read().
|
protected |
Definition at line 99 of file LienCubicKE.H.
Referenced by LienCubicKE::correctNonlinearStress(), and LienCubicKE::read().
|
protected |
Definition at line 100 of file LienCubicKE.H.
Referenced by LienCubicKE::correctNonlinearStress(), and LienCubicKE::read().
|
protected |
Definition at line 102 of file LienCubicKE.H.
Referenced by LienCubicKE::E(), LienCubicKE::fMu(), and LienCubicKE::read().
|
protected |
Definition at line 103 of file LienCubicKE.H.
Referenced by LienCubicKE::E(), LienCubicKE::fMu(), and LienCubicKE::read().
|
protected |
Definition at line 105 of file LienCubicKE.H.
Referenced by LienCubicKE::fMu(), and LienCubicKE::read().
|
protected |
Definition at line 106 of file LienCubicKE.H.
Referenced by LienCubicKE::E(), and LienCubicKE::read().
|
protected |
Definition at line 111 of file LienCubicKE.H.
Referenced by LienCubicKE::correct(), LienCubicKE::correctNonlinearStress(), LienCubicKE::E(), LienCubicKE::f2(), LienCubicKE::fMu(), LienCubicKE::k(), and LienCubicKE::LienCubicKE().
|
protected |
Definition at line 112 of file LienCubicKE.H.
Referenced by LienCubicKE::correct(), LienCubicKE::correctNonlinearStress(), LienCubicKE::E(), LienCubicKE::epsilon(), LienCubicKE::f2(), and LienCubicKE::LienCubicKE().
|
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().