A one-equation (turbulent kinetic energy k
) turbulence closure model for incompressible and compressible geophysical applications.
More...
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 ("kL") | |
Runtime type information. More... | |
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. More... | |
virtual | ~kL ()=default |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. More... | |
tmp< volScalarField > | DkEff () const |
Return the effective diffusivity for k. More... | |
virtual tmp< volScalarField > | k () const |
Return the turbulent kinetic energy field. 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 ()=default |
Destructor. More... | |
virtual tmp< volScalarField > | nut () const |
Return the turbulence viscosity. More... | |
virtual tmp< scalarField > | nut (const label patchi) const |
Return the turbulence viscosity on patch. More... | |
virtual tmp< volSymmTensorField > | R () 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 ()=default |
Destructor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
Return the effective stress tensor based on a given velocity field. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (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 ()=default |
Destructor. More... | |
const dimensionedScalar & | kMin () const noexcept |
Return the lower allowable limit for k (default: SMALL) More... | |
const dimensionedScalar & | epsilonMin () const noexcept |
Return the lower allowable limit for epsilon (default: SMALL) More... | |
const dimensionedScalar & | omegaMin () const noexcept |
Return the lower allowable limit for omega (default: SMALL) More... | |
dimensionedScalar & | kMin () noexcept |
Allow kMin to be changed. More... | |
dimensionedScalar & | epsilonMin () noexcept |
Allow epsilonMin to be changed. More... | |
dimensionedScalar & | omegaMin () noexcept |
Allow omegaMin to be changed. More... | |
virtual const dictionary & | coeffDict () const |
Const access to the coefficients dictionary. More... | |
virtual tmp< volScalarField > | nuEff () const |
Return the effective viscosity. More... | |
virtual tmp< scalarField > | nuEff (const label patchi) const |
Return the effective viscosity on patch. More... | |
virtual tmp< volScalarField > | omega () const |
Return the specific dissipation rate. More... | |
Public Member Functions inherited from RASModelBase | |
ClassName ("RASModelBase") | |
Runtime type information. More... | |
RASModelBase () noexcept=default | |
Constructor. More... | |
virtual | ~RASModelBase ()=default |
Destructor. More... | |
Protected Member Functions | |
virtual void | correctNut () |
Correct the turbulence viscosity. More... | |
virtual tmp< fvScalarMatrix > | kSource () const |
Add explicit source for turbulent kinetic energy. More... | |
Protected Member Functions inherited from RASModel< BasicTurbulenceModel > | |
virtual void | printCoeffs (const word &type) |
Print model coefficients. More... | |
RASModel (const RASModel &)=delete | |
No copy construct. More... | |
void | operator= (const RASModel &)=delete |
No copy assignment. More... | |
Protected Attributes | |
dimensionedScalar | kappa_ |
von Karman constant More... | |
dimensionedScalar | sigmak_ |
Empirical model coefficient. More... | |
dimensionedScalar | beta_ |
Thermal expansion coefficient [1/K]. More... | |
dimensionedScalar | Cmu0_ |
Empirical model coefficient. More... | |
dimensionedScalar | Lmax_ |
Maximum mixing-length scalar [m]. More... | |
dimensionedScalar | CbStable_ |
Stable stratification constant. More... | |
dimensionedScalar | CbUnstable_ |
Unstable stratification constant. More... | |
volScalarField | k_ |
Turbulent kinetic energy [m2/s2]. More... | |
volScalarField | L_ |
Characteristic length scale [m]. More... | |
volScalarField | Rt_ |
Turbulent Richardson number [-]. More... | |
const uniformDimensionedVectorField & | g_ |
Gravitational acceleration [m2/s2]. More... | |
const volScalarField & | y_ |
Wall distance. More... | |
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< RASModel > | New (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... | |
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
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] |
readFields
function object.typedef BasicTurbulenceModel::alphaField alphaField |
typedef BasicTurbulenceModel::transportModel transportModel |
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.
Definition at line 241 of file kL.C.
References Foam::bound(), and Foam::type().
|
virtualdefault |
Destructor.
|
protectedvirtual |
Correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 216 of file kL.C.
References optionList::correct(), options::New(), and Foam::sqrt().
|
protectedvirtual |
Add explicit source for turbulent kinetic energy.
Definition at line 227 of file kL.C.
References Foam::dimTime, Foam::dimVolume, and tmp< T >::New().
TypeName | ( | "kL< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 383 of file kL.C.
References Foam::read().
|
inline |
Return the effective diffusivity for k.
Definition at line 414 of file kL.H.
References tmp< T >::New(), nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inlinevirtual |
Return the turbulent kinetic energy field.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 426 of file kL.H.
References kL< BasicTurbulenceModel >::k_.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 403 of file kL.C.
References Foam::fvc::absolute(), Foam::constant::atomic::alpha, Foam::bound(), tmp< T >::clear(), correct, tmp< T >::cref(), Foam::fvm::ddt(), Foam::devSymm(), Foam::fvm::div(), Foam::fvc::div(), divU, epsilon, fvOptions, Foam::constant::universal::G, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), options::New(), nut, phi, tmp< T >::ref(), rho, Foam::solve(), Foam::sqr(), Foam::fvm::SuSp(), Foam::T(), and U.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Turbulent kinetic energy [m2/s2].
Definition at line 328 of file kL.H.
Referenced by kL< BasicTurbulenceModel >::k().
|
protected |
|
protected |
|
protected |
|
protected |
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.