qZeta Class Reference

Gibson and Dafa'Alla's q-zeta two-equation low-Re turbulence model for incompressible flows. More...

Inheritance diagram for qZeta:
Collaboration diagram for qZeta:

Public Member Functions

 TypeName ("qZeta")
 Runtime type information. More...
 
 qZeta (const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
 Construct from components. More...
 
virtual ~qZeta ()
 Destructor. More...
 
virtual bool read ()
 Read RASProperties dictionary. More...
 
const dimensionedScalarqMin () const
 Return the lower allowable limit for q (default: small) More...
 
const dimensionedScalarzetaMin () const
 Return the lower allowable limit for zeta (default: small) More...
 
dimensionedScalarqMin ()
 Allow qMin to be changed. More...
 
dimensionedScalarzetaMin ()
 Allow zetaMin to be changed. More...
 
tmp< volScalarFieldDqEff () const
 Return the effective diffusivity for q. More...
 
tmp< volScalarFieldDzetaEff () const
 Return the effective diffusivity for epsilon. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 
virtual const volScalarFieldq () const
 
virtual const volScalarFieldzeta () const
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
- Public Member Functions inherited from eddyViscosity< incompressible::RASModel >
 eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
virtual ~eddyViscosity ()
 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< volSymmTensorFieldsigma () const
 Return the Reynolds stress tensor [m^2/s^2]. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
- Public Member Functions inherited from linearViscousStress< incompressible::RASModel >
 linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
virtual ~linearViscousStress ()
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevTau () const
 Return the effective stress tensor. More...
 
virtual tmp< fvVectorMatrixdivDevTau (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevTau (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 

Protected Member Functions

tmp< volScalarFieldfMu () const
 
tmp< volScalarFieldf2 () const
 
virtual void correctNut ()
 

Protected Attributes

dimensionedScalar Cmu_
 
dimensionedScalar C1_
 
dimensionedScalar C2_
 
dimensionedScalar sigmaZeta_
 
Switch anisotropic_
 
dimensionedScalar qMin_
 Lower limit of q. More...
 
dimensionedScalar zetaMin_
 Lower limit of zeta. More...
 
volScalarField k_
 
volScalarField epsilon_
 
volScalarField q_
 
volScalarField zeta_
 
- Protected Attributes inherited from eddyViscosity< incompressible::RASModel >
volScalarField nut_
 

Additional Inherited Members

- Public Types inherited from eddyViscosity< incompressible::RASModel >
typedef incompressible::RASModel ::alphaField alphaField
 
typedef incompressible::RASModel ::rhoField rhoField
 
typedef incompressible::RASModel ::transportModel transportModel
 
- Public Types inherited from linearViscousStress< incompressible::RASModel >
typedef incompressible::RASModel ::alphaField alphaField
 
typedef incompressible::RASModel ::rhoField rhoField
 
typedef incompressible::RASModel ::transportModel transportModel
 

Detailed Description

Gibson and Dafa'Alla's q-zeta two-equation low-Re turbulence model for incompressible flows.

This turbulence model is described in:

    Dafa'Alla, A.A., Juntasaro, E. & Gibson, M.M. (1996).
    Calculation of oscillating boundary layers with the
    q-zeta turbulence model.
    Engineering Turbulence Modelling and Experiments 3:
    Proceedings of the Third International Symposium,
    Crete, Greece, May 27-29, 141.

which is a development of the original q-zeta model described in:

    Gibson, M. M., & Dafa'Alla, A. A. (1995).
    Two-equation model for turbulent wall flow.
    AIAA journal, 33(8), 1514-1518.
Source files

Definition at line 70 of file qZeta.H.

Constructor & Destructor Documentation

◆ qZeta()

qZeta ( const geometricOneField alpha,
const geometricOneField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport,
const word type = typeName 
)

Construct from components.

Definition at line 80 of file qZeta.C.

References Foam::bound(), qZeta::zeta_, and qZeta::zetaMin_.

Referenced by qZeta::correctNut().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ~qZeta()

virtual ~qZeta ( )
inlinevirtual

Destructor.

Definition at line 130 of file qZeta.H.

References qZeta::read().

Here is the call graph for this function:

Member Function Documentation

◆ fMu()

tmp< volScalarField > fMu ( ) const
protected

Definition at line 46 of file qZeta.C.

References qZeta::anisotropic_, Foam::exp(), qZeta::k_, Foam::pow3(), qZeta::q_, Foam::sqr(), and qZeta::zeta_.

Referenced by qZeta::correctNut().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ f2()

tmp< volScalarField > f2 ( ) const
protected

Definition at line 63 of file qZeta.C.

References Foam::exp(), qZeta::k_, qZeta::q_, Foam::sqr(), and qZeta::zeta_.

Referenced by qZeta::correct().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ correctNut()

void correctNut ( )
protectedvirtual

Implements eddyViscosity< incompressible::RASModel >.

Definition at line 70 of file qZeta.C.

References qZeta::Cmu_, GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), qZeta::epsilon_, qZeta::fMu(), qZeta::k_, eddyViscosity< incompressible::RASModel >::nut_, qZeta::qZeta(), and Foam::sqr().

Referenced by qZeta::correct().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ TypeName()

TypeName ( "qZeta"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Read RASProperties dictionary.

Implements eddyViscosity< incompressible::RASModel >.

Definition at line 215 of file qZeta.C.

References qZeta::anisotropic_, qZeta::C1_, qZeta::C2_, qZeta::Cmu_, qZeta::qMin_, Switch::readIfPresent(), dimensioned< Type >::readIfPresent(), qZeta::sigmaZeta_, and qZeta::zetaMin_.

Referenced by qZeta::~qZeta().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ qMin() [1/2]

const dimensionedScalar& qMin ( ) const
inline

Return the lower allowable limit for q (default: small)

Definition at line 140 of file qZeta.H.

References qZeta::qMin_.

◆ zetaMin() [1/2]

const dimensionedScalar& zetaMin ( ) const
inline

Return the lower allowable limit for zeta (default: small)

Definition at line 146 of file qZeta.H.

References qZeta::zetaMin_.

◆ qMin() [2/2]

dimensionedScalar& qMin ( )
inline

Allow qMin to be changed.

Definition at line 152 of file qZeta.H.

References qZeta::qMin_.

◆ zetaMin() [2/2]

dimensionedScalar& zetaMin ( )
inline

Allow zetaMin to be changed.

Definition at line 158 of file qZeta.H.

References qZeta::zetaMin_.

◆ DqEff()

tmp<volScalarField> DqEff ( ) const
inline

Return the effective diffusivity for q.

Definition at line 164 of file qZeta.H.

References GeometricField< scalar, fvPatchField, volMesh >::New(), and eddyViscosity< incompressible::RASModel >::nut_.

Referenced by qZeta::correct().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ DzetaEff()

tmp<volScalarField> DzetaEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 174 of file qZeta.H.

References GeometricField< scalar, fvPatchField, volMesh >::New(), and eddyViscosity< incompressible::RASModel >::nut_.

Referenced by qZeta::correct().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< incompressible::RASModel >.

Definition at line 184 of file qZeta.H.

References qZeta::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 190 of file qZeta.H.

References qZeta::epsilon_.

◆ q()

virtual const volScalarField& q ( ) const
inlinevirtual

Definition at line 195 of file qZeta.H.

References qZeta::q_.

◆ zeta()

virtual const volScalarField& zeta ( ) const
inlinevirtual

Definition at line 200 of file qZeta.H.

References qZeta::correct(), and qZeta::zeta_.

Here is the call graph for this function:

◆ correct()

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 81 of file qZeta.H.

Referenced by qZeta::correctNut(), and qZeta::read().

◆ C1_

dimensionedScalar C1_
protected

Definition at line 82 of file qZeta.H.

Referenced by qZeta::correct(), and qZeta::read().

◆ C2_

dimensionedScalar C2_
protected

Definition at line 83 of file qZeta.H.

Referenced by qZeta::correct(), and qZeta::read().

◆ sigmaZeta_

dimensionedScalar sigmaZeta_
protected

Definition at line 84 of file qZeta.H.

Referenced by qZeta::read().

◆ anisotropic_

Switch anisotropic_
protected

Definition at line 85 of file qZeta.H.

Referenced by qZeta::fMu(), and qZeta::read().

◆ qMin_

dimensionedScalar qMin_
protected

Lower limit of q.

Definition at line 88 of file qZeta.H.

Referenced by qZeta::correct(), qZeta::qMin(), and qZeta::read().

◆ zetaMin_

dimensionedScalar zetaMin_
protected

Lower limit of zeta.

Definition at line 91 of file qZeta.H.

Referenced by qZeta::correct(), qZeta::qZeta(), qZeta::read(), and qZeta::zetaMin().

◆ k_

volScalarField k_
protected

Definition at line 95 of file qZeta.H.

Referenced by qZeta::correct(), qZeta::correctNut(), qZeta::f2(), qZeta::fMu(), and qZeta::k().

◆ epsilon_

volScalarField epsilon_
protected

Definition at line 96 of file qZeta.H.

Referenced by qZeta::correct(), qZeta::correctNut(), and qZeta::epsilon().

◆ q_

volScalarField q_
protected

Definition at line 98 of file qZeta.H.

Referenced by qZeta::correct(), qZeta::f2(), qZeta::fMu(), and qZeta::q().

◆ zeta_

volScalarField zeta_
protected

Definition at line 99 of file qZeta.H.

Referenced by qZeta::correct(), qZeta::f2(), qZeta::fMu(), qZeta::qZeta(), and qZeta::zeta().


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