One-equation SGS model for the gas-phase in a two-phase system supporting phase-inversion. More...
#include <continuousGasKEqn.H>


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< volScalarField > | k () const |
| Return SGS kinetic energy. | |
| tmp< volScalarField > | DkEff () 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< volScalarField > | nut () const |
| Return the turbulence viscosity. | |
| virtual tmp< volSymmTensorField > | R () 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< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
| Return the effective stress tensor based on a given velocity field. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. | |
Protected Member Functions | |
| const turbulenceModel & | liquidTurbulence () const |
| Return the turbulence model for the liquid phase. | |
| tmp< volScalarField > | phaseTransferCoeff () const |
| virtual tmp< fvScalarMatrix > | kSource () 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_ |
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;
}
Definition at line 72 of file continuousGasKEqn.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 116 of file continuousGasKEqn.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 117 of file continuousGasKEqn.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 118 of file continuousGasKEqn.H.
| 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.

|
virtualdefault |
Destructor.
|
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().


|
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().


|
protectedvirtual |
Reimplemented from kEqn< BasicTurbulenceModel >.
Definition at line 147 of file continuousGasKEqn.C.
References kEqn< BasicTurbulenceModel >::k_, liquidTurbulence(), phaseTransferCoeff(), and Foam::fvm::Sp().

| TypeName | ( | "continuousGasKEqn< BasicTurbulenceModel >" | ) |
Runtime type information.
References alpha, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, phi, turbulenceModel::propertiesName, rho, and U.
|
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().

|
protected |
Definition at line 100 of file continuousGasKEqn.H.
Referenced by continuousGasKEqn(), phaseTransferCoeff(), and read().