Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
kOmegaSSTSAS< BasicTurbulenceModel > Class Template Reference

Scale-adaptive URAS model based on the k-omega-SST RAS model. More...

Inheritance diagram for kOmegaSSTSAS< BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for kOmegaSSTSAS< BasicTurbulenceModel >:
Collaboration graph
[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 ("kOmegaSSTSAS")
 
 kOmegaSSTSAS (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)
 
virtual ~kOmegaSSTSAS ()=default
 
virtual bool read ()
 
const volScalarFielddelta () const
 
- Public Member Functions inherited from kOmegaSST< BasicTurbulenceModel >
 TypeName ("kOmegaSST")
 
 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)
 
virtual ~kOmegaSST ()=default
 
- 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)
 
virtual ~kOmegaSSTBase ()=default
 
tmp< volScalarFieldDkEff (const volScalarField &F1) const
 
tmp< volScalarFieldDomegaEff (const volScalarField &F1) const
 
virtual tmp< volScalarFieldk () const
 
virtual tmp< volScalarFieldomega () const
 
virtual void correct ()
 
- 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)
 
virtual ~eddyViscosity ()=default
 
virtual tmp< volScalarFieldnut () const
 
virtual tmp< scalarFieldnut (const label patchi) const
 
virtual tmp< volSymmTensorFieldR () const
 
virtual void validate ()
 
- 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)
 
virtual ~linearViscousStress ()=default
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 
- Public Member Functions inherited from RASModel< BasicTurbulenceModel >
 TypeName ("RAS")
 
 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)
 
virtual ~RASModel ()=default
 
const dimensionedScalarkMin () const
 
const dimensionedScalarepsilonMin () const
 
const dimensionedScalaromegaMin () const
 
dimensionedScalarkMin ()
 
dimensionedScalarepsilonMin ()
 
dimensionedScalaromegaMin ()
 
virtual const dictionarycoeffDict () const
 
virtual tmp< volScalarFieldnuEff () const
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 
virtual tmp< volScalarFieldepsilon () const
 
virtual tmp< volScalarFieldomega () const
 

Protected Member Functions

virtual tmp< fvScalarMatrixQsas (const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
 
- 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< 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
 
virtual tmp< volScalarField::InternalepsilonByk (const volScalarField &F1, const volTensorField &gradU) const
 
virtual tmp< volScalarField::InternalGbyNu (const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
 
virtual tmp< fvScalarMatrixkSource () const
 
virtual tmp< fvScalarMatrixomegaSource () const
 
- Protected Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 
 RASModel (const RASModel &)=delete
 
void operator= (const RASModel &)=delete
 

Protected Attributes

dimensionedScalar Cs_
 
dimensionedScalar kappa_
 
dimensionedScalar zeta2_
 
dimensionedScalar sigmaPhi_
 
dimensionedScalar C_
 
autoPtr< Foam::LESdeltadelta_
 
- 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_
 
const volScalarFieldy_
 
volScalarField k_
 
volScalarField omega_
 
Switch decayControl_
 
dimensionedScalar kInf_
 
dimensionedScalar omegaInf_
 
- Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_
 
- Protected Attributes inherited from RASModel< BasicTurbulenceModel >
dictionary RASDict_
 
Switch turbulence_
 
Switch printCoeffs_
 
dictionary coeffDict_
 
dimensionedScalar kMin_
 
dimensionedScalar epsilonMin_
 
dimensionedScalar omegaMin_
 

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)
 

Detailed Description

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

Scale-adaptive URAS model based on the k-omega-SST RAS model.

References:

    Egorov, Y., & Menter F.R. (2008).
    Development and Application of SST-SAS Model in the DESIDER Project.
    Advances in Hybrid RANS-LES Modelling,
    Notes on Num. Fluid Mech. And Multidisciplinary Design,
    Volume 97, 261-270.

The model coefficients are

    kOmegaSSTSASCoeffs
    {
        // Default SST coefficients
        alphaK1     0.85;
        alphaK2     1.0;
        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.0;
        c1          10.0;
        F3          no;

        // Default SAS coefficients
        Cs          0.11;
        kappa       0.41;
        zeta2       3.51;
        sigmaPhi    2.0/3.0;
        C           2;

        // Delta must be specified for SAS e.g.
        delta cubeRootVol;

        cubeRootVolCoeffs
        {}
    }
Source files

Definition at line 97 of file kOmegaSSTSAS.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 142 of file kOmegaSSTSAS.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 143 of file kOmegaSSTSAS.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 144 of file kOmegaSSTSAS.H.

Constructor & Destructor Documentation

◆ kOmegaSSTSAS()

kOmegaSSTSAS ( 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 
)

Definition at line 93 of file kOmegaSSTSAS.C.

References Foam::type().

Here is the call graph for this function:

◆ ~kOmegaSSTSAS()

virtual ~kOmegaSSTSAS ( )
virtualdefault

Member Function Documentation

◆ Qsas()

tmp< fvScalarMatrix > Qsas ( const volScalarField::Internal S2,
const volScalarField::Internal gamma,
const volScalarField::Internal beta 
) const
protectedvirtual

◆ TypeName()

TypeName ( "kOmegaSSTSAS< BasicTurbulenceModel >"  )

◆ read()

bool read
virtual

◆ delta()

const volScalarField& delta ( ) const
inline

Definition at line 177 of file kOmegaSSTSAS.H.

References kOmegaSSTSAS< BasicTurbulenceModel >::delta_.

Member Data Documentation

◆ Cs_

dimensionedScalar Cs_
protected

Definition at line 116 of file kOmegaSSTSAS.H.

◆ kappa_

dimensionedScalar kappa_
protected

Definition at line 117 of file kOmegaSSTSAS.H.

◆ zeta2_

dimensionedScalar zeta2_
protected

Definition at line 118 of file kOmegaSSTSAS.H.

◆ sigmaPhi_

dimensionedScalar sigmaPhi_
protected

Definition at line 119 of file kOmegaSSTSAS.H.

◆ C_

dimensionedScalar C_
protected

Definition at line 120 of file kOmegaSSTSAS.H.

◆ delta_

autoPtr<Foam::LESdelta> delta_
protected

Definition at line 126 of file kOmegaSSTSAS.H.

Referenced by kOmegaSSTSAS< BasicTurbulenceModel >::delta().


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