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

A one-equation (turbulent kinetic energy k) turbulence closure model for incompressible and compressible geophysical applications. More...

#include <kL.H>

Inheritance diagram for kL< BasicTurbulenceModel >:
Collaboration diagram for kL< 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 ("kL")
 Runtime type information.
 kL (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 ~kL ()=default
 Destructor.
virtual bool read ()
 Re-read model coefficients if they have changed.
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k.
virtual tmp< volScalarFieldk () const
 Return the turbulent kinetic energy field.
virtual void correct ()
 Solve the turbulence equations and correct the turbulence 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 ()
 Correct the turbulence viscosity.
virtual tmp< fvScalarMatrixkSource () const
 Add explicit source for turbulent kinetic energy.

Protected Attributes

dimensionedScalar kappa_
 von Karman constant
dimensionedScalar sigmak_
 Empirical model coefficient.
dimensionedScalar beta_
 Thermal expansion coefficient [1/K].
dimensionedScalar Cmu0_
 Empirical model coefficient.
dimensionedScalar Lmax_
 Maximum mixing-length scalar [m].
dimensionedScalar CbStable_
 Stable stratification constant.
dimensionedScalar CbUnstable_
 Unstable stratification constant.
volScalarField k_
 Turbulent kinetic energy [m2/s2].
volScalarField L_
 Characteristic length scale [m].
volScalarField Rt_
 Turbulent Richardson number [-].
const uniformDimensionedVectorFieldg_
 Gravitational acceleration [m2/s2].
const volScalarFieldy_
 Wall distance.
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_

Detailed Description

template<class BasicTurbulenceModel>
class Foam::RASModels::kL< BasicTurbulenceModel >

A one-equation (turbulent kinetic energy k) turbulence closure model for incompressible and compressible geophysical applications.

Turbulent kinetic energy (k) is computed with a transport equation and the turbulent length scale (L) is computed with an algebraic expression which depends on the local stratification.

References:

    Standard model (tag:A):
        Axell, L. B., & Liungman, O. (2001).
        A one-equation turbulence model for geophysical applications:
        comparison with data and the k−ε model.
        Environmental Fluid Mechanics, 1(1), 71-106.
        DOI:10.1023/A:1011560202388

    Canopy-related models (tag:W):
        Wilson, J. D., & Flesch, T. K. (1999).
        Wind and remnant tree sway in forest cutblocks.
        III. A windflow model to diagnose spatial variation.
        Agricultural and Forest Meteorology, 93(4), 259-282.
        DOI:10.1016/S0168-1923(98)00121-X
Usage
Example by using constant/turbulenceProperties:
RAS
{
    // Mandatory entries
    RASModel        kL;

    // Optional entries
    kLCoeffs
    {
        kappa       <scalar>;
        sigmak      <scalar>;
        beta        <scalar>;
        Cmu0        <scalar>;
        Lmax        <scalar>;
        CbStable    <scalar>;
        CbUnstable  <scalar>;
    }

    // Inherited entries
    ...
}

where the entries mean:

Property Description Type Reqd Deflt
RASModel Type name: kL word yes -
kappa von Karman constant scalar no 0.41
sigmak Empirical model coefficient scalar no 1.0
beta Thermal expansion coefficient [1/K] scalar no 3.3e-3
Cmu0 Empirical model coefficient scalar no 0.556
Lmax Maximum mixing-length scale [m] scalar no GREAT
CbStable Stable stratification constant scalar no 0.25
CbUnstable Unstable stratification constant scalar no 0.35

The inherited entries are elaborated in:


Input fields (mandatory)

k :Turbulent kinetic energy [m2/s2]
T :Potential temperature [K]


Input fields (optional)

canopyHeight :Canopy height [m]
Cd :Plant canopy drag coefficient [-]
LAD :Leaf area density [1/m]
Rt :Turbulent Richardson number [-]
L :Characteristic length scale [m]
Note
  • Optional input fields can/should be input by using readFields function object.
Source files

Definition at line 218 of file kL.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 369 of file kL.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 370 of file kL.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 371 of file kL.H.

Constructor & Destructor Documentation

◆ kL()

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

◆ ~kL()

template<class BasicTurbulenceModel>
virtual ~kL ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

template<class BasicTurbulenceModel>
void correctNut ( )
protectedvirtual

Correct the turbulence viscosity.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 205 of file kL.C.

References optionList::correct(), k_, L_, options::New(), eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_, and Foam::sqrt().

Referenced by correct().

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

◆ kSource()

template<class BasicTurbulenceModel>
tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Add explicit source for turbulent kinetic energy.

Definition at line 216 of file kL.C.

References Foam::dimTime, Foam::dimVolume, k_, and tmp< T >::New().

Referenced by correct().

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

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "kL< 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 372 of file kL.C.

References beta_, CbStable_, CbUnstable_, Cmu0_, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), kappa_, Lmax_, read(), and sigmak_.

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.

Definition at line 414 of file kL.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 426 of file kL.H.

◆ correct()

Member Data Documentation

◆ kappa_

template<class BasicTurbulenceModel>
dimensionedScalar kappa_
protected

von Karman constant

Definition at line 290 of file kL.H.

Referenced by kL(), and read().

◆ sigmak_

template<class BasicTurbulenceModel>
dimensionedScalar sigmak_
protected

Empirical model coefficient.

Definition at line 295 of file kL.H.

Referenced by kL(), and read().

◆ beta_

template<class BasicTurbulenceModel>
dimensionedScalar beta_
protected

Thermal expansion coefficient [1/K].

Definition at line 300 of file kL.H.

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

◆ Cmu0_

template<class BasicTurbulenceModel>
dimensionedScalar Cmu0_
protected

Empirical model coefficient.

Definition at line 305 of file kL.H.

Referenced by kL(), and read().

◆ Lmax_

template<class BasicTurbulenceModel>
dimensionedScalar Lmax_
protected

Maximum mixing-length scalar [m].

Definition at line 310 of file kL.H.

Referenced by kL(), and read().

◆ CbStable_

template<class BasicTurbulenceModel>
dimensionedScalar CbStable_
protected

Stable stratification constant.

Definition at line 315 of file kL.H.

Referenced by kL(), and read().

◆ CbUnstable_

template<class BasicTurbulenceModel>
dimensionedScalar CbUnstable_
protected

Unstable stratification constant.

Definition at line 320 of file kL.H.

Referenced by kL(), and read().

◆ k_

template<class BasicTurbulenceModel>
volScalarField k_
protected

Turbulent kinetic energy [m2/s2].

Definition at line 328 of file kL.H.

Referenced by correct(), correctNut(), kL(), and kSource().

◆ L_

template<class BasicTurbulenceModel>
volScalarField L_
protected

Characteristic length scale [m].

Definition at line 333 of file kL.H.

Referenced by correctNut(), and kL().

◆ Rt_

template<class BasicTurbulenceModel>
volScalarField Rt_
protected

Turbulent Richardson number [-].

Definition at line 338 of file kL.H.

Referenced by correct(), and kL().

◆ g_

template<class BasicTurbulenceModel>
const uniformDimensionedVectorField& g_
protected

Gravitational acceleration [m2/s2].

Definition at line 343 of file kL.H.

Referenced by correct(), and kL().

◆ y_

template<class BasicTurbulenceModel>
const volScalarField& y_
protected

Wall distance.

Note: different to wall distance in parent RASModel which is for near-wall cells only

Definition at line 351 of file kL.H.

Referenced by kL().


The documentation for this class was generated from the following files:
  • src/atmosphericModels/turbulenceModels/RAS/kL/kL.H
  • src/atmosphericModels/turbulenceModels/RAS/kL/kL.C