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

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

Inheritance diagram for kEpsilonPhitF< BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for kEpsilonPhitF< BasicTurbulenceModel >:
Collaboration graph
[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")
 
 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)
 
virtual ~kEpsilonPhitF ()=default
 
virtual bool read ()
 
tmp< volScalarFieldDkEff () const
 
tmp< volScalarFieldDepsilonEff () const
 
tmp< volScalarFieldDphitEff () const
 
virtual tmp< volScalarFieldk () const
 
virtual tmp< volScalarFieldepsilon () const
 
virtual tmp< volScalarFieldphit () const
 
virtual tmp< volScalarFieldf () 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< volScalarFieldk () const=0
 
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< volScalarFieldomega () const
 

Protected Member Functions

virtual void correctNut ()
 
tmp< volScalarFieldTs () const
 
tmp< volScalarFieldLs () 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

Switch includeNu_
 
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_
 
volScalarField epsilon_
 
volScalarField phit_
 
volScalarField f_
 
volScalarField T_
 
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_
 
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::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.


Input fields

k : Turbulent kinetic energy [m2/s2]
epsilon : Turbulent kinetic energy dissipation rate [m2/s3]
phit : Normalised wall-normal fluctuating velocity scale [-]
f : Elliptic relaxation factor [1/s]

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
        {
            includeNu   true;    // include nu in (LUU: Eq. 17), see Notes
            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.

Including nu in DphitEff even though it is not present in (LUU:Eq. 17) provided higher level of resemblance to benchmarks for the tests considered, particularly for the peak skin friction (yet, pressure-related predictions were unaffected). Users can switch off nu in DphitEff by using includeNu entry in kEpsilonPhitFCoeffs as shown above in order to follow the reference paper thereat. includeNu is left true by default. See GitLab issue #1560.

Source files
See also
kEpsilon.C

Definition at line 127 of file kEpsilonPhitF.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 203 of file kEpsilonPhitF.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 204 of file kEpsilonPhitF.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 205 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 
)

Definition at line 93 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

Member Function Documentation

◆ correctNut()

void correctNut
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 35 of file kEpsilonPhitF.C.

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

Here is the call graph for this function:

◆ Ts()

tmp< volScalarField > Ts
protected

Definition at line 47 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
protected

Definition at line 67 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 >"  )

◆ read()

bool read
virtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 341 of file kEpsilonPhitF.C.

References Foam::read().

Here is the call graph for this function:

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

Definition at line 238 of file kEpsilonPhitF.H.

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

◆ DepsilonEff()

tmp<volScalarField> DepsilonEff ( ) const
inline

Definition at line 251 of file kEpsilonPhitF.H.

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

◆ DphitEff()

tmp<volScalarField> DphitEff ( ) const
inline

Definition at line 264 of file kEpsilonPhitF.H.

References kEpsilonPhitF< BasicTurbulenceModel >::includeNu_, tmp::New(), nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.

Here is the call graph for this function:

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Definition at line 278 of file kEpsilonPhitF.H.

References kEpsilonPhitF< BasicTurbulenceModel >::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Reimplemented from RASModel< BasicTurbulenceModel >.

Definition at line 284 of file kEpsilonPhitF.H.

References kEpsilonPhitF< BasicTurbulenceModel >::epsilon_.

◆ phit()

virtual tmp<volScalarField> phit ( ) const
inlinevirtual

Definition at line 290 of file kEpsilonPhitF.H.

References kEpsilonPhitF< BasicTurbulenceModel >::phit_.

◆ f()

virtual tmp<volScalarField> f ( ) const
inlinevirtual

Definition at line 296 of file kEpsilonPhitF.H.

References kEpsilonPhitF< BasicTurbulenceModel >::f_.

◆ correct()

void correct
virtual

Member Data Documentation

◆ includeNu_

Switch includeNu_
protected

Definition at line 144 of file kEpsilonPhitF.H.

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

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 148 of file kEpsilonPhitF.H.

◆ Ceps1a_

dimensionedScalar Ceps1a_
protected

Definition at line 149 of file kEpsilonPhitF.H.

◆ Ceps1b_

dimensionedScalar Ceps1b_
protected

Definition at line 150 of file kEpsilonPhitF.H.

◆ Ceps1c_

dimensionedScalar Ceps1c_
protected

Definition at line 151 of file kEpsilonPhitF.H.

◆ Ceps2_

dimensionedScalar Ceps2_
protected

Definition at line 152 of file kEpsilonPhitF.H.

◆ Cf1_

dimensionedScalar Cf1_
protected

Definition at line 153 of file kEpsilonPhitF.H.

◆ Cf2_

dimensionedScalar Cf2_
protected

Definition at line 154 of file kEpsilonPhitF.H.

◆ CL_

dimensionedScalar CL_
protected

Definition at line 155 of file kEpsilonPhitF.H.

◆ Ceta_

dimensionedScalar Ceta_
protected

Definition at line 156 of file kEpsilonPhitF.H.

◆ CT_

dimensionedScalar CT_
protected

Definition at line 157 of file kEpsilonPhitF.H.

◆ sigmaK_

dimensionedScalar sigmaK_
protected

Definition at line 158 of file kEpsilonPhitF.H.

◆ sigmaEps_

dimensionedScalar sigmaEps_
protected

Definition at line 159 of file kEpsilonPhitF.H.

◆ sigmaPhit_

dimensionedScalar sigmaPhit_
protected

Definition at line 160 of file kEpsilonPhitF.H.

◆ k_

volScalarField k_
protected

Definition at line 166 of file kEpsilonPhitF.H.

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

◆ epsilon_

volScalarField epsilon_
protected

Definition at line 169 of file kEpsilonPhitF.H.

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

◆ phit_

volScalarField phit_
protected

Definition at line 172 of file kEpsilonPhitF.H.

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

◆ f_

volScalarField f_
protected

Definition at line 175 of file kEpsilonPhitF.H.

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

◆ T_

volScalarField T_
protected

Definition at line 178 of file kEpsilonPhitF.H.

◆ phitMin_

dimensionedScalar phitMin_
protected

Definition at line 183 of file kEpsilonPhitF.H.

◆ fMin_

dimensionedScalar fMin_
protected

Definition at line 184 of file kEpsilonPhitF.H.

◆ TMin_

dimensionedScalar TMin_
protected

Definition at line 185 of file kEpsilonPhitF.H.

◆ L2Min_

dimensionedScalar L2Min_
protected

Definition at line 186 of file kEpsilonPhitF.H.


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