Loading...
Searching...
No Matches
kOmegaSSTSato< BasicTurbulenceModel > Class Template Reference

Implementation of the k-omega-SST turbulence model for dispersed bubbly flows with Sato (1981) bubble induced turbulent viscosity model. More...

#include <kOmegaSSTSato.H>

Inheritance diagram for kOmegaSSTSato< BasicTurbulenceModel >:
Collaboration diagram for kOmegaSSTSato< BasicTurbulenceModel >:

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 eddyViscosity< RASModel< BasicTurbulenceModel > >::alphaField alphaField
typedef eddyViscosity< RASModel< BasicTurbulenceModel > >::rhoField rhoField
typedef eddyViscosity< RASModel< 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 ("kOmegaSSTSato")
 Runtime type information.
 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.
virtual ~kOmegaSSTSato ()=default
 Destructor.
virtual bool read ()
 Read model coefficients if they have changed.
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity.
Public Member Functions inherited from kOmegaSST< BasicTurbulenceModel >
 TypeName ("kOmegaSST")
 Runtime type information.
 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.
virtual ~kOmegaSST ()=default
 Destructor.
Public Member Functions inherited from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >
 kOmegaSSTBase (const kOmegaSSTBase &)=delete
 No copy construct.
void operator= (const kOmegaSSTBase &)=delete
 No copy assignment.
virtual ~kOmegaSSTBase ()=default
 Destructor.
tmp< volScalarFieldDkEff (const volScalarField &F1) const
 Return the effective diffusivity for k.
tmp< volScalarFieldDomegaEff (const volScalarField &F1) const
 Return the effective diffusivity for omega.
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy.
virtual tmp< volScalarFieldomega () const
 Return the turbulence kinetic energy dissipation rate.
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< volScalarFieldnut () const
 Return the turbulence viscosity.
virtual tmp< volSymmTensorFieldR () 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< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field.
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation.
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation.

Protected Member Functions

virtual void correctNut (const volScalarField &S2)
Protected Member Functions inherited from kOmegaSST< BasicTurbulenceModel >
virtual void correctNut ()
Protected Member Functions inherited from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >
void setDecayControl (const dictionary &dict)
 Set decay control with kInf and omegaInf.
virtual tmp< volScalarFieldF1 (const volScalarField &CDkOmega) const
virtual tmp< volScalarFieldF2 () const
virtual tmp< volScalarFieldF3 () const
virtual tmp< volScalarFieldF23 () const
tmp< volScalarFieldblend (const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
 Return the blended field.
tmp< volScalarFieldalphaK (const volScalarField &F1) const
tmp< volScalarFieldalphaOmega (const volScalarField &F1) const
tmp< volScalarField::Internalbeta (const volScalarField::Internal &F1) const
tmp< volScalarField::Internalgamma (const volScalarField::Internal &F1) const
virtual tmp< volScalarFieldS2 (const volTensorField &gradU) const
 Return square of strain rate.
virtual tmp< volScalarField::InternalPk (const volScalarField::Internal &G) const
 Return k production rate.
virtual tmp< volScalarField::InternalepsilonByk (const volScalarField &F1, const volTensorField &gradU) const
 Return epsilon/k which for standard RAS is betaStar*omega.
virtual tmp< volScalarField::InternalGbyNu0 (const volTensorField &gradU, const volScalarField &S2) const
 Return (G/nu)_0.
virtual tmp< volScalarField::InternalGbyNu (const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
 Return G/nu.
virtual tmp< fvScalarMatrixkSource () const
virtual tmp< fvScalarMatrixomegaSource () const
virtual tmp< fvScalarMatrixQsas (const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const

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.
const volScalarFieldy_
 Wall distance.
volScalarField k_
 Turbulent kinetic energy field [m^2/s^2].
volScalarField omega_
 Specific dissipation rate field [1/s].
Switch decayControl_
 Flag to include the decay control.
dimensionedScalar kInf_
dimensionedScalar omegaInf_
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_

Detailed Description

template<class BasicTurbulenceModel>
class Foam::RASModels::kOmegaSSTSato< BasicTurbulenceModel >

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;
    }
Source files
Source files

Definition at line 121 of file kOmegaSSTSato.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 172 of file kOmegaSSTSato.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 173 of file kOmegaSSTSato.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 174 of file kOmegaSSTSato.H.

Constructor & Destructor Documentation

◆ kOmegaSSTSato()

template<class BasicTurbulenceModel>
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 36 of file kOmegaSSTSato.C.

References alpha, Cmub_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, phi, rho, and U.

◆ ~kOmegaSSTSato()

template<class BasicTurbulenceModel>
virtual ~kOmegaSSTSato ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "kOmegaSSTSato< BasicTurbulenceModel >" )

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

Read model coefficients if they have changed.

Reimplemented from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >.

Definition at line 82 of file kOmegaSSTSato.C.

References Cmub_, and linearViscousStress< BasicTurbulenceModel >::read().

Here is the call graph for this function:

◆ correct()

template<class BasicTurbulenceModel>
void correct ( )
virtual

Solve the turbulence equations and correct the turbulence viscosity.

Reimplemented from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >.

Definition at line 160 of file kOmegaSSTSato.C.

References linearViscousStress< BasicTurbulenceModel >::correct().

Here is the call graph for this function:

Member Data Documentation

◆ Cmub_

template<class BasicTurbulenceModel>
dimensionedScalar Cmub_
protected

Definition at line 162 of file kOmegaSSTSato.H.

Referenced by correctNut(), kOmegaSSTSato(), and read().


The documentation for this class was generated from the following files: