Implementation of the k-omega-SST turbulence model for dispersed bubbly flows with Sato (1981) bubble induced turbulent viscosity model. More...
Public Types | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
Public Types inherited from kOmegaSST< BasicTurbulenceModel > | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
Public Types inherited from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > > | |
typedef BasicEddyViscosityModel::alphaField | alphaField |
typedef BasicEddyViscosityModel::rhoField | rhoField |
typedef BasicEddyViscosityModel::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 ("kOmegaSSTSato") | |
Runtime type information. More... | |
kOmegaSSTSato (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 | ~kOmegaSSTSato ()=default |
Destructor. More... | |
virtual bool | read () |
Read model coefficients if they have changed. More... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
Public Member Functions inherited from kOmegaSST< BasicTurbulenceModel > | |
TypeName ("kOmegaSST") | |
Runtime type information. More... | |
kOmegaSST (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 | ~kOmegaSST ()=default |
Destructor. More... | |
Public Member Functions inherited from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > > | |
kOmegaSSTBase (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 | ~kOmegaSSTBase ()=default |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. More... | |
tmp< volScalarField > | DkEff (const volScalarField &F1) const |
Return the effective diffusivity for k. More... | |
tmp< volScalarField > | DomegaEff (const volScalarField &F1) const |
Return the effective diffusivity for omega. More... | |
virtual tmp< volScalarField > | k () const |
Return the turbulence kinetic energy. More... | |
virtual tmp< volScalarField > | omega () 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 ()=default |
Destructor. More... | |
virtual bool | read ()=0 |
Re-read model coefficients if they have changed. 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< volScalarField > | k () const=0 |
Return the turbulence kinetic energy. More... | |
virtual tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
virtual void | correct ()=0 |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
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. More... | |
virtual | ~linearViscousStress ()=default |
Destructor. More... | |
virtual bool | read ()=0 |
Re-read model coefficients if they have changed. More... | |
virtual tmp< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
Return the effective stress tensor based on a given velocity field. 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... | |
virtual void | correct ()=0 |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
Protected Attributes | |
dimensionedScalar | Cmub_ |
Protected Attributes inherited from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > > | |
dimensionedScalar | alphaK1_ |
dimensionedScalar | alphaK2_ |
dimensionedScalar | alphaOmega1_ |
dimensionedScalar | alphaOmega2_ |
dimensionedScalar | gamma1_ |
dimensionedScalar | gamma2_ |
dimensionedScalar | beta1_ |
dimensionedScalar | beta2_ |
dimensionedScalar | betaStar_ |
dimensionedScalar | a1_ |
dimensionedScalar | b1_ |
dimensionedScalar | c1_ |
Switch | F3_ |
Flag to include the F3 term. More... | |
const volScalarField & | y_ |
Wall distance. More... | |
volScalarField | k_ |
volScalarField | omega_ |
Switch | decayControl_ |
Flag to include the decay control. More... | |
dimensionedScalar | kInf_ |
dimensionedScalar | omegaInf_ |
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
volScalarField | nut_ |
Implementation of the k-omega-SST turbulence model for dispersed bubbly flows with Sato (1981) bubble induced turbulent viscosity model.
Bubble induced turbulent viscosity model described in:
Sato, Y., Sadatomi, M. "Momentum and heat transfer in two-phase bubble flow - I, Theory" International Journal of Multiphase FLow 7, pp. 167-177, 1981.
Turbulence model described in:
Menter, F., Esch, T. "Elements of Industrial Heat Transfer Prediction" 16th Brazilian Congress of Mechanical Engineering (COBEM), Nov. 2001
with the addition of the optional F3 term for rough walls from
Hellsten, A. "Some Improvements in Menter’s k-omega-SST turbulence model" 29th AIAA Fluid Dynamics Conference, AIAA-98-2554, June 1998.
Note that this implementation is written in terms of alpha diffusion coefficients rather than the more traditional sigma (alpha = 1/sigma) so that the blending can be applied to all coefficuients in a consistent manner. The paper suggests that sigma is blended but this would not be consistent with the blending of the k-epsilon and k-omega models.
Also note that the error in the last term of equation (2) relating to sigma has been corrected.
Wall-functions are applied in this implementation by using equations (14) to specify the near-wall omega as appropriate.
The blending functions (15) and (16) are not currently used because of the uncertainty in their origin, range of applicability and that as y+ becomes sufficiently small blending u_tau in this manner is clearly nonsense.
The default model coefficients correspond to the following:
kOmegaSSTCoeffs { alphaK1 0.85034; alphaK2 1.0; alphaOmega1 0.5; alphaOmega2 0.85616; Prt 1.0; // only for compressible beta1 0.075; beta2 0.0828; betaStar 0.09; gamma1 0.5532; gamma2 0.4403; a1 0.31; b1 1.0; c1 10.0; F3 no; Cmub 0.6; }
Definition at line 123 of file kOmegaSSTSato.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 168 of file kOmegaSSTSato.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 169 of file kOmegaSSTSato.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 170 of file kOmegaSSTSato.H.
kOmegaSSTSato | ( | 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 43 of file kOmegaSSTSato.C.
References Foam::type().
|
virtualdefault |
Destructor.
|
protectedvirtual |
Reimplemented from kOmegaSST< BasicTurbulenceModel >.
Definition at line 135 of file kOmegaSSTSato.C.
References TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >::alpha(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::exp(), Foam::mag(), Foam::max(), Time::New(), nu, Foam::pow(), Foam::sqr(), Foam::sqrt(), TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >::transport(), turbulenceModel::U(), and yPlus.
TypeName | ( | "kOmegaSSTSato< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Reimplemented from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >.
Definition at line 89 of file kOmegaSSTSato.C.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Reimplemented from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >.
Definition at line 167 of file kOmegaSSTSato.C.
References kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct().
|
protected |
Definition at line 158 of file kOmegaSSTSato.H.