Lien cubic non-linear low-Reynolds k-epsilon turbulence models for incompressible flows. More...
#include <incompressibleRASModels/LienCubicKELowRe.H>
Lien cubic non-linear low-Reynolds k-epsilon turbulence models for incompressible flows.
Lien, F.S., Chen, W.L., Leschziner, M.A., "Low-Reynolds-number eddy-viscosity modeling based on non-linear stress-strain/vorticity relations" Engineering Turbulence Modelling and Experiments 3 (Edited by Rodi, W. and Bergeles, G.), 91-100. 1996. Elsevier Science Publishers. Etemad, S., et al., "Turbulent flow and heat transfer in a square-sectioned U bend" Progress in compuational fluid dynamics 6, 89-100. 2006.
Definition at line 69 of file LienCubicKELowRe.H.
Public Member Functions | |
TypeName ("LienCubicKELowRe") | |
Runtime type information.
| |
LienCubicKELowRe (const volVectorField &U, const surfaceScalarField &phi, transportModel &transport) | |
Construct from components.
| |
virtual | ~LienCubicKELowRe () |
Destructor.
| |
virtual tmp< volScalarField > | nut () const |
Return the turbulence viscosity.
| |
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 tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor.
| |
virtual tmp< volSymmTensorField > | devReff () const |
Return the effective stress tensor including the laminar stress.
| |
virtual tmp< fvVectorMatrix > | divDevReff (volVectorField &U) const |
Return the source term for the momentum equation.
| |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity.
| |
virtual bool | read () |
Read RASProperties dictionary.
|
LienCubicKELowRe | ( | const volVectorField & | U, |
const surfaceScalarField & | phi, | ||
transportModel & | transport | ||
) |
Construct from components.
Definition at line 49 of file LienCubicKELowRe.C.
References Foam::exp(), Foam::max(), and Foam::sqr().
virtual ~LienCubicKELowRe | ( | ) | [inline, virtual]
|
Destructor.
Definition at line 134 of file LienCubicKELowRe.H.
TypeName | ( | "LienCubicKELowRe" | ) |
Runtime type information.
virtual tmp<volScalarField> nut | ( | ) | const [inline, virtual]
|
Return the turbulence viscosity.
Implements incompressible::RASModel.
Definition at line 141 of file LienCubicKELowRe.H.
tmp<volScalarField> DkEff | ( | ) | const [inline]
|
Return the effective diffusivity for k.
Definition at line 147 of file LienCubicKELowRe.H.
References incompressible::turbulenceModel::nu().
Referenced by LienCubicKELowRe::correct().
tmp<volScalarField> DepsilonEff | ( | ) | const [inline]
|
Return the effective diffusivity for epsilon.
Definition at line 156 of file LienCubicKELowRe.H.
References incompressible::turbulenceModel::nu().
Referenced by LienCubicKELowRe::correct().
virtual tmp<volScalarField> k | ( | ) | const [inline, virtual]
|
Return the turbulence kinetic energy.
Implements incompressible::RASModel.
Definition at line 165 of file LienCubicKELowRe.H.
virtual tmp<volScalarField> epsilon | ( | ) | const [inline, virtual]
|
Return the turbulence kinetic energy dissipation rate.
Implements incompressible::RASModel.
Definition at line 171 of file LienCubicKELowRe.H.
tmp< volSymmTensorField > R | ( | ) | const [virtual]
|
Return the Reynolds stress tensor.
Implements incompressible::RASModel.
Definition at line 303 of file LienCubicKELowRe.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), Foam::I, incompressible::turbulenceModel::mesh_, IOobject::NO_READ, IOobject::NO_WRITE, incompressible::turbulenceModel::runTime_, Time::timeName(), and Foam::twoSymm().
tmp< volSymmTensorField > devReff | ( | ) | const [virtual]
|
Return the effective stress tensor including the laminar stress.
Implements incompressible::RASModel.
Definition at line 324 of file LienCubicKELowRe.C.
References Foam::dev(), Foam::fvc::grad(), incompressible::turbulenceModel::mesh_, IOobject::NO_READ, IOobject::NO_WRITE, incompressible::RASModel::nuEff(), incompressible::turbulenceModel::runTime_, Time::timeName(), Foam::twoSymm(), and incompressible::turbulenceModel::U_.
tmp< fvVectorMatrix > divDevReff | ( | volVectorField & | U ) | const [virtual]
|
Return the source term for the momentum equation.
Implements incompressible::RASModel.
Definition at line 344 of file LienCubicKELowRe.C.
References Foam::dev(), Foam::fvc::div(), Foam::fvc::grad(), Foam::fvm::laplacian(), incompressible::RASModel::nuEff(), and Foam::T().
void correct | ( | ) | [virtual]
|
Solve the turbulence equations and correct the turbulence viscosity.
Implements incompressible::RASModel.
Definition at line 384 of file LienCubicKELowRe.C.
References Foam::bound(), polyMesh::changing(), wallDist::correct(), Foam::fvm::ddt(), LienCubicKELowRe::DepsilonEff(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::fvm::div(), LienCubicKELowRe::DkEff(), incompressible::RASModel::epsilon0_, Foam::exp(), Foam::fvc::grad(), incompressible::RASModel::k0_, Foam::fvm::laplacian(), Foam::magSqr(), Foam::max(), incompressible::turbulenceModel::mesh_, Foam::min(), incompressible::turbulenceModel::nu(), incompressible::turbulenceModel::phi_, Foam::pow(), Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), Foam::symm(), Foam::T(), GeometricField< Type, PatchField, GeoMesh >::T(), incompressible::RASModel::turbulence_, and incompressible::turbulenceModel::U_.
bool read | ( | ) | [virtual]
|
Read RASProperties dictionary.
Implements incompressible::RASModel.
Definition at line 355 of file LienCubicKELowRe.C.
References incompressible::RASModel::coeffDict(), and dimensioned< Type >::readIfPresent().