Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
continuousGasKEpsilon< BasicTurbulenceModel > Class Template Reference

k-epsilon model for the gas-phase in a two-phase system supporting phase-inversion. More...

Inheritance diagram for continuousGasKEpsilon< BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for continuousGasKEpsilon< BasicTurbulenceModel >:
Collaboration graph
[legend]

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 RASModel< BasicTurbulenceModel > ::alphaField alphaField
 
typedef RASModel< BasicTurbulenceModel > ::rhoField rhoField
 
typedef RASModel< BasicTurbulenceModel > ::transportModel transportModel
 
- Public Types inherited from linearViscousStress< RASModel< BasicTurbulenceModel > >
typedef RASModel< BasicTurbulenceModel > ::alphaField alphaField
 
typedef RASModel< BasicTurbulenceModel > ::rhoField rhoField
 
typedef RASModel< BasicTurbulenceModel > ::transportModel transportModel
 
- Public Types inherited from RASModel< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 

Public Member Functions

 TypeName ("continuousGasKEpsilon")
 Runtime type information. More...
 
 continuousGasKEpsilon (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. More...
 
virtual ~continuousGasKEpsilon ()
 Destructor. More...
 
virtual bool read ()
 Re-read model coefficients if they have changed. More...
 
const turbulenceModelliquidTurbulence () const
 Return the turbulence model for the liquid phase. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< volScalarFieldrhoEff () const
 Return the effective density for the stress. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
- Public Member Functions inherited from kEpsilon< BasicTurbulenceModel >
 TypeName ("kEpsilon")
 Runtime type information. More...
 
 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. More...
 
virtual ~kEpsilon ()
 Destructor. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
- 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. More...
 
virtual ~eddyViscosity ()
 Destructor. More...
 
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity. More...
 
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
- Public Member Functions inherited from linearViscousStress< RASModel< 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. More...
 
virtual ~linearViscousStress ()
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
- Public Member Functions inherited from RASModel< BasicTurbulenceModel >
 TypeName ("RAS")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, RASModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName))
 
 RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 Construct from components. More...
 
virtual ~RASModel ()
 Destructor. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: small) More...
 
const dimensionedScalarepsilonMin () const
 Return the lower allowable limit for epsilon (default: small) More...
 
const dimensionedScalaromegaMin () const
 Return the lower allowable limit for omega (default: small) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
dimensionedScalarepsilonMin ()
 Allow epsilonMin to be changed. More...
 
dimensionedScalaromegaMin ()
 Allow omegaMin to be changed. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 

Protected Member Functions

virtual void correctNut ()
 
tmp< volScalarFieldphaseTransferCoeff () const
 
virtual tmp< fvScalarMatrixkSource () const
 
virtual tmp< fvScalarMatrixepsilonSource () const
 
- Protected Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 

Protected Attributes

dimensionedScalar alphaInversion_
 
- Protected Attributes inherited from kEpsilon< BasicTurbulenceModel >
dimensionedScalar Cmu_
 
dimensionedScalar C1_
 
dimensionedScalar C2_
 
dimensionedScalar C3_
 
dimensionedScalar sigmak_
 
dimensionedScalar sigmaEps_
 
volScalarField k_
 
volScalarField epsilon_
 
- Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_
 
- Protected Attributes inherited from RASModel< BasicTurbulenceModel >
dictionary RASDict_
 RAS coefficients dictionary. More...
 
Switch turbulence_
 Turbulence on/off flag. More...
 
Switch printCoeffs_
 Flag to print the model coeffs at run-time. More...
 
dictionary coeffDict_
 Model coefficients dictionary. More...
 
dimensionedScalar kMin_
 Lower limit of k. More...
 
dimensionedScalar epsilonMin_
 Lower limit of epsilon. More...
 
dimensionedScalar omegaMin_
 Lower limit for omega. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from RASModel< BasicTurbulenceModel >
static autoPtr< RASModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
 Return a reference to the selected RAS model. More...
 

Detailed Description

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

k-epsilon 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 and epsilon equations 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

    continuousGasKEpsilonCoeffs
    {
        Cmu             0.09;
        C1              1.44;
        C2              1.92;
        C3              -0.33;
        sigmak          1.0;
        sigmaEps        1.3;
        alphaInversion  0.7;
    }
Source files

Definition at line 75 of file continuousGasKEpsilon.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 112 of file continuousGasKEpsilon.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 113 of file continuousGasKEpsilon.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 114 of file continuousGasKEpsilon.H.

Constructor & Destructor Documentation

◆ continuousGasKEpsilon()

continuousGasKEpsilon ( 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 43 of file continuousGasKEpsilon.C.

◆ ~continuousGasKEpsilon()

virtual ~continuousGasKEpsilon ( )
inlinevirtual

Member Function Documentation

◆ correctNut()

void correctNut ( )
protectedvirtual

◆ phaseTransferCoeff()

tmp< volScalarField > phaseTransferCoeff ( ) const
protected

Definition at line 229 of file continuousGasKEpsilon.C.

References TimeState::deltaT(), turbulenceModel::epsilon(), turbulenceModel::k(), Foam::max(), Foam::min(), and IOobject::time().

Here is the call graph for this function:

◆ kSource()

tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Reimplemented from kEpsilon< BasicTurbulenceModel >.

Definition at line 252 of file continuousGasKEpsilon.C.

References turbulenceModel::k(), and Foam::fvm::Sp().

Here is the call graph for this function:

◆ epsilonSource()

tmp< fvScalarMatrix > epsilonSource ( ) const
protectedvirtual

Reimplemented from kEpsilon< BasicTurbulenceModel >.

Definition at line 265 of file continuousGasKEpsilon.C.

References turbulenceModel::epsilon(), and Foam::fvm::Sp().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "continuousGasKEpsilon< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Re-read model coefficients if they have changed.

Reimplemented from kEpsilon< BasicTurbulenceModel >.

Definition at line 101 of file continuousGasKEpsilon.C.

Referenced by continuousGasKEpsilon< BasicTurbulenceModel >::~continuousGasKEpsilon().

Here is the caller graph for this function:

◆ liquidTurbulence()

const turbulenceModel & liquidTurbulence ( ) const

Return the turbulence model for the liquid phase.

Definition at line 149 of file continuousGasKEpsilon.C.

References IOobject::db(), fluid, IOobject::groupName(), objectRegistry::lookupObject(), twoPhaseSystem::otherPhase(), and turbulenceModel::propertiesName.

Referenced by continuousGasKEpsilon< BasicTurbulenceModel >::~continuousGasKEpsilon().

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

◆ nuEff()

tmp< Foam::volScalarField > nuEff ( ) const
virtual

Return the effective viscosity.

Reimplemented from RASModel< BasicTurbulenceModel >.

Definition at line 177 of file continuousGasKEpsilon.C.

References IOobject::groupName(), Foam::max(), Foam::min(), and nu.

Referenced by continuousGasKEpsilon< BasicTurbulenceModel >::~continuousGasKEpsilon().

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

◆ rhoEff()

tmp< Foam::volScalarField > rhoEff ( ) const
virtual

Return the effective density for the stress.

Definition at line 207 of file continuousGasKEpsilon.C.

References virtualMassModel::Cvm(), fluid, IOobject::groupName(), twoPhaseSystem::lookupSubModel(), and twoPhaseSystem::otherPhase().

Referenced by continuousGasKEpsilon< BasicTurbulenceModel >::~continuousGasKEpsilon().

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

◆ R()

tmp< volSymmTensorField > R ( ) const
virtual

Return the Reynolds stress tensor.

Reimplemented from eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 278 of file continuousGasKEpsilon.C.

References Foam::dev(), Foam::fvc::grad(), IOobject::groupName(), Foam::I, k, IOobject::NO_READ, IOobject::NO_WRITE, and Foam::twoSymm().

Referenced by continuousGasKEpsilon< BasicTurbulenceModel >::~continuousGasKEpsilon().

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

Member Data Documentation

◆ alphaInversion_

dimensionedScalar alphaInversion_
protected

Definition at line 99 of file continuousGasKEpsilon.H.


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