k-epsilon model for the gas-phase in a two-phase system supporting phase-inversion. More...
#include <continuousGasKEpsilon.H>


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 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 ("continuousGasKEpsilon") | |
| Runtime type information. | |
| 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. | |
| virtual | ~continuousGasKEpsilon ()=default |
| Destructor. | |
| virtual bool | read () |
| Re-read model coefficients if they have changed. | |
| const turbulenceModel & | liquidTurbulence () const |
| Return the turbulence model for the liquid phase. | |
| virtual tmp< volScalarField > | nuEff () const |
| Return the effective viscosity. | |
| virtual tmp< volScalarField > | rhoEff () const |
| Return the effective density for the stress. | |
| virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. | |
| Public Member Functions inherited from kEpsilon< BasicTurbulenceModel > | |
| TypeName ("kEpsilon") | |
| Runtime type information. | |
| 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. | |
| virtual | ~kEpsilon ()=default |
| Destructor. | |
| tmp< volScalarField > | DkEff () const |
| Return the effective diffusivity for k. | |
| tmp< volScalarField > | DepsilonEff () const |
| Return the effective diffusivity for epsilon. | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy. | |
| virtual tmp< volScalarField > | epsilon () const |
| Return the turbulence kinetic energy dissipation rate. | |
| virtual void | correct () |
| Solve the turbulence equations and correct the turbulence viscosity. | |
| 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. | |
| virtual | ~eddyViscosity ()=default |
| Destructor. | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity. | |
| 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 | |
| virtual void | correctNut () |
| tmp< volScalarField > | phaseTransferCoeff () const |
| virtual tmp< fvScalarMatrix > | kSource () const |
| virtual tmp< fvScalarMatrix > | epsilonSource () const |
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_ |
| bool | twoLayerTreatment_ |
| dimensionedScalar | ReyStar_ |
| dimensionedScalar | ReyFactor_ |
| std::unique_ptr< volScalarField::Internal > | lEpsPtr_ |
| std::unique_ptr< volScalarField > | lambdaEpsPtr_ |
| Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| volScalarField | nut_ |
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;
}
Definition at line 77 of file continuousGasKEpsilon.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 120 of file continuousGasKEpsilon.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 121 of file continuousGasKEpsilon.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 122 of file continuousGasKEpsilon.H.
| 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 38 of file continuousGasKEpsilon.C.
References alpha, alphaInversion_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_, phi, rho, timeName, and U.
|
virtualdefault |
Destructor.
|
protectedvirtual |
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 111 of file continuousGasKEpsilon.C.
References optionList::correct(), kEpsilon< BasicTurbulenceModel >::correctNut(), virtualMassModel::Cvm(), Foam::exp(), fluid, liquidTurbulence(), Foam::min(), options::New(), Foam::refCast(), liquid::rho(), and Foam::sqr().

|
protected |
Definition at line 223 of file continuousGasKEpsilon.C.
References alpha, alphaInversion_, liquidTurbulence(), Foam::max(), Foam::min(), rho, and U.
Referenced by epsilonSource(), and kSource().


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

|
protectedvirtual |
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 259 of file continuousGasKEpsilon.C.
References kEpsilon< BasicTurbulenceModel >::epsilon_, liquidTurbulence(), phaseTransferCoeff(), and Foam::fvm::Sp().

| TypeName | ( | "continuousGasKEpsilon< BasicTurbulenceModel >" | ) |
Runtime type information.
References alpha, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, phi, turbulenceModel::propertiesName, rho, and U.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 97 of file continuousGasKEpsilon.C.
References alphaInversion_, and kEpsilon< BasicTurbulenceModel >::read().

| const turbulenceModel & liquidTurbulence | ( | ) | const |
Return the turbulence model for the liquid phase.
Definition at line 143 of file continuousGasKEpsilon.C.
References fluid, IOobject::groupName(), turbulenceModel::propertiesName, Foam::refCast(), and U.
Referenced by correctNut(), epsilonSource(), kSource(), and phaseTransferCoeff().


|
virtual |
Return the effective viscosity.
Definition at line 171 of file continuousGasKEpsilon.C.
References alphaInversion_, IOobject::groupName(), Foam::max(), Foam::min(), nu, rho, and rhoEff().

|
virtual |
Return the effective density for the stress.
Definition at line 201 of file continuousGasKEpsilon.C.
References virtualMassModel::Cvm(), fluid, IOobject::groupName(), Foam::refCast(), and liquid::rho().
Referenced by nuEff().


|
virtual |
Return the Reynolds stress tensor.
Reimplemented from eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 272 of file continuousGasKEpsilon.C.
References Foam::devTwoSymm(), Foam::fvc::grad(), IOobject::groupName(), Foam::I, k, IOobjectOption::NO_READ, and IOobjectOption::NO_WRITE.

|
protected |
Definition at line 107 of file continuousGasKEpsilon.H.
Referenced by continuousGasKEpsilon(), nuEff(), phaseTransferCoeff(), and read().