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

Lien and Kalitzin's v2-f turbulence model for incompressible and compressible flows, with a limit imposed on the turbulent viscosity given by Davidson et al. More...

Inheritance diagram for v2f< BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for v2f< 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 ("v2f")
 Runtime type information. More...
 
 v2f (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 ~v2f ()
 Destructor. More...
 
virtual bool read ()
 Read RASProperties dictionary. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
tmp< volScalarFieldDepsilonEff () 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 tmp< volScalarFieldv2 () const
 Return turbulence stress normal to streamlines. More...
 
virtual tmp< volScalarFieldf () const
 Return the damping function. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence 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 ()
 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 ()
 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 ()
 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 v2fBase
 TypeName ("v2fBase")
 Runtime type information. More...
 
 v2fBase ()
 
virtual ~v2fBase ()
 Destructor. More...
 

Protected Member Functions

virtual void correctNut ()
 
tmp< volScalarFieldTs () const
 Return time scale, Ts. More...
 
tmp< volScalarFieldLs () const
 Return length scale, Ls. More...
 
- Protected Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 

Protected Attributes

dimensionedScalar Cmu_
 
dimensionedScalar CmuKEps_
 
dimensionedScalar C1_
 
dimensionedScalar C2_
 
dimensionedScalar CL_
 
dimensionedScalar Ceta_
 
dimensionedScalar Ceps2_
 
dimensionedScalar Ceps3_
 
dimensionedScalar sigmaK_
 
dimensionedScalar sigmaEps_
 
volScalarField k_
 Turbulence kinetic energy. More...
 
volScalarField epsilon_
 Turbulence dissipation. More...
 
volScalarField v2_
 Turbulence stress normal to streamlines. More...
 
volScalarField f_
 Damping function. More...
 
dimensionedScalar v2Min_
 
dimensionedScalar fMin_
 
- 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::v2f< BasicTurbulenceModel >

Lien and Kalitzin's v2-f turbulence model for incompressible and compressible flows, with a limit imposed on the turbulent viscosity given by Davidson et al.

The model solves for turbulence kinetic energy k and turbulence dissipation rate epsilon, with additional equations for the turbulence stress normal to streamlines, v2, and elliptic damping function, f.

The variant implemented employs N=6, such that f=0 on walls.

Wall boundary conditions are:

k = kLowReWallFunction epsilon = epsilonWallFunction v2 = v2WallFunction f = fWallFunction

These are applicable to both low- and high-Reynolds number flows.

Inlet values can be approximated by:

v2 = 2/3 k f = zero-gradient

References:

    Lien, F. S., & Kalitzin, G. (2001).
    Computations of transonic flow with the v2f turbulence model.
    International Journal of Heat and Fluid Flow, 22(1), 53-61.

    Davidson, L., Nielsen, P., & Sveningsson, A. (2003).
    Modifications of the v2-f model for computing the flow in a
    3D wall jet.
    Turbulence, Heat and Mass Transfer, 4, 577-584

The default model coefficients are

    v2fCoeffs
    {
        Cmu         0.22;
        CmuKEps     0.09;
        C1          1.4;
        C2          0.3;
        CL          0.23;
        Ceta        70;
        Ceps2       1.9;
        Ceps3       -0.33;
        sigmaEps    1.3;
        sigmaK      1;
    }
Note
If the kLowReWallFunction is employed, a velocity variant of the turbulent viscosity wall function should be used, e.g. nutUWallFunction. Turbulence k variants (nutk...) for this case will not behave correctly.
See also
Foam::RASModels::v2fBase Foam::RASModels::kEpsilon Foam::kLowReWallFunctionFvPatchScalarField Foam::epsilonWallFunctionFvPatchScalarField Foam::v2WallFunctionFvPatchScalarField Foam::fWallFunctionFvPatchScalarField
Source files

Definition at line 117 of file v2f.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 175 of file v2f.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 176 of file v2f.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 177 of file v2f.H.

Constructor & Destructor Documentation

◆ v2f()

v2f ( 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 70 of file v2f.C.

References Foam::bound().

Referenced by v2f< BasicTurbulenceModel >::correctNut().

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

◆ ~v2f()

virtual ~v2f ( )
inlinevirtual

Destructor.

Definition at line 201 of file v2f.H.

References v2f< BasicTurbulenceModel >::read().

Here is the call graph for this function:

Member Function Documentation

◆ correctNut()

void correctNut ( )
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 56 of file v2f.C.

References optionList::correct(), Foam::min(), options::New(), Foam::sqr(), and v2f< BasicTurbulenceModel >::v2f().

Here is the call graph for this function:

◆ Ts()

tmp< volScalarField > Ts ( ) const
protected

Return time scale, Ts.

Definition at line 40 of file v2f.C.

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

Here is the call graph for this function:

◆ Ls()

tmp< volScalarField > Ls ( ) const
protected

Return length scale, Ls.

Definition at line 47 of file v2f.C.

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

Here is the call graph for this function:

◆ TypeName()

TypeName ( "v2f< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Read RASProperties dictionary.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 251 of file v2f.C.

References Foam::read().

Referenced by v2f< BasicTurbulenceModel >::~v2f().

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

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 211 of file v2f.H.

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

◆ DepsilonEff()

tmp<volScalarField> DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 224 of file v2f.H.

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

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 237 of file v2f.H.

References v2f< BasicTurbulenceModel >::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 243 of file v2f.H.

References v2f< BasicTurbulenceModel >::epsilon_.

◆ v2()

virtual tmp<volScalarField> v2 ( ) const
inlinevirtual

Return turbulence stress normal to streamlines.

Implements v2fBase.

Definition at line 249 of file v2f.H.

References v2f< BasicTurbulenceModel >::v2_.

◆ f()

virtual tmp<volScalarField> f ( ) const
inlinevirtual

Return the damping function.

Implements v2fBase.

Definition at line 255 of file v2f.H.

References v2f< BasicTurbulenceModel >::correct(), and v2f< BasicTurbulenceModel >::f_.

Here is the call graph for this function:

◆ correct()

void correct ( )
virtual

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 129 of file v2f.H.

◆ CmuKEps_

dimensionedScalar CmuKEps_
protected

Definition at line 130 of file v2f.H.

◆ C1_

dimensionedScalar C1_
protected

Definition at line 131 of file v2f.H.

◆ C2_

dimensionedScalar C2_
protected

Definition at line 132 of file v2f.H.

◆ CL_

dimensionedScalar CL_
protected

Definition at line 133 of file v2f.H.

◆ Ceta_

dimensionedScalar Ceta_
protected

Definition at line 134 of file v2f.H.

◆ Ceps2_

dimensionedScalar Ceps2_
protected

Definition at line 135 of file v2f.H.

◆ Ceps3_

dimensionedScalar Ceps3_
protected

Definition at line 136 of file v2f.H.

◆ sigmaK_

dimensionedScalar sigmaK_
protected

Definition at line 137 of file v2f.H.

◆ sigmaEps_

dimensionedScalar sigmaEps_
protected

Definition at line 138 of file v2f.H.

◆ k_

volScalarField k_
protected

Turbulence kinetic energy.

Definition at line 144 of file v2f.H.

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

◆ epsilon_

volScalarField epsilon_
protected

Turbulence dissipation.

Definition at line 147 of file v2f.H.

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

◆ v2_

volScalarField v2_
protected

Turbulence stress normal to streamlines.

Definition at line 150 of file v2f.H.

Referenced by v2f< BasicTurbulenceModel >::v2().

◆ f_

volScalarField f_
protected

Damping function.

Definition at line 153 of file v2f.H.

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

◆ v2Min_

dimensionedScalar v2Min_
protected

Definition at line 158 of file v2f.H.

◆ fMin_

dimensionedScalar fMin_
protected

Definition at line 159 of file v2f.H.


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