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


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 turbulenceModel & | liquidTurbulence () const |
| Return the turbulence model for the liquid phase. More... | |
| virtual tmp< volScalarField > | nuEff () const |
| Return the effective viscosity. More... | |
| virtual tmp< volScalarField > | rhoEff () const |
| Return the effective density for the stress. More... | |
| virtual tmp< volSymmTensorField > | R () 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< volScalarField > | DkEff () const |
| Return the effective diffusivity for k. More... | |
| tmp< volScalarField > | DepsilonEff () const |
| Return the effective diffusivity for epsilon. More... | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy. More... | |
| virtual tmp< volScalarField > | epsilon () 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... | |
| volScalarField & | evNut () |
| Return non-const access to the turbulence viscosity. More... | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity. More... | |
| virtual tmp< scalarField > | nut (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< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. More... | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. More... | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (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 dimensionedScalar & | kMin () const |
| Return the lower allowable limit for k (default: SMALL) More... | |
| const dimensionedScalar & | epsilonMin () const |
| Return the lower allowable limit for epsilon (default: SMALL) More... | |
| const dimensionedScalar & | omegaMin () const |
| Return the lower allowable limit for omega (default: SMALL) More... | |
| dimensionedScalar & | kMin () |
| Allow kMin to be changed. More... | |
| dimensionedScalar & | epsilonMin () |
| Allow epsilonMin to be changed. More... | |
| dimensionedScalar & | omegaMin () |
| Allow omegaMin to be changed. More... | |
| virtual const dictionary & | coeffDict () const |
| Const access to the coefficients dictionary. More... | |
| virtual tmp< scalarField > | nuEff (const label patchi) const |
| Return the effective viscosity on patch. More... | |
Protected Member Functions | |
| virtual void | correctNut () |
| tmp< volScalarField > | phaseTransferCoeff () const |
| virtual tmp< fvScalarMatrix > | kSource () const |
| virtual tmp< fvScalarMatrix > | epsilonSource () 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< RASModel > | New (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... | |
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 78 of file continuousGasKEpsilon.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 115 of file continuousGasKEpsilon.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 116 of file continuousGasKEpsilon.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 117 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 42 of file continuousGasKEpsilon.C.
|
inlinevirtual |
Destructor.
Definition at line 141 of file continuousGasKEpsilon.H.
References continuousGasKEpsilon< BasicTurbulenceModel >::liquidTurbulence(), continuousGasKEpsilon< BasicTurbulenceModel >::nuEff(), continuousGasKEpsilon< BasicTurbulenceModel >::R(), continuousGasKEpsilon< BasicTurbulenceModel >::read(), and continuousGasKEpsilon< BasicTurbulenceModel >::rhoEff().

|
protectedvirtual |
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 116 of file continuousGasKEpsilon.C.
References optionList::correct(), kEpsilon< BasicTurbulenceModel >::correctNut(), virtualMassModel::Cvm(), turbulenceModel::epsilon(), Foam::exp(), fluid, turbulenceModel::k(), Foam::min(), options::New(), turbulenceModel::nut(), twoPhaseSystem::otherPhase(), Foam::sqr(), and twoPhaseSystem::virtualMass().

|
protected |
Definition at line 222 of file continuousGasKEpsilon.C.
References TimeState::deltaT(), turbulenceModel::epsilon(), turbulenceModel::k(), Foam::max(), Foam::min(), and IOobject::time().

|
protectedvirtual |
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 245 of file continuousGasKEpsilon.C.
References turbulenceModel::k(), and Foam::fvm::Sp().

|
protectedvirtual |
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 258 of file continuousGasKEpsilon.C.
References turbulenceModel::epsilon(), and Foam::fvm::Sp().

| TypeName | ( | "continuousGasKEpsilon< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented from kEpsilon< BasicTurbulenceModel >.
Definition at line 100 of file continuousGasKEpsilon.C.
Referenced by continuousGasKEpsilon< BasicTurbulenceModel >::~continuousGasKEpsilon().

| const turbulenceModel & liquidTurbulence | ( | ) | const |
Return the turbulence model for the liquid phase.
Definition at line 145 of file continuousGasKEpsilon.C.
References IOobject::db(), fluid, IOobject::groupName(), objectRegistry::lookupObject(), twoPhaseSystem::otherPhase(), and turbulenceModel::propertiesName.
Referenced by continuousGasKEpsilon< BasicTurbulenceModel >::~continuousGasKEpsilon().


|
virtual |
Return the effective viscosity.
Reimplemented from RASModel< BasicTurbulenceModel >.
Definition at line 173 of file continuousGasKEpsilon.C.
References IOobject::groupName(), Foam::max(), Foam::min(), and nu.
Referenced by continuousGasKEpsilon< BasicTurbulenceModel >::~continuousGasKEpsilon().


|
virtual |
Return the effective density for the stress.
Definition at line 203 of file continuousGasKEpsilon.C.
References virtualMassModel::Cvm(), fluid, IOobject::groupName(), twoPhaseSystem::otherPhase(), and twoPhaseSystem::virtualMass().
Referenced by continuousGasKEpsilon< BasicTurbulenceModel >::~continuousGasKEpsilon().


|
virtual |
Return the Reynolds stress tensor.
Reimplemented from eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 271 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().


|
protected |
Definition at line 102 of file continuousGasKEpsilon.H.
1.8.11