Dynamic one equation eddy-viscosity model. More...
Public Types | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
typedef BasicMomentumTransportModel::transportModel | transportModel |
Public Types inherited from LESeddyViscosity< BasicMomentumTransportModel > | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
typedef BasicMomentumTransportModel::transportModel | transportModel |
Public Types inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > > | |
typedef LESModel< BasicMomentumTransportModel > ::alphaField | alphaField |
typedef LESModel< BasicMomentumTransportModel > ::rhoField | rhoField |
typedef LESModel< BasicMomentumTransportModel > ::transportModel | transportModel |
Public Types inherited from linearViscousStress< LESModel< BasicMomentumTransportModel > > | |
typedef LESModel< BasicMomentumTransportModel > ::alphaField | alphaField |
typedef LESModel< BasicMomentumTransportModel > ::rhoField | rhoField |
typedef LESModel< BasicMomentumTransportModel > ::transportModel | transportModel |
Public Types inherited from LESModel< BasicMomentumTransportModel > | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
typedef BasicMomentumTransportModel::transportModel | transportModel |
Public Member Functions | |
TypeName ("dynamicKEqn") | |
Runtime type information. More... | |
dynamicKEqn (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName) | |
Construct from components. More... | |
dynamicKEqn (const dynamicKEqn &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~dynamicKEqn () |
Destructor. More... | |
virtual bool | read () |
Read model coefficients if they have changed. More... | |
virtual tmp< volScalarField > | k () const |
Return SGS kinetic energy. More... | |
virtual tmp< volScalarField > | epsilon () const |
Return sub-grid disipation rate. More... | |
tmp< volScalarField > | DkEff () const |
Return the effective diffusivity for k. More... | |
virtual void | correct () |
Correct Eddy-Viscosity and related properties. More... | |
void | operator= (const dynamicKEqn &)=delete |
Disallow default bitwise assignment. More... | |
Public Member Functions inherited from LESeddyViscosity< BasicMomentumTransportModel > | |
LESeddyViscosity (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) | |
Construct from components. More... | |
LESeddyViscosity (const LESeddyViscosity &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~LESeddyViscosity () |
Destructor. More... | |
void | operator= (const LESeddyViscosity &)=delete |
Disallow default bitwise assignment. More... | |
Public Member Functions inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > > | |
eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) | |
Construct from components. More... | |
virtual | ~eddyViscosity () |
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 > | sigma () const |
Return the Reynolds stress tensor [m^2/s^2]. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
Public Member Functions inherited from linearViscousStress< LESModel< BasicMomentumTransportModel > > | |
linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) | |
Construct from components. More... | |
virtual | ~linearViscousStress () |
Destructor. More... | |
virtual tmp< volSymmTensorField > | devTau () const |
Return the effective stress tensor. More... | |
virtual tmp< fvVectorMatrix > | divDevTau (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevTau (const volScalarField &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
Public Member Functions inherited from LESModel< BasicMomentumTransportModel > | |
TypeName ("LES") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, LESModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport),(alpha, rho, U, alphaRhoPhi, phi, transport)) | |
LESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) | |
Construct from components. More... | |
LESModel (const LESModel &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~LESModel () |
Destructor. More... | |
virtual const dictionary & | coeffDict () const |
Const access to the coefficients dictionary. More... | |
const dimensionedScalar & | kMin () const |
Return the lower allowable limit for k (default: small) More... | |
dimensionedScalar & | kMin () |
Allow kMin to be changed. More... | |
const volScalarField & | delta () const |
Access function to filter width. 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... | |
void | operator= (const LESModel &)=delete |
Disallow default bitwise assignment. More... | |
Protected Member Functions | |
tmp< volScalarField > | KK () const |
Return the KK field obtained by filtering the velocity field U. More... | |
volScalarField | Ck (const volSymmTensorField &D, const volScalarField &KK) const |
Calculate Ck by filtering the velocity field U. More... | |
volScalarField | Ce (const volSymmTensorField &D, const volScalarField &KK) const |
Calculate Ce by filtering the velocity field U. More... | |
volScalarField | Ce () const |
void | correctNut (const volSymmTensorField &D, const volScalarField &KK) |
Update sub-grid eddy-viscosity. More... | |
virtual void | correctNut () |
virtual tmp< fvScalarMatrix > | kSource () const |
Protected Member Functions inherited from LESModel< BasicMomentumTransportModel > | |
virtual void | printCoeffs (const word &type) |
Print model coefficients. More... | |
Protected Attributes | |
volScalarField | k_ |
simpleFilter | simpleFilter_ |
autoPtr< LESfilter > | filterPtr_ |
LESfilter & | filter_ |
Protected Attributes inherited from LESeddyViscosity< BasicMomentumTransportModel > | |
dimensionedScalar | Ce_ |
Protected Attributes inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > > | |
volScalarField | nut_ |
Protected Attributes inherited from LESModel< BasicMomentumTransportModel > | |
dictionary | LESDict_ |
LES 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... | |
autoPtr< Foam::LESdelta > | delta_ |
Run-time selectable delta model. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from LESModel< BasicMomentumTransportModel > | |
static autoPtr< LESModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) |
Return a reference to the selected LES model. More... | |
Dynamic one equation eddy-viscosity model.
Eddy viscosity SGS model using a modeled balance equation to simulate the behaviour of k in which a dynamic procedure is applied to evaluate the coefficients.
Reference:
Kim, W and Menon, S. (1995). A new dynamic one-equation subgrid-scale model for large eddy simulation. In 33rd Aerospace Sciences Meeting and Exhibit, Reno, NV, 1995.
There are no default model coefficients but the filter used for KK must be supplied, e.g.
dynamicKEqnCoeffs { filter simple; }
Definition at line 73 of file dynamicKEqn.H.
typedef BasicMomentumTransportModel::alphaField alphaField |
Definition at line 128 of file dynamicKEqn.H.
typedef BasicMomentumTransportModel::rhoField rhoField |
Definition at line 129 of file dynamicKEqn.H.
typedef BasicMomentumTransportModel::transportModel transportModel |
Definition at line 130 of file dynamicKEqn.H.
dynamicKEqn | ( | const alphaField & | alpha, |
const rhoField & | rho, | ||
const volVectorField & | U, | ||
const surfaceScalarField & | alphaRhoPhi, | ||
const surfaceScalarField & | phi, | ||
const transportModel & | transport, | ||
const word & | type = typeName |
||
) |
Construct from components.
Definition at line 147 of file dynamicKEqn.C.
References Foam::bound().
Referenced by dynamicKEqn< BasicMomentumTransportModel >::kSource().
|
delete |
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 156 of file dynamicKEqn.H.
References dynamicKEqn< BasicMomentumTransportModel >::read().
|
protected |
Return the KK field obtained by filtering the velocity field U.
Definition at line 41 of file dynamicKEqn.C.
References dynamicKEqn< BasicMomentumTransportModel >::Ck(), Foam::dimVelocity, Foam::magSqr(), Foam::max(), and Foam::sqr().
|
protected |
Calculate Ck by filtering the velocity field U.
Definition at line 53 of file dynamicKEqn.C.
References dynamicKEqn< BasicMomentumTransportModel >::Ce(), delta, Foam::dev(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::mag(), Foam::magSqr(), Foam::sqr(), and Foam::sqrt().
Referenced by dynamicKEqn< BasicMomentumTransportModel >::KK().
|
protected |
Calculate Ce by filtering the velocity field U.
Definition at line 84 of file dynamicKEqn.C.
References delta, Foam::mag(), Foam::magSqr(), and Foam::pow().
|
protected |
Definition at line 101 of file dynamicKEqn.C.
References dynamicKEqn< BasicMomentumTransportModel >::correctNut(), Foam::dev(), Foam::fvc::grad(), and Foam::symm().
Referenced by dynamicKEqn< BasicMomentumTransportModel >::Ck().
|
protected |
Update sub-grid eddy-viscosity.
Definition at line 110 of file dynamicKEqn.C.
References delta, dictionary::New(), and Foam::sqrt().
|
protectedvirtual |
Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.
Definition at line 122 of file dynamicKEqn.C.
References Foam::fvc::grad(), and Foam::symm().
Referenced by dynamicKEqn< BasicMomentumTransportModel >::Ce().
|
protectedvirtual |
Definition at line 129 of file dynamicKEqn.C.
References Foam::dimTime, Foam::dimVolume, and dynamicKEqn< BasicMomentumTransportModel >::dynamicKEqn().
TypeName | ( | "dynamicKEqn< BasicMomentumTransportModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Reimplemented from LESeddyViscosity< BasicMomentumTransportModel >.
Definition at line 197 of file dynamicKEqn.C.
Referenced by dynamicKEqn< BasicMomentumTransportModel >::~dynamicKEqn().
|
inlinevirtual |
Return SGS kinetic energy.
Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.
Definition at line 166 of file dynamicKEqn.H.
References dynamicKEqn< BasicMomentumTransportModel >::epsilon(), and dynamicKEqn< BasicMomentumTransportModel >::k_.
|
virtual |
Return sub-grid disipation rate.
Reimplemented from LESeddyViscosity< BasicMomentumTransportModel >.
Definition at line 213 of file dynamicKEqn.C.
References delta, IOobject::groupName(), k, GeometricField< scalar, fvPatchField, volMesh >::New(), and Foam::sqrt().
Referenced by dynamicKEqn< BasicMomentumTransportModel >::k().
|
inline |
Return the effective diffusivity for k.
Definition at line 175 of file dynamicKEqn.H.
References dynamicKEqn< BasicMomentumTransportModel >::correct(), GeometricField< scalar, fvPatchField, volMesh >::New(), eddyViscosity< LESModel< BasicMomentumTransportModel > >::nut_, and dynamicKEqn< BasicMomentumTransportModel >::operator=().
|
virtual |
Correct Eddy-Viscosity and related properties.
Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.
Definition at line 224 of file dynamicKEqn.C.
References Foam::fvc::absolute(), Foam::bound(), tmp< T >::clear(), fvConstraints::constrain(), eddyViscosity< LESModel< BasicMomentumTransportModel > >::correct(), Foam::fvm::ddt(), delta, Foam::dev(), Foam::fvm::div(), Foam::fvc::div(), divU, fvConstraints, fvModels, Foam::constant::universal::G, Foam::fvc::grad(), Foam::fvm::laplacian(), dictionary::New(), nut, phi, Foam::solve(), Foam::fvm::Sp(), Foam::sqrt(), Foam::fvm::SuSp(), and Foam::symm().
Referenced by dynamicKEqn< BasicMomentumTransportModel >::DkEff().
|
delete |
Disallow default bitwise assignment.
Referenced by dynamicKEqn< BasicMomentumTransportModel >::DkEff().
|
protected |
Definition at line 83 of file dynamicKEqn.H.
Referenced by dynamicKEqn< BasicMomentumTransportModel >::k().
|
protected |
Definition at line 88 of file dynamicKEqn.H.
Definition at line 89 of file dynamicKEqn.H.
|
protected |
Definition at line 90 of file dynamicKEqn.H.