The Wall-adapting local eddy-viscosity (WALE) SGS model. More...
#include <WALE.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from LESeddyViscosity< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from eddyViscosity< LESModel< 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 ("WALE") | |
| Runtime type information. | |
| WALE (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 | ~WALE ()=default |
| Destructor. | |
| virtual bool | read () |
| Read model coefficients if they have changed. | |
| virtual tmp< volScalarField > | k () const |
| Return SGS kinetic energy. | |
| virtual void | correct () |
| Correct Eddy-Viscosity and related properties. | |
| Public Member Functions inherited from LESeddyViscosity< BasicTurbulenceModel > | |
| LESeddyViscosity (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) | |
| Construct from components. | |
| virtual | ~LESeddyViscosity ()=default |
| Destructor. | |
| Public Member Functions inherited from eddyViscosity< LESModel< 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< volScalarField > | nut () const |
| Return the turbulence viscosity. | |
| virtual tmp< volSymmTensorField > | R () 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< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
| Return the effective stress tensor based on a given velocity field. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. | |
Protected Member Functions | |
| tmp< volSymmTensorField > | Sd (const volTensorField &gradU) const |
| Return the deviatoric symmetric part of the square of the given. | |
| tmp< volScalarField > | k (const volTensorField &gradU) const |
| Return SGS kinetic energy. | |
| virtual void | correctNut () |
| Update the SGS eddy-viscosity. | |
Protected Attributes | |
| dimensionedScalar | Ck_ |
| dimensionedScalar | Cw_ |
| Protected Attributes inherited from eddyViscosity< LESModel< BasicTurbulenceModel > > | |
| volScalarField | nut_ |
The Wall-adapting local eddy-viscosity (WALE) SGS model.
Reference:
Nicoud, F., & Ducros, F. (1999).
Subgrid-scale stress modelling based on the square of the velocity
gradient tensor.
Flow, Turbulence and Combustion, 62(3), 183-200.
The default model coefficients are
WALECoeffs
{
Ck 0.094;
Ce 1.048;e
Cw 0.325;
}
| typedef BasicTurbulenceModel::alphaField alphaField |
| typedef BasicTurbulenceModel::rhoField rhoField |
| typedef BasicTurbulenceModel::transportModel transportModel |
| WALE | ( | 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 ) |
|
virtualdefault |
Destructor.
|
protected |
Return the deviatoric symmetric part of the square of the given.
velocity gradient field
Definition at line 35 of file WALE.C.
References Foam::devSymm().
Referenced by k().


|
protected |
Return SGS kinetic energy.
calculated from the given velocity gradient
Definition at line 45 of file WALE.C.
References Ck_, Cw_, delta, IOobject::groupName(), Foam::magSqr(), Foam::pow(), Foam::pow3(), Sd(), Foam::sqr(), and Foam::symm().

|
protectedvirtual |
Update the SGS eddy-viscosity.
Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.
Definition at line 85 of file WALE.C.
References Ck_, optionList::correct(), delta, Foam::fvc::grad(), k, options::New(), eddyViscosity< LESModel< BasicTurbulenceModel > >::nut_, and Foam::sqrt().
Referenced by correct().


| TypeName | ( | "WALE< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Reimplemented from LESeddyViscosity< BasicTurbulenceModel >.
Definition at line 152 of file WALE.C.
References Ck_, Cw_, and LESeddyViscosity< BasicTurbulenceModel >::read().

|
inlinevirtual |
Return SGS kinetic energy.
Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.
Definition at line 169 of file WALE.H.
Referenced by WALE< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >::k().

|
virtual |
Correct Eddy-Viscosity and related properties.
Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.
Definition at line 167 of file WALE.C.
References linearViscousStress< BasicTurbulenceModel >::correct(), and correctNut().

|
protected |
|
protected |