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

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

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

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

Public Member Functions

 TypeName ("continuousGasKEqn")
 Runtime type information. More...
 
 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. More...
 
virtual ~continuousGasKEqn ()
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
- Public Member Functions inherited from kEqn< BasicTurbulenceModel >
 TypeName ("kEqn")
 Runtime type information. More...
 
 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. More...
 
virtual ~kEqn ()
 Destructor. More...
 
virtual tmp< volScalarFieldk () const
 Return SGS kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return sub-grid disipation rate. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
virtual void correct ()
 Correct eddy-Viscosity and related properties. More...
 
- 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. More...
 
virtual ~LESeddyViscosity ()
 Destructor. More...
 
- 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. 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 tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
- Public Member Functions inherited from linearViscousStress< LESModel< 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 LESModel< BasicTurbulenceModel >
 TypeName ("LES")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, LESModel, 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))
 
 LESModel (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 ~LESModel ()
 Destructor. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: small) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
const volScalarFielddelta () const
 Access function to filter width. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 

Protected Member Functions

const turbulenceModelliquidTurbulence () const
 Return the turbulence model for the liquid phase. More...
 
tmp< volScalarFieldphaseTransferCoeff () const
 
virtual tmp< fvScalarMatrixkSource () const
 
- Protected Member Functions inherited from kEqn< BasicTurbulenceModel >
virtual void correctNut ()
 
- Protected Member Functions inherited from LESModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 

Protected Attributes

dimensionedScalar alphaInversion_
 
- Protected Attributes inherited from kEqn< BasicTurbulenceModel >
volScalarField k_
 
dimensionedScalar Ck_
 
- Protected Attributes inherited from LESeddyViscosity< BasicTurbulenceModel >
dimensionedScalar Ce_
 
- Protected Attributes inherited from eddyViscosity< LESModel< BasicTurbulenceModel > >
volScalarField nut_
 
- Protected Attributes inherited from LESModel< BasicTurbulenceModel >
dictionary LESDict_
 LES 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...
 
autoPtr< Foam::LESdeltadelta_
 Run-time selectable delta model. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from LESModel< BasicTurbulenceModel >
static autoPtr< LESModelNew (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 LES model. More...
 

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 70 of file continuousGasKEqn.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 106 of file continuousGasKEqn.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 107 of file continuousGasKEqn.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 108 of file continuousGasKEqn.H.

Constructor & Destructor Documentation

◆ continuousGasKEqn()

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

◆ ~continuousGasKEqn()

virtual ~continuousGasKEqn ( )
inlinevirtual

Destructor.

Definition at line 132 of file continuousGasKEqn.H.

References continuousGasKEqn< BasicTurbulenceModel >::read().

Here is the call graph for this function:

Member Function Documentation

◆ liquidTurbulence()

const turbulenceModel & liquidTurbulence ( ) const
protected

Return the turbulence model for the liquid phase.

Definition at line 102 of file continuousGasKEqn.C.

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

Here is the call graph for this function:

◆ phaseTransferCoeff()

tmp< volScalarField > phaseTransferCoeff ( ) const
protected

Definition at line 130 of file continuousGasKEqn.C.

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

Here is the call graph for this function:

◆ kSource()

tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Reimplemented from kEqn< BasicTurbulenceModel >.

Definition at line 153 of file continuousGasKEqn.C.

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

Here is the call graph for this function:

◆ TypeName()

TypeName ( "continuousGasKEqn< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Read model coefficients if they have changed.

Reimplemented from kEqn< BasicTurbulenceModel >.

Definition at line 85 of file continuousGasKEqn.C.

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

Here is the caller graph for this function:

Member Data Documentation

◆ alphaInversion_

dimensionedScalar alphaInversion_
protected

Definition at line 92 of file continuousGasKEqn.H.


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