v2f< BasicMomentumTransportModel > 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< BasicMomentumTransportModel >:
Collaboration diagram for v2f< BasicMomentumTransportModel >:

Public Types

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

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 viscosity &viscosity, 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< volScalarFieldomega () const
 Return the turbulence specific 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< BasicMomentumTransportModel > >
 eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
 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< BasicMomentumTransportModel >
 linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
 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...
 
- Public Member Functions inherited from v2fBase
 TypeName ("v2fBase")
 Runtime type information. More...
 
 v2fBase ()
 
virtual ~v2fBase ()
 Destructor. More...
 

Protected Member Functions

tmp< volScalarFieldboundEpsilon ()
 Bound epsilon and return Cmu*sqr(k) for nut. More...
 
virtual void correctNut ()
 Correct the eddy-viscosity nut. More...
 
tmp< volScalarFieldTs () const
 Return time scale, Ts. More...
 
tmp< volScalarFieldLs () const
 Return length scale, Ls. 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< BasicMomentumTransportModel > >
volScalarField nut_
 

Detailed Description

template<class BasicMomentumTransportModel>
class Foam::RASModels::v2f< BasicMomentumTransportModel >

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 BasicMomentumTransportModel::alphaField alphaField

Definition at line 179 of file v2f.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 180 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 viscosity viscosity,
const word type = typeName 
)

◆ ~v2f()

virtual ~v2f ( )
inlinevirtual

Destructor.

Definition at line 203 of file v2f.H.

Member Function Documentation

◆ boundEpsilon()

tmp< volScalarField > boundEpsilon
protected

Bound epsilon and return Cmu*sqr(k) for nut.

Definition at line 41 of file v2f.C.

References Foam::max(), and Foam::sqr().

Referenced by v2f< BasicMomentumTransportModel >::v2f().

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

◆ correctNut()

void correctNut
protectedvirtual

Correct the eddy-viscosity nut.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 66 of file v2f.C.

References Foam::min(), and dictionary::New().

Here is the call graph for this function:

◆ Ts()

tmp< volScalarField > Ts
protected

Return time scale, Ts.

Definition at line 50 of file v2f.C.

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

Here is the call graph for this function:

◆ Ls()

tmp< volScalarField > Ls
protected

Return length scale, Ls.

Definition at line 57 of file v2f.C.

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

Here is the call graph for this function:

◆ TypeName()

TypeName ( "v2f< BasicMomentumTransportModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Read RASProperties dictionary.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 257 of file v2f.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 v2f.H.

References GeometricField< Type, PatchField, GeoMesh >::New(), and eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_.

Here is the call graph for this function:

◆ DepsilonEff()

tmp<volScalarField> DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 223 of file v2f.H.

References GeometricField< Type, PatchField, GeoMesh >::New(), and eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_.

Here is the call graph for this function:

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 233 of file v2f.H.

References v2f< BasicMomentumTransportModel >::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 239 of file v2f.H.

References v2f< BasicMomentumTransportModel >::epsilon_.

◆ omega()

virtual tmp<volScalarField> omega ( ) const
inlinevirtual

Return the turbulence specific dissipation rate.

Definition at line 245 of file v2f.H.

References v2f< BasicMomentumTransportModel >::Cmu_, v2f< BasicMomentumTransportModel >::epsilon_, v2f< BasicMomentumTransportModel >::k_, and GeometricField< Type, PatchField, GeoMesh >::New().

Here is the call graph for this function:

◆ v2()

virtual tmp<volScalarField> v2 ( ) const
inlinevirtual

Return turbulence stress normal to streamlines.

Implements v2fBase.

Definition at line 255 of file v2f.H.

References v2f< BasicMomentumTransportModel >::v2_.

◆ f()

virtual tmp<volScalarField> f ( ) const
inlinevirtual

Return the damping function.

Implements v2fBase.

Definition at line 261 of file v2f.H.

References v2f< BasicMomentumTransportModel >::f_.

◆ correct()

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 129 of file v2f.H.

Referenced by v2f< BasicMomentumTransportModel >::omega().

◆ 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_

◆ epsilon_

volScalarField epsilon_
protected

Turbulence dissipation.

Definition at line 147 of file v2f.H.

Referenced by v2f< BasicMomentumTransportModel >::epsilon(), and v2f< BasicMomentumTransportModel >::omega().

◆ v2_

volScalarField v2_
protected

Turbulence stress normal to streamlines.

Definition at line 150 of file v2f.H.

Referenced by v2f< BasicMomentumTransportModel >::v2(), and v2f< BasicMomentumTransportModel >::v2f().

◆ f_

volScalarField f_
protected

Damping function.

Definition at line 153 of file v2f.H.

Referenced by v2f< BasicMomentumTransportModel >::f(), and v2f< BasicMomentumTransportModel >::v2f().

◆ v2Min_

dimensionedScalar v2Min_
protected

Definition at line 158 of file v2f.H.

Referenced by v2f< BasicMomentumTransportModel >::v2f().

◆ fMin_

dimensionedScalar fMin_
protected

Definition at line 159 of file v2f.H.

Referenced by v2f< BasicMomentumTransportModel >::v2f().


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