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

Mixture k-epsilon turbulence model for two-phase gas-liquid systems. More...

#include <mixtureKEpsilon.H>

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

Public Types

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 ("mixtureKEpsilon")
 Runtime type information.
 mixtureKEpsilon (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 ~mixtureKEpsilon ()=default
 Destructor.
virtual bool read ()
 Re-read model coefficients if they have changed.
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy.
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate.
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity.
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

wordList epsilonBoundaryTypes (const volScalarField &epsilon) const
void correctInletOutlet (volScalarField &vsf, const volScalarField &refVsf) const
void initMixtureFields ()
virtual void correctNut ()
tmp< volScalarFieldCt2 () const
tmp< volScalarFieldrholEff () const
tmp< volScalarFieldrhogEff () const
tmp< volScalarFieldrhom () const
tmp< volScalarFieldmix (const volScalarField &fc, const volScalarField &fd) const
tmp< volScalarFieldmixU (const volScalarField &fc, const volScalarField &fd) const
tmp< surfaceScalarFieldmixFlux (const surfaceScalarField &fc, const surfaceScalarField &fd) const
tmp< volScalarFieldbubbleG () const
virtual tmp< fvScalarMatrixkSource () const
virtual tmp< fvScalarMatrixepsilonSource () const
tmp< volScalarFieldDkEff (const volScalarField &nutm) const
 Return the effective diffusivity for k.
tmp< volScalarFieldDepsilonEff (const volScalarField &nutm) const
 Return the effective diffusivity for epsilon.

Protected Attributes

dimensionedScalar Cmu_
dimensionedScalar C1_
dimensionedScalar C2_
dimensionedScalar C3_
dimensionedScalar Cp_
dimensionedScalar sigmak_
dimensionedScalar sigmaEps_
volScalarField k_
volScalarField epsilon_
autoPtr< volScalarFieldCt2_
autoPtr< volScalarFieldrhom_
autoPtr< volScalarFieldkm_
autoPtr< volScalarFieldepsilonm_
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_

Detailed Description

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

Mixture k-epsilon turbulence model for two-phase gas-liquid systems.

The basic structure of the model is based on:

    Behzadi, A., Issa, R. I., & Rusche, H. (2004).
    Modelling of dispersed bubble and droplet flow at high phase fractions.
    Chemical Engineering Science, 59(4), 759-770.

but an effective density for the gas is used in the averaging and an alternative model for bubble-generated turbulence from:

    Lahey Jr, R. T. (2005).
    The simulation of multidimensional multiphase flows.
    Nuclear Engineering and Design, 235(10), 1043-1060.

The default model coefficients are

    mixtureKEpsilonCoeffs
    {
        Cmu         0.09;
        C1          1.44;
        C2          1.92;
        C3          C2;
        sigmak      1.0;
        sigmaEps    1.3;
    }
Source files

Definition at line 82 of file mixtureKEpsilon.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 210 of file mixtureKEpsilon.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 211 of file mixtureKEpsilon.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 212 of file mixtureKEpsilon.H.

Constructor & Destructor Documentation

◆ mixtureKEpsilon()

template<class BasicTurbulenceModel>
mixtureKEpsilon ( 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 42 of file mixtureKEpsilon.C.

References alpha, Foam::bound(), C1_, C2_, C3_, Cmu_, Cp_, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), epsilon_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, k_, phi, rho, sigmaEps_, sigmak_, timeName, and U.

Here is the call graph for this function:

◆ ~mixtureKEpsilon()

template<class BasicTurbulenceModel>
virtual ~mixtureKEpsilon ( )
virtualdefault

Destructor.

Member Function Documentation

◆ epsilonBoundaryTypes()

template<class BasicTurbulenceModel>
wordList epsilonBoundaryTypes ( const volScalarField & epsilon) const
protected

Definition at line 168 of file mixtureKEpsilon.C.

References epsilon(), forAll, Foam::isA(), and GeometricBoundaryField< Type, PatchField, GeoMesh >::types().

Referenced by initMixtureFields().

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

◆ correctInletOutlet()

template<class BasicTurbulenceModel>
void correctInletOutlet ( volScalarField & vsf,
const volScalarField & refVsf ) const
protected

Definition at line 190 of file mixtureKEpsilon.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), forAll, Foam::isA(), and Foam::refCast().

Referenced by initMixtureFields().

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

◆ initMixtureFields()

template<class BasicTurbulenceModel>
void initMixtureFields ( )
protected

◆ correctNut()

template<class BasicTurbulenceModel>
void correctNut ( )
protectedvirtual

◆ Ct2()

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

Definition at line 371 of file mixtureKEpsilon.C.

References beta(), Cmu_, Foam::exp(), fluid, Foam::mag(), Foam::refCast(), liquid::rho(), Foam::sqr(), and Foam::sqrt().

Referenced by correct(), and initMixtureFields().

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

◆ rholEff()

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

Definition at line 398 of file mixtureKEpsilon.C.

References fluid, and Foam::refCast().

Referenced by correct(), mix(), mixFlux(), mixU(), and rhom().

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

◆ rhogEff()

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

Definition at line 407 of file mixtureKEpsilon.C.

References virtualMassModel::Cvm(), fluid, and Foam::refCast().

Referenced by correct(), mix(), mixFlux(), mixU(), and rhom().

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

◆ rhom()

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

Definition at line 420 of file mixtureKEpsilon.C.

References alphal, rhogEff(), and rholEff().

Referenced by correct(), and initMixtureFields().

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

◆ mix()

template<class BasicTurbulenceModel>
tmp< volScalarField > mix ( const volScalarField & fc,
const volScalarField & fd ) const
protected

Definition at line 430 of file mixtureKEpsilon.C.

References alphal, rhogEff(), rholEff(), and rhom_.

Referenced by correct(), and initMixtureFields().

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

◆ mixU()

template<class BasicTurbulenceModel>
tmp< volScalarField > mixU ( const volScalarField & fc,
const volScalarField & fd ) const
protected

Definition at line 444 of file mixtureKEpsilon.C.

References alphal, Ct2_, rhogEff(), and rholEff().

Referenced by correct().

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

◆ mixFlux()

template<class BasicTurbulenceModel>
tmp< surfaceScalarField > mixFlux ( const surfaceScalarField & fc,
const surfaceScalarField & fd ) const
protected

Definition at line 460 of file mixtureKEpsilon.C.

References alphal, Ct2_, Foam::fvc::interpolate(), rhogEff(), and rholEff().

Referenced by correct().

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

◆ bubbleG()

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

Definition at line 482 of file mixtureKEpsilon.C.

References bubbleG(), Cp_, drag, fluid, Foam::mag(), Foam::pow(), Foam::pow3(), Foam::refCast(), and liquid::rho().

Referenced by bubbleG(), 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

Definition at line 520 of file mixtureKEpsilon.C.

References bubbleG(), km_, rhom_, and Foam::fvm::Su().

Referenced by correct().

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

◆ epsilonSource()

template<class BasicTurbulenceModel>
tmp< fvScalarMatrix > epsilonSource ( ) const
protectedvirtual

Definition at line 527 of file mixtureKEpsilon.C.

References bubbleG(), C3_, epsilonm_, km_, rhom_, and Foam::fvm::Su().

Referenced by correct().

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

◆ DkEff()

template<class BasicTurbulenceModel>
tmp< volScalarField > DkEff ( const volScalarField & nutm) const
inlineprotected

Return the effective diffusivity for k.

Definition at line 181 of file mixtureKEpsilon.H.

References sigmak_.

Referenced by correct().

Here is the caller graph for this function:

◆ DepsilonEff()

template<class BasicTurbulenceModel>
tmp< volScalarField > DepsilonEff ( const volScalarField & nutm) const
inlineprotected

Return the effective diffusivity for epsilon.

Definition at line 196 of file mixtureKEpsilon.H.

References sigmaEps_.

Referenced by correct().

Here is the caller graph for this function:

◆ TypeName()

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

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 309 of file mixtureKEpsilon.C.

References C1_, C2_, C3_, Cmu_, Cp_, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), read(), sigmaEps_, and sigmak_.

Referenced by read().

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

◆ k()

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 255 of file mixtureKEpsilon.H.

References k_.

◆ epsilon()

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 263 of file mixtureKEpsilon.H.

References epsilon_.

Referenced by epsilonBoundaryTypes().

Here is the caller graph for this function:

◆ correct()

Member Data Documentation

◆ Cmu_

template<class BasicTurbulenceModel>
dimensionedScalar Cmu_
protected

Definition at line 115 of file mixtureKEpsilon.H.

Referenced by correctNut(), Ct2(), mixtureKEpsilon(), and read().

◆ C1_

template<class BasicTurbulenceModel>
dimensionedScalar C1_
protected

Definition at line 116 of file mixtureKEpsilon.H.

Referenced by correct(), mixtureKEpsilon(), and read().

◆ C2_

template<class BasicTurbulenceModel>
dimensionedScalar C2_
protected

Definition at line 117 of file mixtureKEpsilon.H.

Referenced by correct(), mixtureKEpsilon(), and read().

◆ C3_

template<class BasicTurbulenceModel>
dimensionedScalar C3_
protected

Definition at line 118 of file mixtureKEpsilon.H.

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

◆ Cp_

template<class BasicTurbulenceModel>
dimensionedScalar Cp_
protected

Definition at line 119 of file mixtureKEpsilon.H.

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

◆ sigmak_

template<class BasicTurbulenceModel>
dimensionedScalar sigmak_
protected

Definition at line 120 of file mixtureKEpsilon.H.

Referenced by DkEff(), mixtureKEpsilon(), and read().

◆ sigmaEps_

template<class BasicTurbulenceModel>
dimensionedScalar sigmaEps_
protected

Definition at line 121 of file mixtureKEpsilon.H.

Referenced by DepsilonEff(), mixtureKEpsilon(), and read().

◆ k_

template<class BasicTurbulenceModel>
volScalarField k_
protected

Definition at line 125 of file mixtureKEpsilon.H.

Referenced by correct(), correctNut(), initMixtureFields(), k(), and mixtureKEpsilon().

◆ epsilon_

template<class BasicTurbulenceModel>
volScalarField epsilon_
protected

Definition at line 126 of file mixtureKEpsilon.H.

Referenced by correct(), correctNut(), epsilon(), initMixtureFields(), and mixtureKEpsilon().

◆ Ct2_

template<class BasicTurbulenceModel>
autoPtr<volScalarField> Ct2_
protected

Definition at line 130 of file mixtureKEpsilon.H.

Referenced by correct(), initMixtureFields(), mixFlux(), and mixU().

◆ rhom_

template<class BasicTurbulenceModel>
autoPtr<volScalarField> rhom_
protected

Definition at line 131 of file mixtureKEpsilon.H.

Referenced by correct(), epsilonSource(), initMixtureFields(), kSource(), and mix().

◆ km_

template<class BasicTurbulenceModel>
autoPtr<volScalarField> km_
protected

Definition at line 132 of file mixtureKEpsilon.H.

Referenced by correct(), epsilonSource(), initMixtureFields(), and kSource().

◆ epsilonm_

template<class BasicTurbulenceModel>
autoPtr<volScalarField> epsilonm_
protected

Definition at line 133 of file mixtureKEpsilon.H.

Referenced by correct(), epsilonSource(), and initMixtureFields().


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