The Smagorinsky SGS model including bubble-generated turbulence. More...


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Types inherited from Smagorinsky< BasicTurbulenceModel > | |
| 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 LESModel< BasicTurbulenceModel >::alphaField | alphaField |
| typedef LESModel< BasicTurbulenceModel >::rhoField | rhoField |
| typedef LESModel< BasicTurbulenceModel >::transportModel | transportModel |
Public Types inherited from linearViscousStress< LESModel< BasicTurbulenceModel > > | |
| typedef LESModel< BasicTurbulenceModel >::alphaField | alphaField |
| typedef LESModel< BasicTurbulenceModel >::rhoField | rhoField |
| typedef LESModel< BasicTurbulenceModel >::transportModel | transportModel |
Public Types inherited from LESModel< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("SmagorinskyZhang") | |
| Runtime type information. More... | |
| SmagorinskyZhang (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 | ~SmagorinskyZhang () |
| Destructor. More... | |
| virtual bool | read () |
| Read model coefficients if they have changed. More... | |
Public Member Functions inherited from Smagorinsky< BasicTurbulenceModel > | |
| TypeName ("Smagorinsky") | |
| Runtime type information. More... | |
| Smagorinsky (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 | ~Smagorinsky () |
| Destructor. More... | |
| virtual tmp< volScalarField > | k () const |
| Return SGS kinetic energy. More... | |
| virtual tmp< volScalarField > | epsilon () const |
| Return sub-grid disipation rate. More... | |
| virtual void | correct () |
| Correct Eddy-Viscosity and related properties. More... | |
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. More... | |
| virtual | ~LESeddyViscosity () |
| Destructor. More... | |
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. More... | |
| virtual | ~eddyViscosity () |
| Destructor. More... | |
| volScalarField & | evNut () |
| Return non-const access to the turbulence viscosity. 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< LESModel< 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 () |
| Destructor. More... | |
| virtual tmp< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. 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 LESModel< BasicTurbulenceModel > | |
| 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, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)) | |
| LESModel (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 | ~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... | |
Protected Member Functions | |
| virtual void | correctNut () |
| Update the SGS eddy viscosity. More... | |
Protected Member Functions inherited from Smagorinsky< BasicTurbulenceModel > | |
| tmp< volScalarField > | k (const tmp< volTensorField > &gradU) const |
| Return SGS kinetic energy. More... | |
Protected Member Functions inherited from LESModel< BasicTurbulenceModel > | |
| virtual void | printCoeffs (const word &type) |
| Print model coefficients. More... | |
Protected Attributes | |
| dimensionedScalar | Cmub_ |
Protected Attributes inherited from Smagorinsky< BasicTurbulenceModel > | |
| dimensionedScalar | Ck_ |
Protected Attributes inherited from LESeddyViscosity< BasicTurbulenceModel > | |
| dimensionedScalar | Ce_ |
Protected Attributes inherited from eddyViscosity< LESModel< BasicTurbulenceModel > > | |
| volScalarField | nut_ |
Protected Attributes inherited from LESModel< BasicTurbulenceModel > | |
| 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< BasicTurbulenceModel > | |
| static autoPtr< LESModel > | 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 LES model. More... | |
The Smagorinsky SGS model including bubble-generated turbulence.
Zhang, D., Deen, N. G., & Kuipers, J. A. M. (2006).
Numerical simulation of the dynamic flow behavior in a bubble column:
a study of closures for turbulence and interface forces.
Chemical Engineering Science, 61(23), 7593-7608.
The default model coefficients are
SmagorinskyZhangCoeffs
{
Ck 0.094;
Ce 1.048;
Cmub 0.6;
}
Definition at line 73 of file SmagorinskyZhang.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 115 of file SmagorinskyZhang.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 116 of file SmagorinskyZhang.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 117 of file SmagorinskyZhang.H.
| SmagorinskyZhang | ( | 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 40 of file SmagorinskyZhang.C.
|
inlinevirtual |
Destructor.
Definition at line 141 of file SmagorinskyZhang.H.
References SmagorinskyZhang< BasicTurbulenceModel >::read().

|
protectedvirtual |
Update the SGS eddy viscosity.
Reimplemented from Smagorinsky< BasicTurbulenceModel >.
Definition at line 133 of file SmagorinskyZhang.C.
References TurbulenceModel< volScalarField, volScalarField, compressibleTurbulenceModel, TransportModel >::alpha(), optionList::correct(), delta, Foam::fvc::grad(), k, Foam::mag(), options::New(), Foam::sqrt(), TurbulenceModel< volScalarField, volScalarField, compressibleTurbulenceModel, TransportModel >::transport(), and turbulenceModel::U().

| TypeName | ( | "SmagorinskyZhang< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Reimplemented from Smagorinsky< BasicTurbulenceModel >.
Definition at line 85 of file SmagorinskyZhang.C.
References IOobject::db(), fluid, IOobject::groupName(), objectRegistry::lookupObject(), twoPhaseSystem::otherPhase(), and turbulenceModel::propertiesName.
Referenced by SmagorinskyZhang< BasicTurbulenceModel >::~SmagorinskyZhang().


|
protected |
Definition at line 105 of file SmagorinskyZhang.H.
1.8.11