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...

Inheritance diagram for kOmegaSSTSato< BasicTurbulenceModel >:
[legend]
Collaboration diagram for kOmegaSSTSato< BasicTurbulenceModel >:
[legend]

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 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 ("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...
 
tmp< volScalarFieldDkEff (const volScalarField &F1) const
 Return the effective diffusivity for k. More...
 
tmp< volScalarFieldDomegaEff (const volScalarField &F1) const
 Return the effective diffusivity for omega. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy. More...
 
virtual tmp< volScalarFieldomega () const
 Return the turbulence kinetic energy dissipation rate. 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 tmp< volScalarFieldnut () const
 Return the turbulence viscosity. More...
 
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. 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 ()=default
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (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 ()=default
 Destructor. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: SMALL) More...
 
const dimensionedScalarepsilonMin () const
 Return the lower allowable limit for epsilon (default: SMALL) More...
 
const dimensionedScalaromegaMin () const
 Return the lower allowable limit for omega (default: SMALL) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
dimensionedScalarepsilonMin ()
 Allow epsilonMin to be changed. More...
 
dimensionedScalaromegaMin ()
 Allow omegaMin to be changed. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 

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)
 
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
 
tmp< volScalarField::Internalblend (const volScalarField::Internal &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
 
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< volScalarField::InternalPk (const volScalarField::Internal &G) const
 Return k production rate. More...
 
virtual tmp< volScalarField::InternalepsilonByk (const volScalarField &F1, const volTensorField &gradU) const
 Return epsilon/k which for standard RAS is betaStar*omega. More...
 
virtual tmp< volScalarField::InternalGbyNu (const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
 Return G/nu. More...
 
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 Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 
 RASModel (const RASModel &)=delete
 No copy construct. More...
 
void operator= (const RASModel &)=delete
 No copy assignment. 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 volScalarFieldy_
 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_
 
- 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< RASModelNew (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...
 

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 123 of file kOmegaSSTSato.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 168 of file kOmegaSSTSato.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 169 of file kOmegaSSTSato.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 170 of file kOmegaSSTSato.H.

Constructor & Destructor Documentation

◆ kOmegaSSTSato()

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 44 of file kOmegaSSTSato.C.

References Foam::type().

Here is the call graph for this function:

◆ ~kOmegaSSTSato()

virtual ~kOmegaSSTSato ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

void correctNut ( const volScalarField S2)
protectedvirtual

Reimplemented from kOmegaSST< BasicTurbulenceModel >.

Definition at line 136 of file kOmegaSSTSato.C.

References optionList::correct(), Foam::exp(), Foam::mag(), Foam::max(), options::New(), nu, Foam::pow(), Foam::sqr(), Foam::sqrt(), and yPlus.

Here is the call graph for this function:

◆ TypeName()

TypeName ( "kOmegaSSTSato< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Read model coefficients if they have changed.

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

Definition at line 89 of file kOmegaSSTSato.C.

◆ correct()

void correct ( )
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 kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct().

Here is the call graph for this function:

Member Data Documentation

◆ Cmub_

dimensionedScalar Cmub_
protected

Definition at line 158 of file kOmegaSSTSato.H.


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