Mixture k-epsilon turbulence model for two-phase gas-liquid systems. More...
Public Types | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
![]() | |
typedef RASModel< BasicMomentumTransportModel > ::alphaField | alphaField |
typedef RASModel< BasicMomentumTransportModel > ::rhoField | rhoField |
![]() | |
typedef RASModel< BasicMomentumTransportModel > ::alphaField | alphaField |
typedef RASModel< BasicMomentumTransportModel > ::rhoField | rhoField |
![]() | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Member Functions | |
TypeName ("mixtureKEpsilon") | |
Runtime type information. More... | |
mixtureKEpsilon (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName) | |
Construct from components. More... | |
mixtureKEpsilon (const mixtureKEpsilon &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~mixtureKEpsilon () |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. 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 tmp< volScalarField > | omega () const |
Return the turbulence specific dissipation rate. More... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
void | operator= (const mixtureKEpsilon &)=delete |
Disallow default bitwise assignment. More... | |
![]() | |
eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
Construct from components. More... | |
virtual | ~eddyViscosity () |
Destructor. 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 tmp< volSymmTensorField > | sigma () const |
Return the Reynolds stress tensor [m^2/s^2]. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
![]() | |
linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
Construct from components. More... | |
virtual | ~linearViscousStress () |
Destructor. More... | |
virtual tmp< volSymmTensorField > | devTau () const |
Return the effective stress tensor. More... | |
virtual tmp< fvVectorMatrix > | divDevTau (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevTau (const volScalarField &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
![]() | |
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 viscosity &viscosity),(alpha, rho, U, alphaRhoPhi, phi, viscosity)) | |
RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
Construct from components. More... | |
RASModel (const RASModel &)=delete | |
Disallow default bitwise copy construction. 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< volScalarField > | nu () const |
Return the laminar viscosity. More... | |
virtual tmp< scalarField > | nu (const label patchi) const |
Return the laminar viscosity on patchi. More... | |
virtual tmp< volScalarField > | nuEff () const |
Return the effective viscosity. More... | |
virtual tmp< scalarField > | nuEff (const label patchi) const |
Return the effective viscosity on patch. More... | |
void | operator= (const RASModel &)=delete |
Disallow default bitwise assignment. More... | |
Protected Member Functions | |
void | initMixtureFields () |
virtual void | correctNut () |
tmp< volScalarField > | Ct2 () const |
tmp< volScalarField > | rholEff () const |
tmp< volScalarField > | rhogEff () const |
tmp< volScalarField > | rhom () const |
tmp< volScalarField > | mix (const volScalarField &fc, const volScalarField &fd) const |
tmp< volScalarField > | mixU (const volScalarField &fc, const volScalarField &fd) const |
tmp< surfaceScalarField > | mixFlux (const surfaceScalarField &fc, const surfaceScalarField &fd) const |
tmp< volScalarField > | bubbleG () const |
virtual tmp< fvScalarMatrix > | kSource () const |
virtual tmp< fvScalarMatrix > | epsilonSource () const |
tmp< volScalarField > | DkEff (const volScalarField &nutm) const |
Return the effective diffusivity for k. More... | |
tmp< volScalarField > | DepsilonEff (const volScalarField &nutm) const |
Return the effective diffusivity for epsilon. More... | |
![]() | |
virtual void | printCoeffs (const word &type) |
Print model coefficients. More... | |
Additional Inherited Members | |
![]() | |
static autoPtr< RASModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) |
Return a reference to the selected RAS model. More... | |
Mixture k-epsilon turbulence model for two-phase gas-liquid systems.
The basic structure of the model is based on:
Behzadi, A., Issa, R. I., & Rusche, H. (2004). Modelling of dispersed bubble and droplet flow at high phase fractions. Chemical Engineering Science, 59(4), 759-770.
but an effective density for the gas is used in the averaging and an alternative model for bubble-generated turbulence from:
Lahey Jr, R. T. (2005). The simulation of multidimensional multiphase flows. Nuclear Engineering and Design, 235(10), 1043-1060.
The default model coefficients are
mixtureKEpsilonCoeffs { Cmu 0.09; C1 1.44; C2 1.92; C3 C2; Cp 0.25; // Bubble generated turbulence alphap 0.3; // Gas phase fraction below which // bubble generated turbulence is included sigmak 1.0; sigmaEps 1.3; }
Definition at line 83 of file mixtureKEpsilon.H.
typedef BasicMomentumTransportModel::alphaField alphaField |
Definition at line 183 of file mixtureKEpsilon.H.
typedef BasicMomentumTransportModel::rhoField rhoField |
Definition at line 184 of file mixtureKEpsilon.H.
mixtureKEpsilon | ( | const alphaField & | alpha, |
const rhoField & | rho, | ||
const volVectorField & | U, | ||
const surfaceScalarField & | alphaRhoPhi, | ||
const surfaceScalarField & | phi, | ||
const viscosity & | viscosity, | ||
const word & | type = typeName |
||
) |
Construct from components.
Definition at line 49 of file mixtureKEpsilon.C.
References IOobject::AUTO_WRITE, Foam::bound(), Foam::dimDensity, Foam::dimless, phaseModel::index(), IOobject::MUST_READ, and IOobject::READ_IF_PRESENT.
|
delete |
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 210 of file mixtureKEpsilon.H.
References mixtureKEpsilon< BasicMomentumTransportModel >::read().
|
protected |
Definition at line 253 of file mixtureKEpsilon.C.
|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 291 of file mixtureKEpsilon.C.
References IOobject::db(), fluid, phaseModel::fluid(), IOobject::groupName(), objectRegistry::lookupObject(), phaseModel::name(), dictionary::New(), and Foam::sqr().
|
protected |
Definition at line 334 of file mixtureKEpsilon.C.
References drag, Foam::exp(), fluid, phaseModel::fluid(), dispersedDragModel::K(), Foam::mag(), phaseModel::rho(), Foam::sqr(), Foam::sqrt(), and U.
|
protected |
Definition at line 365 of file mixtureKEpsilon.C.
References phaseModel::rho().
|
protected |
Definition at line 374 of file mixtureKEpsilon.C.
References dispersedVirtualMassModel::Cvm(), fluid, phaseModel::fluid(), and phaseModel::rho().
|
protected |
Definition at line 390 of file mixtureKEpsilon.C.
References alphal, and mixtureKEpsilon< BasicMomentumTransportModel >::mix().
|
protected |
Definition at line 401 of file mixtureKEpsilon.C.
References alphal, and mixtureKEpsilon< BasicMomentumTransportModel >::mixU().
Referenced by mixtureKEpsilon< BasicMomentumTransportModel >::rhom().
|
protected |
Definition at line 415 of file mixtureKEpsilon.C.
References alphal, and mixtureKEpsilon< BasicMomentumTransportModel >::mixFlux().
Referenced by mixtureKEpsilon< BasicMomentumTransportModel >::mix().
|
protected |
Definition at line 431 of file mixtureKEpsilon.C.
References alphal, and Foam::fvc::interpolate().
Referenced by mixtureKEpsilon< BasicMomentumTransportModel >::mixU().
|
protected |
Definition at line 453 of file mixtureKEpsilon.C.
References dispersedDragModel::CdRe(), phaseModel::d(), drag, fluid, phaseModel::fluid(), Foam::mag(), fluidThermo::nu(), Foam::pos(), Foam::pow(), Foam::pow3(), phaseModel::rho(), phaseModel::thermo(), and U.
|
protectedvirtual |
Definition at line 494 of file mixtureKEpsilon.C.
References Foam::fvm::Su().
|
protectedvirtual |
Definition at line 502 of file mixtureKEpsilon.C.
References Foam::fvm::Su().
|
inlineprotected |
Return the effective diffusivity for k.
Definition at line 162 of file mixtureKEpsilon.H.
References GeometricField< scalar, fvPatchField, volMesh >::New().
|
inlineprotected |
Return the effective diffusivity for epsilon.
Definition at line 172 of file mixtureKEpsilon.H.
References GeometricField< scalar, fvPatchField, volMesh >::New().
TypeName | ( | "mixtureKEpsilon< BasicMomentumTransportModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 269 of file mixtureKEpsilon.C.
References Foam::read().
Referenced by mixtureKEpsilon< BasicMomentumTransportModel >::~mixtureKEpsilon().
|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 220 of file mixtureKEpsilon.H.
References mixtureKEpsilon< BasicMomentumTransportModel >::k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 226 of file mixtureKEpsilon.H.
References mixtureKEpsilon< BasicMomentumTransportModel >::epsilon_.
|
inlinevirtual |
Return the turbulence specific dissipation rate.
Definition at line 232 of file mixtureKEpsilon.H.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), mixtureKEpsilon< BasicMomentumTransportModel >::correct(), GeometricField< scalar, fvPatchField, volMesh >::New(), mixtureKEpsilon< BasicMomentumTransportModel >::operator=(), and GeometricBoundaryField< Type, PatchField, GeoMesh >::types().
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 509 of file mixtureKEpsilon.C.
References Foam::fvc::absolute(), alphal, Foam::bound(), GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), tmp< T >::clear(), fvConstraints::constrain(), correct, GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), Foam::fvc::div(), mixtureKEpsilon< BasicMomentumTransportModel >::epsilon_, fvConstraints, fvModels, Foam::fvc::grad(), phaseModel::index(), mixtureKEpsilon< BasicMomentumTransportModel >::k_, Foam::fvm::laplacian(), dictionary::New(), RASModel< BasicMomentumTransportModel >::nu(), eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_, phi, phig(), tmp< T >::ref(), Foam::solve(), Foam::fvm::Sp(), Foam::fvm::SuSp(), Foam::twoSymm(), and GeometricBoundaryField< Type, PatchField, GeoMesh >::updateCoeffs().
Referenced by mixtureKEpsilon< BasicMomentumTransportModel >::omega().
|
delete |
Disallow default bitwise assignment.
Referenced by mixtureKEpsilon< BasicMomentumTransportModel >::omega().
|
protected |
Definition at line 105 of file mixtureKEpsilon.H.
|
protected |
Definition at line 106 of file mixtureKEpsilon.H.
|
protected |
Definition at line 107 of file mixtureKEpsilon.H.
|
protected |
Definition at line 108 of file mixtureKEpsilon.H.
|
protected |
Definition at line 109 of file mixtureKEpsilon.H.
|
protected |
Definition at line 110 of file mixtureKEpsilon.H.
|
protected |
Definition at line 111 of file mixtureKEpsilon.H.
|
protected |
Definition at line 112 of file mixtureKEpsilon.H.
|
protected |
Definition at line 116 of file mixtureKEpsilon.H.
Referenced by mixtureKEpsilon< BasicMomentumTransportModel >::correct(), and mixtureKEpsilon< BasicMomentumTransportModel >::k().
|
protected |
Definition at line 117 of file mixtureKEpsilon.H.
Referenced by mixtureKEpsilon< BasicMomentumTransportModel >::correct(), and mixtureKEpsilon< BasicMomentumTransportModel >::epsilon().
|
protected |
Definition at line 121 of file mixtureKEpsilon.H.
|
protected |
Definition at line 122 of file mixtureKEpsilon.H.
|
protected |
Definition at line 123 of file mixtureKEpsilon.H.
|
protected |
Definition at line 124 of file mixtureKEpsilon.H.