Loading...
Searching...
No Matches
kEpsilonPhitF< BasicTurbulenceModel > Class Template Reference

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

#include <kEpsilonPhitF.H>

Inheritance diagram for kEpsilonPhitF< BasicTurbulenceModel >:
Collaboration diagram for kEpsilonPhitF< BasicTurbulenceModel >:

Public Types

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

Public Member Functions

 TypeName ("kEpsilonPhitF")
 Runtime type information.
 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.
virtual ~kEpsilonPhitF ()=default
 Destructor.
virtual bool read ()
 Re-read model coefficients if they have changed.
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k (LUU:Eq. 3).
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon (LUU:Eq. 4).
tmp< volScalarFieldDphitEff () const
 Return the effective diffusivity for phit (LUU:Eq. 17).
virtual tmp< volScalarFieldk () const
 Return the turbulent kinetic energy field.
virtual tmp< volScalarFieldepsilon () const
 Return the turbulent kinetic energy dissipation rate field.
virtual tmp< volScalarFieldphit () const
 Return the normalised wall-normal fluctuating velocity scale field.
virtual tmp< volScalarFieldf () const
 Return the elliptic relaxation factor field.
virtual void correct ()
 Solve the transport equations and correct the turbulent viscosity.
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.
virtual ~eddyViscosity ()=default
 Destructor.
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity.
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor.
virtual void validate ()
 Validate the turbulence fields after construction.
Public Member Functions inherited from linearViscousStress< 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.
virtual ~linearViscousStress ()=default
 Destructor.
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field.
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation.
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation.

Protected Member Functions

virtual void correctNut ()
 Update nut with the latest available k, phit, and T.
tmp< volScalarFieldTs () const
 Return the turbulent time scale, T.
tmp< volScalarFieldLs () const
 Return the turbulent length scale, L.

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_
 Turbulent kinetic energy [m2/s2].
volScalarField epsilon_
 Turbulent kinetic energy dissipation rate [m2/s3].
volScalarField phit_
 Normalised wall-normal fluctuating velocity scale [-].
volScalarField f_
 Elliptic relaxation factor [1/s].
volScalarField T_
 Turbulent time scale [s].
dimensionedScalar phitMin_
dimensionedScalar fMin_
dimensionedScalar TMin_
dimensionedScalar L2Min_
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_

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

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 223 of file kEpsilonPhitF.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 224 of file kEpsilonPhitF.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 225 of file kEpsilonPhitF.H.

Constructor & Destructor Documentation

◆ kEpsilonPhitF()

template<class BasicTurbulenceModel>
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 )

◆ ~kEpsilonPhitF()

template<class BasicTurbulenceModel>
virtual ~kEpsilonPhitF ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

template<class BasicTurbulenceModel>
void correctNut ( )
protectedvirtual

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

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 35 of file kEpsilonPhitF.C.

References Cmu_, optionList::correct(), k_, options::New(), eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_, phit_, and T_.

Referenced by correct().

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

◆ Ts()

template<class BasicTurbulenceModel>
tmp< volScalarField > Ts ( ) const
protected

Return the turbulent time scale, T.

Definition at line 47 of file kEpsilonPhitF.C.

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

Referenced by correct().

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

◆ Ls()

template<class BasicTurbulenceModel>
tmp< volScalarField > Ls ( ) const
protected

Return the turbulent length scale, L.

Definition at line 67 of file kEpsilonPhitF.C.

References Ceta_, CL_, epsilon_, k_, Foam::max(), nu, Foam::pow(), Foam::pow025(), Foam::pow3(), and Foam::Zero.

Referenced by correct().

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

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "kEpsilonPhitF< BasicTurbulenceModel >" )

Runtime type information.

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 341 of file kEpsilonPhitF.C.

References Ceps1a_, Ceps1b_, Ceps1c_, Ceps2_, Ceta_, Cf1_, Cf2_, CL_, Cmu_, CT_, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), includeNu_, read(), sigmaEps_, sigmaK_, and sigmaPhit_.

Referenced by read().

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

◆ DkEff()

template<class BasicTurbulenceModel>
tmp< volScalarField > DkEff ( ) const
inline

Return the effective diffusivity for k (LUU:Eq. 3).

Definition at line 268 of file kEpsilonPhitF.H.

Referenced by correct().

Here is the caller graph for this function:

◆ DepsilonEff()

template<class BasicTurbulenceModel>
tmp< volScalarField > DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon (LUU:Eq. 4).

Definition at line 283 of file kEpsilonPhitF.H.

Referenced by correct().

Here is the caller graph for this function:

◆ DphitEff()

template<class BasicTurbulenceModel>
tmp< volScalarField > DphitEff ( ) const
inline

Return the effective diffusivity for phit (LUU:Eq. 17).

Definition at line 298 of file kEpsilonPhitF.H.

Referenced by correct().

Here is the caller graph for this function:

◆ k()

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > k ( ) const
inlinevirtual

Return the turbulent kinetic energy field.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 314 of file kEpsilonPhitF.H.

◆ epsilon()

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > epsilon ( ) const
inlinevirtual

Return the turbulent kinetic energy dissipation rate field.

Definition at line 322 of file kEpsilonPhitF.H.

◆ phit()

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > phit ( ) const
inlinevirtual

Return the normalised wall-normal fluctuating velocity scale field.

Definition at line 330 of file kEpsilonPhitF.H.

◆ f()

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > f ( ) const
inlinevirtual

Return the elliptic relaxation factor field.

Definition at line 338 of file kEpsilonPhitF.H.

◆ correct()

Member Data Documentation

◆ includeNu_

template<class BasicTurbulenceModel>
Switch includeNu_
protected

Definition at line 148 of file kEpsilonPhitF.H.

Referenced by kEpsilonPhitF(), and read().

◆ Cmu_

template<class BasicTurbulenceModel>
dimensionedScalar Cmu_
protected

Definition at line 152 of file kEpsilonPhitF.H.

Referenced by correctNut(), kEpsilonPhitF(), and read().

◆ Ceps1a_

template<class BasicTurbulenceModel>
dimensionedScalar Ceps1a_
protected

Definition at line 153 of file kEpsilonPhitF.H.

Referenced by correct(), kEpsilonPhitF(), and read().

◆ Ceps1b_

template<class BasicTurbulenceModel>
dimensionedScalar Ceps1b_
protected

Definition at line 154 of file kEpsilonPhitF.H.

Referenced by correct(), kEpsilonPhitF(), and read().

◆ Ceps1c_

template<class BasicTurbulenceModel>
dimensionedScalar Ceps1c_
protected

Definition at line 155 of file kEpsilonPhitF.H.

Referenced by correct(), kEpsilonPhitF(), and read().

◆ Ceps2_

template<class BasicTurbulenceModel>
dimensionedScalar Ceps2_
protected

Definition at line 156 of file kEpsilonPhitF.H.

Referenced by correct(), kEpsilonPhitF(), and read().

◆ Cf1_

template<class BasicTurbulenceModel>
dimensionedScalar Cf1_
protected

Definition at line 157 of file kEpsilonPhitF.H.

Referenced by correct(), kEpsilonPhitF(), and read().

◆ Cf2_

template<class BasicTurbulenceModel>
dimensionedScalar Cf2_
protected

Definition at line 158 of file kEpsilonPhitF.H.

Referenced by correct(), kEpsilonPhitF(), and read().

◆ CL_

template<class BasicTurbulenceModel>
dimensionedScalar CL_
protected

Definition at line 159 of file kEpsilonPhitF.H.

Referenced by kEpsilonPhitF(), Ls(), and read().

◆ Ceta_

template<class BasicTurbulenceModel>
dimensionedScalar Ceta_
protected

Definition at line 160 of file kEpsilonPhitF.H.

Referenced by kEpsilonPhitF(), Ls(), and read().

◆ CT_

template<class BasicTurbulenceModel>
dimensionedScalar CT_
protected

Definition at line 161 of file kEpsilonPhitF.H.

Referenced by kEpsilonPhitF(), read(), and Ts().

◆ sigmaK_

template<class BasicTurbulenceModel>
dimensionedScalar sigmaK_
protected

Definition at line 162 of file kEpsilonPhitF.H.

Referenced by kEpsilonPhitF(), and read().

◆ sigmaEps_

template<class BasicTurbulenceModel>
dimensionedScalar sigmaEps_
protected

Definition at line 163 of file kEpsilonPhitF.H.

Referenced by kEpsilonPhitF(), and read().

◆ sigmaPhit_

template<class BasicTurbulenceModel>
dimensionedScalar sigmaPhit_
protected

Definition at line 164 of file kEpsilonPhitF.H.

Referenced by correct(), kEpsilonPhitF(), and read().

◆ k_

template<class BasicTurbulenceModel>
volScalarField k_
protected

Turbulent kinetic energy [m2/s2].

Definition at line 172 of file kEpsilonPhitF.H.

Referenced by correct(), correctNut(), kEpsilonPhitF(), Ls(), and Ts().

◆ epsilon_

template<class BasicTurbulenceModel>
volScalarField epsilon_
protected

Turbulent kinetic energy dissipation rate [m2/s3].

Definition at line 177 of file kEpsilonPhitF.H.

Referenced by correct(), kEpsilonPhitF(), Ls(), and Ts().

◆ phit_

template<class BasicTurbulenceModel>
volScalarField phit_
protected

Normalised wall-normal fluctuating velocity scale [-].

Definition at line 182 of file kEpsilonPhitF.H.

Referenced by correct(), correctNut(), and kEpsilonPhitF().

◆ f_

template<class BasicTurbulenceModel>
volScalarField f_
protected

Elliptic relaxation factor [1/s].

Definition at line 187 of file kEpsilonPhitF.H.

Referenced by correct(), and kEpsilonPhitF().

◆ T_

template<class BasicTurbulenceModel>
volScalarField T_
protected

Turbulent time scale [s].

Definition at line 192 of file kEpsilonPhitF.H.

Referenced by correct(), correctNut(), and kEpsilonPhitF().

◆ phitMin_

template<class BasicTurbulenceModel>
dimensionedScalar phitMin_
protected

Definition at line 197 of file kEpsilonPhitF.H.

Referenced by correct(), and kEpsilonPhitF().

◆ fMin_

template<class BasicTurbulenceModel>
dimensionedScalar fMin_
protected

Definition at line 198 of file kEpsilonPhitF.H.

Referenced by correct(), and kEpsilonPhitF().

◆ TMin_

template<class BasicTurbulenceModel>
dimensionedScalar TMin_
protected

Definition at line 199 of file kEpsilonPhitF.H.

Referenced by correct(), and kEpsilonPhitF().

◆ L2Min_

template<class BasicTurbulenceModel>
dimensionedScalar L2Min_
protected

Definition at line 200 of file kEpsilonPhitF.H.

Referenced by correct(), and kEpsilonPhitF().


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