OpenFOAM logo
The Open Source CFD Toolbox
  Source Guide OpenCFD Solutions Contact OpenFOAM

PDRkEpsilon Class Reference

Standard k-epsilon turbulence model with additional source terms corresponding to PDR basic drag model (basic.H). More...

Inheritance diagram for PDRkEpsilon:
Collaboration diagram for PDRkEpsilon:

List of all members.


Public Member Functions

 TypeName ("PDRkEpsilon")
 Runtime type information.
 PDRkEpsilon (const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermophysicalModel)
 Construct from components.
virtual ~PDRkEpsilon ()
 Destructor.
tmp< volScalarFieldmut () const
 Return the turbulence viscosity.
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k.
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon.
tmp< volScalarFieldalphaEff () const
 Return the effective turbulent thermal diffusivity.
tmp< volScalarFieldk () const
 Return the turbulence kinetic energy.
tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate.
tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor.
tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor including the laminar stress.
tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation.
void correct ()
 Solve the turbulence equations and correct the turbulence viscosity.
bool read ()
 Read turbulenceProperties dictionary.

Detailed Description

Standard k-epsilon turbulence model with additional source terms corresponding to PDR basic drag model (basic.H).

The default model coefficients correspond to the following:

    kEpsilonCoeffs
    {
        Cmu         0.09;
        C1          1.44;
        C2          1.92;
        C3          -0.33;  // only for compressible
        sigmak      1.0;    // only for compressible
        sigmaEps    1.3;
        Prt         1.0;    // only for compressible
    }

The turbulence source term $ G_{R} $ appears in the $ \kappa-\epsilon $ equation for the generation of turbulence due to interaction with unresolved obstacles.

In the $ \epsilon $ equation $ C_{1} G_{R} $ is added as a source term.

In the $ \kappa $ equation $ G_{R} $ is added as a source term.

Source files

Definition at line 73 of file PDRkEpsilon.H.


Constructor & Destructor Documentation

PDRkEpsilon ( const volScalarField rho,
const volVectorField U,
const surfaceScalarField phi,
const basicThermo thermophysicalModel 
)

Construct from components.

virtual ~PDRkEpsilon (  )  [inline, virtual]

Destructor.

Definition at line 115 of file PDRkEpsilon.H.


Member Function Documentation

TypeName ( "PDRkEpsilon"   ) 

Runtime type information.

tmp<volScalarField> mut (  )  const [inline, virtual]

Return the turbulence viscosity.

Implements RASModel.

Definition at line 121 of file PDRkEpsilon.H.

Referenced by PDRkEpsilon::DkEff().

Here is the caller graph for this function:

tmp<volScalarField> DkEff (  )  const [inline]

Return the effective diffusivity for k.

Definition at line 127 of file PDRkEpsilon.H.

References PDRkEpsilon::mut().

Here is the call graph for this function:

tmp<volScalarField> DepsilonEff (  )  const [inline]

Return the effective diffusivity for epsilon.

Definition at line 136 of file PDRkEpsilon.H.

tmp<volScalarField> alphaEff (  )  const [inline, virtual]

Return the effective turbulent thermal diffusivity.

Implements RASModel.

Definition at line 145 of file PDRkEpsilon.H.

tmp<volScalarField> k (  )  const [inline, virtual]

Return the turbulence kinetic energy.

Implements RASModel.

Definition at line 154 of file PDRkEpsilon.H.

References turbulenceModel::alpha().

Here is the call graph for this function:

tmp<volScalarField> epsilon (  )  const [inline, virtual]

Return the turbulence kinetic energy dissipation rate.

Implements RASModel.

Definition at line 160 of file PDRkEpsilon.H.

tmp<volSymmTensorField> R (  )  const [virtual]

Return the Reynolds stress tensor.

Implements RASModel.

tmp<volSymmTensorField> devRhoReff (  )  const [virtual]

Return the effective stress tensor including the laminar stress.

Implements RASModel.

tmp<fvVectorMatrix> divDevRhoReff ( volVectorField U  )  const [virtual]

Return the source term for the momentum equation.

Implements RASModel.

void correct (  )  [virtual]

Solve the turbulence equations and correct the turbulence viscosity.

Implements RASModel.

bool read (  )  [virtual]

Read turbulenceProperties dictionary.

Implements RASModel.


The documentation for this class was generated from the following file:
  • applications/solvers/combustion/PDRFoam/PDRModels/turbulence/PDRkEpsilon/PDRkEpsilon.H
Copyright © 2000-2009 OpenCFD Ltd