The Smagorinsky SGS model. More...


Public Types | |
| typedef BasicMomentumTransportModel::alphaField | alphaField |
| typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Types inherited from LESeddyViscosity< BasicMomentumTransportModel > | |
| typedef BasicMomentumTransportModel::alphaField | alphaField |
| typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Types inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > > | |
| typedef BasicMomentumTransportModel::alphaField | alphaField |
| typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Types inherited from linearViscousStress< BasicMomentumTransportModel > | |
| typedef BasicMomentumTransportModel::alphaField | alphaField |
| typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Member Functions | |
| TypeName ("Smagorinsky") | |
| Runtime type information. More... | |
| Smagorinsky (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName) | |
| Construct from components. More... | |
| Smagorinsky (const Smagorinsky &)=delete | |
| Disallow default bitwise copy construction. More... | |
| virtual | ~Smagorinsky () |
| Destructor. More... | |
| virtual bool | read () |
| Read model coefficients if they have changed. More... | |
| virtual tmp< volScalarField > | k () const |
| Return SGS kinetic energy. More... | |
| virtual void | correct () |
| Correct Eddy-Viscosity and related properties. More... | |
| void | operator= (const Smagorinsky &)=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 viscosity &viscosity) | |
| Construct from components. More... | |
| LESeddyViscosity (const LESeddyViscosity &)=delete | |
| Disallow default bitwise copy construction. More... | |
| virtual | ~LESeddyViscosity () |
| Destructor. More... | |
| virtual tmp< volScalarField > | epsilon () const |
| Return sub-grid disipation rate. More... | |
| virtual tmp< volScalarField > | omega () const |
| Return the turbulence specific dissipation rate. 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 viscosity &viscosity) | |
| 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 > | R () 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< BasicMomentumTransportModel > | |
| linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
| Construct from components. More... | |
| virtual | ~linearViscousStress () |
| Destructor. More... | |
| virtual tmp< surfaceVectorField > | devTau () const |
| Return the effective surface stress. 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... | |
| template<class RhoFieldType > | |
| Foam::tmp< Foam::fvVectorMatrix > | DivDevTau (const RhoFieldType &rho, volVectorField &U) const |
Protected Member Functions | |
| tmp< volScalarField > | k (const tmp< volTensorField > &gradU) const |
| Return SGS kinetic energy. More... | |
| virtual void | correctNut () |
| Update the SGS eddy viscosity. More... | |
Protected Member Functions inherited from linearViscousStress< BasicMomentumTransportModel > | |
| template<class RhoFieldType > | |
| tmp< fvVectorMatrix > | DivDevTau (const RhoFieldType &rho, volVectorField &U) const |
| Return the source term for the momentum equation. More... | |
Additional Inherited Members | |
Protected Attributes inherited from LESeddyViscosity< BasicMomentumTransportModel > | |
| dimensionedScalar | Ck_ |
| dimensionedScalar | Ce_ |
Protected Attributes inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > > | |
| volScalarField | nut_ |
The Smagorinsky SGS model.
Reference:
Smagorinsky, J. (1963).
General circulation experiments with the primitive equations: I.
The basic experiment*.
Monthly weather review, 91(3), 99-164.
The form of the Smagorinsky model implemented is obtained from the k-equation model assuming local equilibrium which provides estimates of both k and epsilon separate from the sub-grid scale viscosity:
B = 2/3*k*I - 2*nuSgs*dev(D)
where
D = symm(grad(U));
k from D:B + Ce*k^3/2/delta = 0
nuSgs = Ck*sqrt(k)*delta
The default model coefficients are
SmagorinskyCoeffs
{
Ck 0.094;
Ce 1.048;
}
Definition at line 86 of file Smagorinsky.H.
| typedef BasicMomentumTransportModel::alphaField alphaField |
Definition at line 104 of file Smagorinsky.H.
| typedef BasicMomentumTransportModel::rhoField rhoField |
Definition at line 105 of file Smagorinsky.H.
| Smagorinsky | ( | const alphaField & | alpha, |
| const rhoField & | rho, | ||
| const volVectorField & | U, | ||
| const surfaceScalarField & | alphaRhoPhi, | ||
| const surfaceScalarField & | phi, | ||
| const viscosity & | viscosity, | ||
| const word & | type = typeName |
||
| ) |
Construct from components.
Definition at line 73 of file Smagorinsky.C.
|
delete |
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 131 of file Smagorinsky.H.
|
protected |
Return SGS kinetic energy.
calculated from the given velocity gradient
Definition at line 40 of file Smagorinsky.C.
References b, Foam::constant::universal::c, Foam::saturationModels::D, delta, Foam::dev(), GeometricField< Type, GeoMesh, PrimitiveField >::New(), Foam::sqr(), Foam::sqrt(), Foam::symm(), and Foam::tr().

|
protectedvirtual |
Update the SGS eddy viscosity.
Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.
Reimplemented in SmagorinskyZhang< BasicMomentumTransportModel >.
Definition at line 60 of file Smagorinsky.C.
References delta, Foam::fvc::grad(), k, dictionary::New(), and Foam::sqrt().

| TypeName | ( | "Smagorinsky< BasicMomentumTransportModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Reimplemented from LESeddyViscosity< BasicMomentumTransportModel >.
Reimplemented in SmagorinskyZhang< BasicMomentumTransportModel >.
Definition at line 100 of file Smagorinsky.C.
References LESeddyViscosity< BasicMomentumTransportModel >::read().

|
inlinevirtual |
Return SGS kinetic energy.
Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.
Definition at line 141 of file Smagorinsky.H.
References Foam::fvc::grad().

|
virtual |
Correct Eddy-Viscosity and related properties.
Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.
Definition at line 107 of file Smagorinsky.C.
References eddyViscosity< LESModel< BasicMomentumTransportModel > >::correct().

|
delete |
Disallow default bitwise assignment.