Loading...
Searching...
No Matches
LaheyKEpsilon< BasicTurbulenceModel > Class Template Reference

Continuous-phase k-epsilon model including bubble-generated turbulence. More...

#include <LaheyKEpsilon.H>

Inheritance diagram for LaheyKEpsilon< BasicTurbulenceModel >:
Collaboration diagram for LaheyKEpsilon< BasicTurbulenceModel >:

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< volScalarFieldDkEff () const
 Return the effective diffusivity for k.
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon.
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy.
virtual tmp< volScalarFieldepsilon () 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< volScalarFieldnut () const
 Return the turbulence viscosity.
virtual tmp< volSymmTensorFieldR () 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< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field.
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation.
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation.

Protected Member Functions

virtual void correctNut ()
tmp< volScalarFieldbubbleG () const
tmp< volScalarFieldphaseTransferCoeff () const
virtual tmp< fvScalarMatrixkSource () const
virtual tmp< fvScalarMatrixepsilonSource () 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::InternallEpsPtr_
std::unique_ptr< volScalarFieldlambdaEpsPtr_
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_

Detailed Description

template<class BasicTurbulenceModel>
class Foam::RASModels::LaheyKEpsilon< BasicTurbulenceModel >

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;
    }
Source files

Definition at line 76 of file LaheyKEpsilon.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 133 of file LaheyKEpsilon.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 134 of file LaheyKEpsilon.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 135 of file LaheyKEpsilon.H.

Constructor & Destructor Documentation

◆ LaheyKEpsilon()

template<class BasicTurbulenceModel>
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 )

◆ ~LaheyKEpsilon()

template<class BasicTurbulenceModel>
virtual ~LaheyKEpsilon ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

template<class BasicTurbulenceModel>
void correctNut ( )
protectedvirtual

◆ bubbleG()

template<class BasicTurbulenceModel>
tmp< volScalarField > bubbleG ( ) const
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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ phaseTransferCoeff()

template<class BasicTurbulenceModel>
tmp< volScalarField > phaseTransferCoeff ( ) const
protected

Definition at line 211 of file LaheyKEpsilon.C.

References alpha, alphaInversion_, Foam::max(), Foam::min(), rho, and U.

Referenced by epsilonSource(), and kSource().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ kSource()

template<class BasicTurbulenceModel>
tmp< fvScalarMatrix > kSource ( ) const
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().

Here is the call graph for this function:

◆ epsilonSource()

template<class BasicTurbulenceModel>
tmp< fvScalarMatrix > epsilonSource ( ) const
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().

Here is the call graph for this function:

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "LaheyKEpsilon< BasicTurbulenceModel >" )

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
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().

Here is the call graph for this function:

◆ correct()

template<class BasicTurbulenceModel>
void correct ( )
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().

Here is the call graph for this function:

Member Data Documentation

◆ alphaInversion_

template<class BasicTurbulenceModel>
dimensionedScalar alphaInversion_
protected

Definition at line 116 of file LaheyKEpsilon.H.

Referenced by LaheyKEpsilon(), phaseTransferCoeff(), and read().

◆ Cp_

template<class BasicTurbulenceModel>
dimensionedScalar Cp_
protected

Definition at line 117 of file LaheyKEpsilon.H.

Referenced by bubbleG(), LaheyKEpsilon(), and read().

◆ C3_

template<class BasicTurbulenceModel>
dimensionedScalar C3_
protected

Definition at line 118 of file LaheyKEpsilon.H.

Referenced by epsilonSource(), LaheyKEpsilon(), and read().

◆ Cmub_

template<class BasicTurbulenceModel>
dimensionedScalar Cmub_
protected

Definition at line 119 of file LaheyKEpsilon.H.

Referenced by correctNut(), LaheyKEpsilon(), and read().


The documentation for this class was generated from the following files: