kOmegaSSTLM< BasicTurbulenceModel > Class Template Reference

Langtry-Menter 4-equation transitional SST model based on the k-omega-SST RAS model. More...

Inheritance diagram for kOmegaSSTLM< BasicTurbulenceModel >:
[legend]
Collaboration diagram for kOmegaSSTLM< 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 ("kOmegaSSTLM")
 Runtime type information. More...
 
 kOmegaSSTLM (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 ~kOmegaSSTLM ()=default
 Destructor. More...
 
virtual bool read ()
 Re-read model coefficients if they have changed. More...
 
const volScalarFieldReThetat () const
 Access function transition onset momentum-thickness Reynolds number. More...
 
const volScalarFieldgammaInt () const
 Access function to intermittency. More...
 
tmp< volScalarFieldDReThetatEff () const
 Return the effective diffusivity for transition onset. More...
 
tmp< volScalarFieldDgammaIntEff () const
 Return the effective diffusivity for intermittency. 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 tmp< volScalarFieldF1 (const volScalarField &CDkOmega) const
 Modified form of the k-omega SST F1 function. More...
 
virtual tmp< volScalarField::InternalPk (const volScalarField::Internal &G) const
 Modified form of the k-omega SST k production rate. More...
 
virtual tmp< volScalarField::InternalepsilonByk (const volScalarField &F1, const volTensorField &gradU) const
 Modified form of the k-omega SST epsilon/k. More...
 
tmp< volScalarField::InternalFthetat (const volScalarField::Internal &Us, const volScalarField::Internal &Omega, const volScalarField::Internal &nu) const
 Freestream blending-function. More...
 
tmp< volScalarField::InternalReThetac () const
 Empirical correlation for critical Reynolds number where the. More...
 
tmp< volScalarField::InternalFlength (const volScalarField::Internal &nu) const
 Empirical correlation that controls the length of the. More...
 
tmp< volScalarField::InternalFonset (const volScalarField::Internal &Rev, const volScalarField::Internal &ReThetac, const volScalarField::Internal &RT) const
 Transition onset location control function. More...
 
tmp< volScalarField::InternalReThetat0 (const volScalarField::Internal &Us, const volScalarField::Internal &dUsds, const volScalarField::Internal &nu) const
 Return the transition onset momentum-thickness Reynolds number. More...
 
void correctReThetatGammaInt ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
- Protected Member Functions inherited from kOmegaSST< BasicTurbulenceModel >
virtual void correctNut (const volScalarField &S2)
 
virtual void correctNut ()
 
- Protected Member Functions inherited from kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >
void setDecayControl (const dictionary &dict)
 
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::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 ca1_
 
dimensionedScalar ca2_
 
dimensionedScalar ce1_
 
dimensionedScalar ce2_
 
dimensionedScalar cThetat_
 
dimensionedScalar sigmaThetat_
 
scalar lambdaErr_
 Convergence criterion for the lambda/thetat loop. More...
 
label maxLambdaIter_
 Maximum number of iterations to converge the lambda/thetat loop. More...
 
const dimensionedScalar deltaU_
 Stabilization for division by the magnitude of the velocity. More...
 
volScalarField ReThetat_
 Transition onset momentum-thickness Reynolds number. More...
 
volScalarField gammaInt_
 Intermittency. More...
 
volScalarField::Internal gammaIntEff_
 Effective intermittency. More...
 
- 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::kOmegaSSTLM< BasicTurbulenceModel >

Langtry-Menter 4-equation transitional SST model based on the k-omega-SST RAS model.

References:

    Langtry, R. B., & Menter, F. R. (2009).
    Correlation-based transition modeling for unstructured parallelized
    computational fluid dynamics codes.
    AIAA journal, 47(12), 2894-2906.

    Menter, F. R., Langtry, R., & Volker, S. (2006).
    Transition modelling for general purpose CFD codes.
    Flow, turbulence and combustion, 77(1-4), 277-303.

    Langtry, R. B. (2006).
    A correlation-based transition model using local variables for
    unstructured parallelized CFD codes.
    Phd. Thesis, Universität Stuttgart.

The model coefficients are

    kOmegaSSTCoeffs
    {
        // Default SST coefficients
        alphaK1     0.85;
        alphaK2     1;
        alphaOmega1 0.5;
        alphaOmega2 0.856;
        beta1       0.075;
        beta2       0.0828;
        betaStar    0.09;
        gamma1      5/9;
        gamma2      0.44;
        a1          0.31;
        b1          1;
        c1          10;
        F3          no;

        // Default LM coefficients
        ca1         2;
        ca2         0.06;
        ce1         1;
        ce2         50;
        cThetat     0.03;
        sigmaThetat 2;

        lambdaErr   1e-6;
        maxLambdaIter 10;
    }
Source files

Definition at line 107 of file kOmegaSSTLM.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 217 of file kOmegaSSTLM.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 218 of file kOmegaSSTLM.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 219 of file kOmegaSSTLM.H.

Constructor & Destructor Documentation

◆ kOmegaSSTLM()

kOmegaSSTLM ( 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 367 of file kOmegaSSTLM.C.

References Foam::type().

Here is the call graph for this function:

◆ ~kOmegaSSTLM()

virtual ~kOmegaSSTLM ( )
virtualdefault

Destructor.

Member Function Documentation

◆ F1()

tmp< volScalarField > F1 ( const volScalarField CDkOmega) const
protectedvirtual

Modified form of the k-omega SST F1 function.

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

Definition at line 44 of file kOmegaSSTLM.C.

References Foam::exp(), F3, Foam::max(), nu, Foam::pow(), Foam::Ry(), and Foam::sqrt().

Here is the call graph for this function:

◆ Pk()

tmp< volScalarField::Internal > Pk ( const volScalarField::Internal G) const
protectedvirtual

Modified form of the k-omega SST k production rate.

Definition at line 57 of file kOmegaSSTLM.C.

References Foam::constant::universal::G, and kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::Pk().

Here is the call graph for this function:

◆ epsilonByk()

tmp< volScalarField::Internal > epsilonByk ( const volScalarField F1,
const volTensorField gradU 
) const
protectedvirtual

Modified form of the k-omega SST epsilon/k.

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

Definition at line 67 of file kOmegaSSTLM.C.

References kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::epsilonByk(), F1, Foam::max(), and Foam::min().

Here is the call graph for this function:

◆ Fthetat()

tmp< volScalarField::Internal > Fthetat ( const volScalarField::Internal Us,
const volScalarField::Internal Omega,
const volScalarField::Internal nu 
) const
protected

Freestream blending-function.

Definition at line 80 of file kOmegaSSTLM.C.

References delta, Foam::dimLength, Foam::exp(), IOobject::groupName(), Foam::max(), Foam::min(), nu, Foam::pow4(), Foam::sqr(), Us, and y.

Here is the call graph for this function:

◆ ReThetac()

tmp< volScalarField::Internal > ReThetac ( ) const
protected

Empirical correlation for critical Reynolds number where the.

intermittency first starts to increase in the boundary layer

Definition at line 119 of file kOmegaSSTLM.C.

References Foam::dimless, forAll, IOobject::groupName(), Foam::pow3(), Foam::pow4(), and Foam::sqr().

Here is the call graph for this function:

◆ Flength()

tmp< volScalarField::Internal > Flength ( const volScalarField::Internal nu) const
protected

Empirical correlation that controls the length of the.

transition region

Definition at line 160 of file kOmegaSSTLM.C.

References Foam::dimless, Foam::constant::electromagnetic::e, Foam::exp(), forAll, IOobject::groupName(), nu, Foam::pow3(), Foam::sqr(), and y.

Here is the call graph for this function:

◆ Fonset()

tmp< volScalarField::Internal > Fonset ( const volScalarField::Internal Rev,
const volScalarField::Internal ReThetac,
const volScalarField::Internal RT 
) const
protected

Transition onset location control function.

Definition at line 337 of file kOmegaSSTLM.C.

References IOobject::groupName(), Foam::max(), Foam::min(), Foam::pow3(), and Foam::pow4().

Here is the call graph for this function:

◆ ReThetat0()

tmp< volScalarField::Internal > ReThetat0 ( const volScalarField::Internal Us,
const volScalarField::Internal dUsds,
const volScalarField::Internal nu 
) const
protected

Return the transition onset momentum-thickness Reynolds number.

(based on freestream conditions)

Definition at line 223 of file kOmegaSSTLM.C.

References Foam::dimless, Foam::endl(), Foam::exp(), forAll, IOobject::groupName(), k, lambda(), Foam::mag(), Foam::max(), Foam::min(), nu, Foam::pow(), Foam::pow3(), Foam::sqr(), Foam::sqrt(), Us, and WarningInFunction.

Here is the call graph for this function:

◆ correctReThetatGammaInt()

void correctReThetatGammaInt ( )
protected

Solve the turbulence equations and correct the turbulence viscosity.

Definition at line 523 of file kOmegaSSTLM.C.

References Foam::constant::atomic::alpha, Foam::bound(), tmp< T >::clear(), Foam::fvm::ddt(), Foam::fvm::div(), Foam::exp(), fvOptions, Foam::fvc::grad(), k, Foam::fvm::laplacian(), Foam::mag(), Foam::magSqr(), Foam::max(), Foam::min(), options::New(), nu, Foam::pow4(), tmp< T >::ref(), rho, Foam::skew(), Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), Foam::symm(), U, Us, and y.

Here is the call graph for this function:

◆ TypeName()

TypeName ( "kOmegaSSTLM< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Re-read model coefficients if they have changed.

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

Definition at line 502 of file kOmegaSSTLM.C.

◆ ReThetat()

const volScalarField& ReThetat ( ) const
inline

Access function transition onset momentum-thickness Reynolds number.

Definition at line 252 of file kOmegaSSTLM.H.

References kOmegaSSTLM< BasicTurbulenceModel >::ReThetat_.

◆ gammaInt()

const volScalarField& gammaInt ( ) const
inline

Access function to intermittency.

Definition at line 258 of file kOmegaSSTLM.H.

References kOmegaSSTLM< BasicTurbulenceModel >::gammaInt_.

◆ DReThetatEff()

tmp<volScalarField> DReThetatEff ( ) const
inline

Return the effective diffusivity for transition onset.

momentum-thickness Reynolds number

Definition at line 265 of file kOmegaSSTLM.H.

◆ DgammaIntEff()

tmp<volScalarField> DgammaIntEff ( ) const
inline

Return the effective diffusivity for intermittency.

Definition at line 278 of file kOmegaSSTLM.H.

References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.

◆ correct()

void correct ( )
virtual

Solve the turbulence equations and correct the turbulence viscosity.

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

Definition at line 621 of file kOmegaSSTLM.C.

References kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct().

Here is the call graph for this function:

Member Data Documentation

◆ ca1_

dimensionedScalar ca1_
protected

Definition at line 126 of file kOmegaSSTLM.H.

◆ ca2_

dimensionedScalar ca2_
protected

Definition at line 127 of file kOmegaSSTLM.H.

◆ ce1_

dimensionedScalar ce1_
protected

Definition at line 129 of file kOmegaSSTLM.H.

◆ ce2_

dimensionedScalar ce2_
protected

Definition at line 130 of file kOmegaSSTLM.H.

◆ cThetat_

dimensionedScalar cThetat_
protected

Definition at line 132 of file kOmegaSSTLM.H.

◆ sigmaThetat_

dimensionedScalar sigmaThetat_
protected

Definition at line 133 of file kOmegaSSTLM.H.

◆ lambdaErr_

scalar lambdaErr_
protected

Convergence criterion for the lambda/thetat loop.

Definition at line 136 of file kOmegaSSTLM.H.

◆ maxLambdaIter_

label maxLambdaIter_
protected

Maximum number of iterations to converge the lambda/thetat loop.

Definition at line 139 of file kOmegaSSTLM.H.

◆ deltaU_

const dimensionedScalar deltaU_
protected

Stabilization for division by the magnitude of the velocity.

Definition at line 142 of file kOmegaSSTLM.H.

◆ ReThetat_

volScalarField ReThetat_
protected

Transition onset momentum-thickness Reynolds number.

Definition at line 148 of file kOmegaSSTLM.H.

Referenced by kOmegaSSTLM< BasicTurbulenceModel >::ReThetat().

◆ gammaInt_

volScalarField gammaInt_
protected

Intermittency.

Definition at line 151 of file kOmegaSSTLM.H.

Referenced by kOmegaSSTLM< BasicTurbulenceModel >::gammaInt().

◆ gammaIntEff_

volScalarField::Internal gammaIntEff_
protected

Effective intermittency.

Definition at line 154 of file kOmegaSSTLM.H.


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