continuousGasKEpsilon< BasicMomentumTransportModel > Class Template Reference

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

Inheritance diagram for continuousGasKEpsilon< BasicMomentumTransportModel >:
Collaboration diagram for continuousGasKEpsilon< BasicMomentumTransportModel >:

Public Types

typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::transportModel transportModel
 
- Public Types inherited from kEpsilon< BasicMomentumTransportModel >
typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::transportModel transportModel
 
- Public Types inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > >
typedef RASModel< BasicMomentumTransportModel > ::alphaField alphaField
 
typedef RASModel< BasicMomentumTransportModel > ::rhoField rhoField
 
typedef RASModel< BasicMomentumTransportModel > ::transportModel transportModel
 
- Public Types inherited from linearViscousStress< RASModel< BasicMomentumTransportModel > >
typedef RASModel< BasicMomentumTransportModel > ::alphaField alphaField
 
typedef RASModel< BasicMomentumTransportModel > ::rhoField rhoField
 
typedef RASModel< BasicMomentumTransportModel > ::transportModel transportModel
 
- Public Types inherited from RASModel< BasicMomentumTransportModel >
typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::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 &type=typeName)
 Construct from components. More...
 
 continuousGasKEpsilon (const continuousGasKEpsilon &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~continuousGasKEpsilon ()
 Destructor. More...
 
virtual bool read ()
 Re-read model coefficients if they have changed. More...
 
const momentumTransportModelliquidTurbulence () 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< volSymmTensorFieldsigma () const
 Return the stress tensor [m^2/s^2]. More...
 
void operator= (const continuousGasKEpsilon &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from kEpsilon< BasicMomentumTransportModel >
 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 &type=typeName)
 Construct from components. More...
 
 kEpsilon (const kEpsilon &)=delete
 Disallow default bitwise copy construction. 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...
 
void operator= (const kEpsilon &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > >
 eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 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< BasicMomentumTransportModel > >
 linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
virtual ~linearViscousStress ()
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevTau () const
 Return the effective stress tensor. More...
 
virtual tmp< fvVectorMatrixdivDevTau (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevTau (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
- Public Member Functions inherited from RASModel< BasicMomentumTransportModel >
 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),(alpha, rho, U, alphaRhoPhi, phi, transport))
 
 RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
 RASModel (const RASModel &)=delete
 Disallow default bitwise copy construction. 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...
 
void operator= (const RASModel &)=delete
 Disallow default bitwise assignment. More...
 

Protected Member Functions

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

Protected Attributes

dimensionedScalar alphaInversion_
 
- Protected Attributes inherited from kEpsilon< BasicMomentumTransportModel >
dimensionedScalar Cmu_
 
dimensionedScalar C1_
 
dimensionedScalar C2_
 
dimensionedScalar C3_
 
dimensionedScalar sigmak_
 
dimensionedScalar sigmaEps_
 
volScalarField k_
 
volScalarField epsilon_
 
- Protected Attributes inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > >
volScalarField nut_
 
- Protected Attributes inherited from RASModel< BasicMomentumTransportModel >
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< BasicMomentumTransportModel >
static autoPtr< RASModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Return a reference to the selected RAS model. More...
 

Detailed Description

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

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 BasicMomentumTransportModel::alphaField alphaField

Definition at line 105 of file continuousGasKEpsilon.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 106 of file continuousGasKEpsilon.H.

◆ transportModel

typedef BasicMomentumTransportModel::transportModel transportModel

Definition at line 107 of file continuousGasKEpsilon.H.

Constructor & Destructor Documentation

◆ continuousGasKEpsilon() [1/2]

continuousGasKEpsilon ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport,
const word type = typeName 
)

Construct from components.

Definition at line 44 of file continuousGasKEpsilon.C.

◆ continuousGasKEpsilon() [2/2]

continuousGasKEpsilon ( const continuousGasKEpsilon< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~continuousGasKEpsilon()

Member Function Documentation

◆ correctNut()

◆ phaseTransferCoeff()

tmp< volScalarField > phaseTransferCoeff ( ) const
protected

Definition at line 225 of file continuousGasKEpsilon.C.

References TimeState::deltaT(), momentumTransportModel::epsilon(), momentumTransportModel::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< BasicMomentumTransportModel >.

Definition at line 248 of file continuousGasKEpsilon.C.

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

Here is the call graph for this function:

◆ epsilonSource()

tmp< fvScalarMatrix > epsilonSource ( ) const
protectedvirtual

Reimplemented from kEpsilon< BasicMomentumTransportModel >.

Definition at line 261 of file continuousGasKEpsilon.C.

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

Here is the call graph for this function:

◆ TypeName()

TypeName ( "continuousGasKEpsilon< BasicMomentumTransportModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Re-read model coefficients if they have changed.

Reimplemented from kEpsilon< BasicMomentumTransportModel >.

Definition at line 100 of file continuousGasKEpsilon.C.

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

Here is the caller graph for this function:

◆ liquidTurbulence()

const momentumTransportModel & liquidTurbulence ( ) const

Return the turbulence model for the liquid phase.

Definition at line 152 of file continuousGasKEpsilon.C.

References IOobject::db(), fluid, phaseModel::fluid(), IOobject::groupName(), objectRegistry::lookupObject(), and phaseModel::name().

Referenced by continuousGasKEpsilon< BasicMomentumTransportModel >::~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< BasicMomentumTransportModel >.

Definition at line 179 of file continuousGasKEpsilon.C.

References IOobject::groupName(), Foam::max(), Foam::min(), and GeometricField< scalar, fvPatchField, volMesh >::New().

Referenced by continuousGasKEpsilon< BasicMomentumTransportModel >::~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 206 of file continuousGasKEpsilon.C.

References virtualMassModel::Cvm(), fluid, phaseModel::fluid(), IOobject::groupName(), GeometricField< scalar, fvPatchField, volMesh >::New(), and phaseModel::rho().

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

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

◆ sigma()

tmp< volSymmTensorField > sigma ( ) const
virtual

Return the stress tensor [m^2/s^2].

Reimplemented from eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 274 of file continuousGasKEpsilon.C.

References Foam::dev(), Foam::fvc::grad(), IOobject::groupName(), Foam::I, k, GeometricField< symmTensor, fvPatchField, volMesh >::New(), and Foam::twoSymm().

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

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

◆ operator=()

void operator= ( const continuousGasKEpsilon< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise assignment.

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

Here is the caller graph for this function:

Member Data Documentation

◆ alphaInversion_

dimensionedScalar alphaInversion_
protected

Definition at line 92 of file continuousGasKEpsilon.H.


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