|
| 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...
|
|
| 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 bool | read ()=0 |
| Re-read model coefficients if they have changed. 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< volScalarField > | k () const=0 |
| Return the turbulence kinetic energy. More...
|
|
virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. More...
|
|
virtual void | validate () |
| Validate the turbulence fields after construction. More...
|
|
virtual void | correct ()=0 |
| Solve the turbulence equations and correct the turbulence viscosity. More...
|
|
| 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 bool | read ()=0 |
| Re-read model coefficients if they have changed. 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...
|
|
virtual void | correct ()=0 |
| Solve the turbulence equations and correct the turbulence viscosity. More...
|
|
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] |
plantCd | : | Plant canopy drag coefficient [-] |
leafAreaDensity | : | 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 222 of file kL.H.