Continuous-phase k-epsilon model including bubble-generated turbulence. More...
#include <LaheyKEpsilon.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from kEpsilon< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from linearViscousStress< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("LaheyKEpsilon") | |
| Runtime type information. | |
| LaheyKEpsilon (const alphaField &alpha, const rhoField &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 | ~LaheyKEpsilon ()=default |
| Destructor. | |
| virtual bool | read () |
| Read model coefficients if they have changed. | |
| virtual void | correct () |
| Solve the turbulence equations and correct the turbulence viscosity. | |
| Public Member Functions inherited from kEpsilon< BasicTurbulenceModel > | |
| TypeName ("kEpsilon") | |
| Runtime type information. | |
| kEpsilon (const alphaField &alpha, const rhoField &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 | ~kEpsilon ()=default |
| Destructor. | |
| 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. | |
| Public Member Functions inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| 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 tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. | |
| virtual void | validate () |
| Validate the turbulence fields after construction. | |
| Public Member Functions inherited from linearViscousStress< BasicTurbulenceModel > | |
| 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. | |
| virtual tmp< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
| Return the effective stress tensor based on a given velocity field. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. | |
Protected Member Functions | |
| virtual void | correctNut () |
| tmp< volScalarField > | bubbleG () const |
| tmp< volScalarField > | phaseTransferCoeff () const |
| virtual tmp< fvScalarMatrix > | kSource () const |
| virtual tmp< fvScalarMatrix > | epsilonSource () const |
Protected Attributes | |
| dimensionedScalar | alphaInversion_ |
| dimensionedScalar | Cp_ |
| dimensionedScalar | C3_ |
| dimensionedScalar | Cmub_ |
| Protected Attributes inherited from kEpsilon< BasicTurbulenceModel > | |
| dimensionedScalar | Cmu_ |
| dimensionedScalar | C1_ |
| dimensionedScalar | C2_ |
| dimensionedScalar | C3_ |
| dimensionedScalar | sigmak_ |
| dimensionedScalar | sigmaEps_ |
| volScalarField | k_ |
| volScalarField | epsilon_ |
| bool | twoLayerTreatment_ |
| dimensionedScalar | ReyStar_ |
| dimensionedScalar | ReyFactor_ |
| std::unique_ptr< volScalarField::Internal > | lEpsPtr_ |
| std::unique_ptr< volScalarField > | lambdaEpsPtr_ |
| Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| volScalarField | nut_ |
Continuous-phase k-epsilon model including bubble-generated turbulence.
Reference:
Lahey Jr, R. T. (2005).
The simulation of multidimensional multiphase flows.
Nuclear Engineering and Design, 235(10), 1043-1060.
The default model coefficients are
LaheyKEpsilonCoeffs
{
Cmu 0.09;
C1 1.44;
C2 1.92;
C3 -0.33;
sigmak 1.0;
sigmaEps 1.3;
Cp 0.25;
Cmub 0.6;
alphaInversion 0.3;
}
Definition at line 76 of file LaheyKEpsilon.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 133 of file LaheyKEpsilon.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 134 of file LaheyKEpsilon.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 135 of file LaheyKEpsilon.H.
| LaheyKEpsilon | ( | const alphaField & | alpha, |
| const rhoField & | 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 37 of file LaheyKEpsilon.C.
References alpha, alphaInversion_, kEpsilon< BasicTurbulenceModel >::C2_, C3_, Cmub_, Cp_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, phi, rho, and U.
|
virtualdefault |
Destructor.
|
protectedvirtual |
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 162 of file LaheyKEpsilon.C.
References kEpsilon< BasicTurbulenceModel >::Cmu_, Cmub_, optionList::correct(), kEpsilon< BasicTurbulenceModel >::epsilon_, kEpsilon< BasicTurbulenceModel >::k_, Foam::mag(), options::New(), eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_, and Foam::sqr().

|
protected |
Definition at line 180 of file LaheyKEpsilon.C.
References bubbleG(), Cp_, drag, fluid, Foam::mag(), Foam::pow(), Foam::pow3(), and Foam::refCast().
Referenced by bubbleG(), epsilonSource(), and kSource().


|
protected |
Definition at line 211 of file LaheyKEpsilon.C.
References alpha, alphaInversion_, Foam::max(), Foam::min(), rho, and U.
Referenced by epsilonSource(), and kSource().


|
protectedvirtual |
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 229 of file LaheyKEpsilon.C.
References alpha, bubbleG(), kEpsilon< BasicTurbulenceModel >::k_, phaseTransferCoeff(), rho, and Foam::fvm::Sp().

|
protectedvirtual |
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 247 of file LaheyKEpsilon.C.
References alpha, bubbleG(), C3_, kEpsilon< BasicTurbulenceModel >::epsilon_, kEpsilon< BasicTurbulenceModel >::k_, phaseTransferCoeff(), rho, and Foam::fvm::Sp().

| TypeName | ( | "LaheyKEpsilon< BasicTurbulenceModel >" | ) |
Runtime type information.
References alpha, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, phi, turbulenceModel::propertiesName, rho, and U.
|
virtual |
Read model coefficients if they have changed.
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 113 of file LaheyKEpsilon.C.
References alphaInversion_, C3_, Cmub_, Cp_, and kEpsilon< BasicTurbulenceModel >::read().

|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 265 of file LaheyKEpsilon.C.
References kEpsilon< BasicTurbulenceModel >::correct().

|
protected |
Definition at line 116 of file LaheyKEpsilon.H.
Referenced by LaheyKEpsilon(), phaseTransferCoeff(), and read().
|
protected |
Definition at line 117 of file LaheyKEpsilon.H.
Referenced by bubbleG(), LaheyKEpsilon(), and read().
|
protected |
Definition at line 118 of file LaheyKEpsilon.H.
Referenced by epsilonSource(), LaheyKEpsilon(), and read().
|
protected |
Definition at line 119 of file LaheyKEpsilon.H.
Referenced by correctNut(), LaheyKEpsilon(), and read().