kEpsilonPhitF< BasicTurbulenceModel > Class Template Reference

The k-epsilon-phit-f turbulence closure model for incompressible and compressible flows. More...

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

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef 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 ("kEpsilonPhitF")
 Runtime type information. More...
 
 kEpsilonPhitF (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 ~kEpsilonPhitF ()=default
 Destructor. More...
 
virtual bool read ()
 Re-read model coefficients if they have changed. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon. More...
 
tmp< volScalarFieldDphitEff () const
 Return the effective diffusivity for phit. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulent kinetic energy field. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulent kinetic energy dissipation rate field. More...
 
virtual tmp< volScalarFieldphit () const
 Return the normalised wall-normal fluctuating velocity scale field. More...
 
virtual tmp< volScalarFieldf () const
 Return the elliptic relaxation factor field. More...
 
virtual void correct ()
 Solve the transport equations and correct the turbulent viscosity. 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< 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...
 
- Public Member Functions inherited from kEpsilonPhitFBase
 TypeName ("kEpsilonPhitFBase")
 Runtime type information. More...
 
 kEpsilonPhitFBase ()
 
virtual ~kEpsilonPhitFBase ()
 Destructor. More...
 

Protected Member Functions

virtual void correctNut ()
 Update nut with the latest available k, phit, and T. More...
 
tmp< volScalarFieldTs () const
 Return the turbulent time scale, T. More...
 
tmp< volScalarFieldLs () const
 Return the turbulent length scale, L. More...
 
- 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 Cmu_
 
dimensionedScalar Ceps1a_
 
dimensionedScalar Ceps1b_
 
dimensionedScalar Ceps1c_
 
dimensionedScalar Ceps2_
 
dimensionedScalar Cf1_
 
dimensionedScalar Cf2_
 
dimensionedScalar CL_
 
dimensionedScalar Ceta_
 
dimensionedScalar CT_
 
dimensionedScalar sigmaK_
 
dimensionedScalar sigmaEps_
 
dimensionedScalar sigmaPhit_
 
volScalarField k_
 Turbulent kinetic energy [m2/s2]. More...
 
volScalarField epsilon_
 Turbulent kinetic energy dissipation rate [m2/s3]. More...
 
volScalarField phit_
 Normalised wall-normal fluctuating velocity scale [-]. More...
 
volScalarField f_
 Elliptic relaxation factor [1/s]. More...
 
volScalarField T_
 Turbulent time scale [s]. More...
 
dimensionedScalar phitMin_
 
dimensionedScalar fMin_
 
dimensionedScalar TMin_
 
dimensionedScalar L2Min_
 
- 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::kEpsilonPhitF< BasicTurbulenceModel >

The k-epsilon-phit-f turbulence closure model for incompressible and compressible flows.

The model is a three-transport-equation linear-eddy-viscosity turbulence closure model alongside an elliptic relaxation equation:

  • Turbulent kinetic energy, k,
  • Turbulent kinetic energy dissipation rate, epsilon,
  • Normalised wall-normal fluctuating velocity scale, phit,
  • Elliptic relaxation factor, f.

Reference:

    Standard model (Tag:LUU):
        Laurence, D. R., Uribe, J. C., & Utyuzhnikov, S. V. (2005).
        A robust formulation of the v2-f model.
        Flow, Turbulence and Combustion, 73(3-4), 169–185.
        DOI:10.1007/s10494-005-1974-8

The default model coefficients are (LUU:Eqs. 19-20):

    kEpsilonPhitFCoeffs
    {
        Cmu         0.22,    // Turbulent viscosity constant
        Ceps1a      1.4,     // Model constant for epsilon
        Ceps1b      1.0,     // Model constant for epsilon
        Ceps1c      0.05,    // Model constant for epsilon
        Ceps2       1.9,     // Model constant for epsilon
        Cf1         1.4,     // Model constant for f
        Cf2         0.3,     // Model constant for f
        CL          0.25,    // Model constant for L
        Ceta        110.0,   // Model constant for L
        CT          6.0,     // Model constant for T
        sigmaK      1.0,     // Turbulent Prandtl number for k
        sigmaEps    1.3,     // Turbulent Prandtl number for epsilon
        sigmaPhit   1.0,     // Turbulent Prandtl number for phit = sigmaK
    }
Note
The name of the original variable replacing 'v2' is 'phi' (LUU:Eq. 14). However, the name 'phi' preexisted in OpenFOAM; therefore, this name was replaced by 'phit' herein.
Source files
See also
v2f.H

Definition at line 103 of file kEpsilonPhitF.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 178 of file kEpsilonPhitF.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 179 of file kEpsilonPhitF.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 180 of file kEpsilonPhitF.H.

Constructor & Destructor Documentation

◆ kEpsilonPhitF()

kEpsilonPhitF ( 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 100 of file kEpsilonPhitF.C.

References Foam::bound(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::mag(), Foam::nl, and Foam::type().

Here is the call graph for this function:

◆ ~kEpsilonPhitF()

virtual ~kEpsilonPhitF ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

void correctNut ( )
protectedvirtual

Update nut with the latest available k, phit, and T.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 42 of file kEpsilonPhitF.C.

References optionList::correct(), and options::New().

Here is the call graph for this function:

◆ Ts()

tmp< volScalarField > Ts ( ) const
protected

Return the turbulent time scale, T.

Definition at line 54 of file kEpsilonPhitF.C.

References Foam::max(), nu, Foam::sqrt(), and Foam::Zero.

Here is the call graph for this function:

◆ Ls()

tmp< volScalarField > Ls ( ) const
protected

Return the turbulent length scale, L.

Definition at line 74 of file kEpsilonPhitF.C.

References Foam::max(), nu, Foam::pow(), Foam::pow025(), Foam::pow3(), and Foam::Zero.

Here is the call graph for this function:

◆ TypeName()

TypeName ( "kEpsilonPhitF< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 340 of file kEpsilonPhitF.C.

References Foam::read().

Here is the call graph for this function:

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 213 of file kEpsilonPhitF.H.

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

◆ DepsilonEff()

tmp<volScalarField> DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 226 of file kEpsilonPhitF.H.

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

◆ DphitEff()

tmp<volScalarField> DphitEff ( ) const
inline

Return the effective diffusivity for phit.

Definition at line 239 of file kEpsilonPhitF.H.

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

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulent kinetic energy field.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 252 of file kEpsilonPhitF.H.

References kEpsilonPhitF< BasicTurbulenceModel >::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulent kinetic energy dissipation rate field.

Definition at line 258 of file kEpsilonPhitF.H.

References kEpsilonPhitF< BasicTurbulenceModel >::epsilon_.

◆ phit()

virtual tmp<volScalarField> phit ( ) const
inlinevirtual

Return the normalised wall-normal fluctuating velocity scale field.

Implements kEpsilonPhitFBase.

Definition at line 264 of file kEpsilonPhitF.H.

References kEpsilonPhitF< BasicTurbulenceModel >::phit_.

◆ f()

virtual tmp<volScalarField> f ( ) const
inlinevirtual

Return the elliptic relaxation factor field.

Implements kEpsilonPhitFBase.

Definition at line 270 of file kEpsilonPhitF.H.

References kEpsilonPhitF< BasicTurbulenceModel >::f_.

◆ correct()

void correct ( )
virtual

Solve the transport equations and correct the turbulent viscosity.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 366 of file kEpsilonPhitF.C.

References Foam::fvc::absolute(), Foam::constant::atomic::alpha, Foam::bound(), tmp< T >::clear(), correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), Foam::fvc::div(), divU, fvOptions, Foam::constant::universal::G, Foam::fvc::grad(), Foam::fvc::laplacian(), Foam::fvm::laplacian(), options::New(), nu, nut, phi, tmp< T >::ref(), rho, Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), Foam::fvm::SuSp(), Foam::symm(), Foam::type(), and U.

Here is the call graph for this function:

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 123 of file kEpsilonPhitF.H.

◆ Ceps1a_

dimensionedScalar Ceps1a_
protected

Definition at line 124 of file kEpsilonPhitF.H.

◆ Ceps1b_

dimensionedScalar Ceps1b_
protected

Definition at line 125 of file kEpsilonPhitF.H.

◆ Ceps1c_

dimensionedScalar Ceps1c_
protected

Definition at line 126 of file kEpsilonPhitF.H.

◆ Ceps2_

dimensionedScalar Ceps2_
protected

Definition at line 127 of file kEpsilonPhitF.H.

◆ Cf1_

dimensionedScalar Cf1_
protected

Definition at line 128 of file kEpsilonPhitF.H.

◆ Cf2_

dimensionedScalar Cf2_
protected

Definition at line 129 of file kEpsilonPhitF.H.

◆ CL_

dimensionedScalar CL_
protected

Definition at line 130 of file kEpsilonPhitF.H.

◆ Ceta_

dimensionedScalar Ceta_
protected

Definition at line 131 of file kEpsilonPhitF.H.

◆ CT_

dimensionedScalar CT_
protected

Definition at line 132 of file kEpsilonPhitF.H.

◆ sigmaK_

dimensionedScalar sigmaK_
protected

Definition at line 133 of file kEpsilonPhitF.H.

◆ sigmaEps_

dimensionedScalar sigmaEps_
protected

Definition at line 134 of file kEpsilonPhitF.H.

◆ sigmaPhit_

dimensionedScalar sigmaPhit_
protected

Definition at line 135 of file kEpsilonPhitF.H.

◆ k_

volScalarField k_
protected

Turbulent kinetic energy [m2/s2].

Definition at line 141 of file kEpsilonPhitF.H.

Referenced by kEpsilonPhitF< BasicTurbulenceModel >::k().

◆ epsilon_

volScalarField epsilon_
protected

Turbulent kinetic energy dissipation rate [m2/s3].

Definition at line 144 of file kEpsilonPhitF.H.

Referenced by kEpsilonPhitF< BasicTurbulenceModel >::epsilon().

◆ phit_

volScalarField phit_
protected

Normalised wall-normal fluctuating velocity scale [-].

Definition at line 147 of file kEpsilonPhitF.H.

Referenced by kEpsilonPhitF< BasicTurbulenceModel >::phit().

◆ f_

volScalarField f_
protected

Elliptic relaxation factor [1/s].

Definition at line 150 of file kEpsilonPhitF.H.

Referenced by kEpsilonPhitF< BasicTurbulenceModel >::f().

◆ T_

volScalarField T_
protected

Turbulent time scale [s].

Definition at line 153 of file kEpsilonPhitF.H.

◆ phitMin_

dimensionedScalar phitMin_
protected

Definition at line 158 of file kEpsilonPhitF.H.

◆ fMin_

dimensionedScalar fMin_
protected

Definition at line 159 of file kEpsilonPhitF.H.

◆ TMin_

dimensionedScalar TMin_
protected

Definition at line 160 of file kEpsilonPhitF.H.

◆ L2Min_

dimensionedScalar L2Min_
protected

Definition at line 161 of file kEpsilonPhitF.H.


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