Lam and Bremhorst low-Reynolds number k-epsilon turbulence model for incompressible flows. More...
Public Member Functions | |
TypeName ("LamBremhorstKE") | |
Runtime type information. More... | |
LamBremhorstKE (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 | ~LamBremhorstKE () |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. 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... | |
![]() | |
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 tmp< volScalarField > | k () const=0 |
Return the turbulence kinetic energy. More... | |
virtual tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
![]() | |
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... | |
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... | |
Protected Member Functions | |
virtual void | correctNut () |
Protected Attributes | |
dimensionedScalar | Cmu_ |
dimensionedScalar | Ceps1_ |
dimensionedScalar | Ceps2_ |
dimensionedScalar | sigmaEps_ |
volScalarField | k_ |
volScalarField | epsilon_ |
const volScalarField & | y_ |
Wall distance. More... | |
![]() | |
volScalarField | nut_ |
Private Member Functions | |
LamBremhorstKE (const LamBremhorstKE &) | |
LamBremhorstKE & | operator= (const LamBremhorstKE &) |
tmp< volScalarField > | Rt () const |
tmp< volScalarField > | fMu (const volScalarField &Rt) const |
tmp< volScalarField > | f1 (const volScalarField &fMu) const |
tmp< volScalarField > | f2 (const volScalarField &Rt) const |
void | correctNut (const volScalarField &fMu) |
Additional Inherited Members | |
![]() | |
typedef incompressible::RASModel ::alphaField | alphaField |
typedef incompressible::RASModel ::rhoField | rhoField |
typedef incompressible::RASModel ::transportModel | transportModel |
![]() | |
typedef incompressible::RASModel ::alphaField | alphaField |
typedef incompressible::RASModel ::rhoField | rhoField |
typedef incompressible::RASModel ::transportModel | transportModel |
Lam and Bremhorst low-Reynolds number k-epsilon turbulence model for incompressible flows.
This turbulence model is described in:
Lam, C. K. G., & Bremhorst, K. (1981). A modified form of the k-ε model for predicting wall turbulence. Journal of Fluids Engineering, 103(3), 456-460.
Definition at line 64 of file LamBremhorstKE.H.
|
private |
LamBremhorstKE | ( | 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 90 of file LamBremhorstKE.C.
References Foam::bound(), and Foam::type().
|
inlinevirtual |
Destructor.
Definition at line 127 of file LamBremhorstKE.H.
|
private |
|
private |
Definition at line 47 of file LamBremhorstKE.C.
References LamBremhorstKE::epsilon_, LamBremhorstKE::k_, nu, and Foam::sqr().
Referenced by LamBremhorstKE::correct(), LamBremhorstKE::correctNut(), LamBremhorstKE::f2(), and LamBremhorstKE::fMu().
|
private |
Definition at line 53 of file LamBremhorstKE.C.
References Foam::exp(), LamBremhorstKE::k_, nu, LamBremhorstKE::Rt(), Foam::sqr(), Foam::sqrt(), and LamBremhorstKE::y_.
Referenced by LamBremhorstKE::correct(), LamBremhorstKE::correctNut(), and LamBremhorstKE::f1().
|
private |
Definition at line 60 of file LamBremhorstKE.C.
References LamBremhorstKE::fMu(), and Foam::pow3().
Referenced by LamBremhorstKE::correct().
|
private |
Definition at line 66 of file LamBremhorstKE.C.
References Foam::exp(), LamBremhorstKE::Rt(), and Foam::sqr().
Referenced by LamBremhorstKE::correct().
|
private |
Definition at line 72 of file LamBremhorstKE.C.
References LamBremhorstKE::Cmu_, GeometricField::correctBoundaryConditions(), LamBremhorstKE::epsilon_, LamBremhorstKE::fMu(), LamBremhorstKE::k_, eddyViscosity< incompressible::RASModel >::nut_, and Foam::sqr().
|
protectedvirtual |
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 81 of file LamBremhorstKE.C.
References LamBremhorstKE::fMu(), and LamBremhorstKE::Rt().
Referenced by LamBremhorstKE::correct().
TypeName | ( | "LamBremhorstKE" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 190 of file LamBremhorstKE.C.
References LamBremhorstKE::Ceps1_, LamBremhorstKE::Ceps2_, LamBremhorstKE::Cmu_, dimensioned::readIfPresent(), and LamBremhorstKE::sigmaEps_.
|
inline |
Return the effective diffusivity for k.
Definition at line 137 of file LamBremhorstKE.H.
References nu, and eddyViscosity< incompressible::RASModel >::nut_.
Referenced by LamBremhorstKE::correct().
|
inline |
Return the effective diffusivity for epsilon.
Definition at line 146 of file LamBremhorstKE.H.
References nu, eddyViscosity< incompressible::RASModel >::nut_, and LamBremhorstKE::sigmaEps_.
Referenced by LamBremhorstKE::correct().
|
inlinevirtual |
Return the turbulence kinetic energy.
Definition at line 155 of file LamBremhorstKE.H.
References LamBremhorstKE::k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 161 of file LamBremhorstKE.H.
References LamBremhorstKE::epsilon_.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< incompressible::RASModel >.
Definition at line 208 of file LamBremhorstKE.C.
References Foam::bound(), GeometricField::boundaryField(), LamBremhorstKE::Ceps1_, LamBremhorstKE::Ceps2_, tmp::clear(), eddyViscosity< BasicTurbulenceModel >::correct(), LamBremhorstKE::correctNut(), Foam::fvm::ddt(), LamBremhorstKE::DepsilonEff(), Foam::fvm::div(), LamBremhorstKE::DkEff(), LamBremhorstKE::epsilon_, LamBremhorstKE::f1(), LamBremhorstKE::f2(), LamBremhorstKE::fMu(), Foam::constant::universal::G, Foam::fvc::grad(), LamBremhorstKE::k_, Foam::fvm::laplacian(), eddyViscosity< incompressible::RASModel >::nut_, LamBremhorstKE::Rt(), Foam::solve(), Foam::fvm::Sp(), Foam::twoSymm(), and GeometricField::GeometricBoundaryField::updateCoeffs().
|
protected |
Definition at line 85 of file LamBremhorstKE.H.
Referenced by LamBremhorstKE::correctNut(), and LamBremhorstKE::read().
|
protected |
Definition at line 86 of file LamBremhorstKE.H.
Referenced by LamBremhorstKE::correct(), and LamBremhorstKE::read().
|
protected |
Definition at line 87 of file LamBremhorstKE.H.
Referenced by LamBremhorstKE::correct(), and LamBremhorstKE::read().
|
protected |
Definition at line 88 of file LamBremhorstKE.H.
Referenced by LamBremhorstKE::DepsilonEff(), and LamBremhorstKE::read().
|
protected |
Definition at line 90 of file LamBremhorstKE.H.
Referenced by LamBremhorstKE::correct(), LamBremhorstKE::correctNut(), LamBremhorstKE::fMu(), LamBremhorstKE::k(), and LamBremhorstKE::Rt().
|
protected |
Definition at line 91 of file LamBremhorstKE.H.
Referenced by LamBremhorstKE::correct(), LamBremhorstKE::correctNut(), LamBremhorstKE::epsilon(), and LamBremhorstKE::Rt().
|
protected |
Wall distance.
Note: different to wall distance in parent RASModel which is for near-wall cells only
Definition at line 96 of file LamBremhorstKE.H.
Referenced by LamBremhorstKE::fMu().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.