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

One-equation SGS model for the gas-phase in a two-phase system supporting phase-inversion. More...

#include <continuousGasKEqn.H>

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

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel
Public Types inherited from kEqn< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel
Public Types inherited from LESeddyViscosity< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel
Public Types inherited from eddyViscosity< LESModel< 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 ("continuousGasKEqn")
 Runtime type information.
 continuousGasKEqn (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 ~continuousGasKEqn ()=default
 Destructor.
virtual bool read ()
 Read model coefficients if they have changed.
Public Member Functions inherited from kEqn< BasicTurbulenceModel >
 TypeName ("kEqn")
 Runtime type information.
 kEqn (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)
 Constructor from components.
virtual ~kEqn ()=default
 Destructor.
virtual tmp< volScalarFieldk () const
 Return SGS kinetic energy.
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k.
virtual void correct ()
 Correct eddy-Viscosity and related properties.
Public Member Functions inherited from LESeddyViscosity< BasicTurbulenceModel >
 LESeddyViscosity (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
 Construct from components.
virtual ~LESeddyViscosity ()=default
 Destructor.
Public Member Functions inherited from eddyViscosity< LESModel< 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

const turbulenceModelliquidTurbulence () const
 Return the turbulence model for the liquid phase.
tmp< volScalarFieldphaseTransferCoeff () const
virtual tmp< fvScalarMatrixkSource () const
Protected Member Functions inherited from kEqn< BasicTurbulenceModel >
 kEqn (const kEqn &)=delete
 No copy construct.
void operator= (const kEqn &)=delete
 No copy assignment.
virtual void correctNut ()

Protected Attributes

dimensionedScalar alphaInversion_
Protected Attributes inherited from kEqn< BasicTurbulenceModel >
volScalarField k_
dimensionedScalar Ck_
Protected Attributes inherited from eddyViscosity< LESModel< BasicTurbulenceModel > >
volScalarField nut_

Detailed Description

template<class BasicTurbulenceModel>
class Foam::LESModels::continuousGasKEqn< BasicTurbulenceModel >

One-equation SGS model for the gas-phase in a two-phase system supporting phase-inversion.

In the limit that the gas-phase fraction approaches zero a contribution from the other phase is blended into the k-equation up to the phase-fraction of alphaInversion at which point phase-inversion is considered to have occurred and the model reverts to the pure single-phase form.

This model is unpublished and is provided as a stable numerical framework on which a more physical model may be built.

The default model coefficients are

    continuousKEqnCoeffs
    {
        Ck              0.094;
        Ce              1.048;
        alphaInversion  0.7;
    }
Source files

Definition at line 72 of file continuousGasKEqn.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 116 of file continuousGasKEqn.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 117 of file continuousGasKEqn.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 118 of file continuousGasKEqn.H.

Constructor & Destructor Documentation

◆ continuousGasKEqn()

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

References alpha, alphaInversion_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, kEqn< BasicTurbulenceModel >::kEqn(), phi, rho, and U.

Here is the call graph for this function:

◆ ~continuousGasKEqn()

template<class BasicTurbulenceModel>
virtual ~continuousGasKEqn ( )
virtualdefault

Destructor.

Member Function Documentation

◆ liquidTurbulence()

template<class BasicTurbulenceModel>
const turbulenceModel & liquidTurbulence ( ) const
protected

Return the turbulence model for the liquid phase.

Definition at line 96 of file continuousGasKEqn.C.

References fluid, IOobject::groupName(), turbulenceModel::propertiesName, Foam::refCast(), and U.

Referenced by kSource(), and phaseTransferCoeff().

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 124 of file continuousGasKEqn.C.

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

Referenced by 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 kEqn< BasicTurbulenceModel >.

Definition at line 147 of file continuousGasKEqn.C.

References kEqn< BasicTurbulenceModel >::k_, liquidTurbulence(), phaseTransferCoeff(), and Foam::fvm::Sp().

Here is the call graph for this function:

◆ TypeName()

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

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

Read model coefficients if they have changed.

Reimplemented from kEqn< BasicTurbulenceModel >.

Definition at line 81 of file continuousGasKEqn.C.

References alphaInversion_, and kEqn< BasicTurbulenceModel >::read().

Here is the call graph for this function:

Member Data Documentation

◆ alphaInversion_

template<class BasicTurbulenceModel>
dimensionedScalar alphaInversion_
protected

Definition at line 100 of file continuousGasKEqn.H.

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


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