Differential SGS Stress Equation Model for incompressible and compressible flows. More...
#include <DeardorffDiffStress.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from ReynoldsStress< LESModel< BasicTurbulenceModel > > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("DeardorffDiffStress") | |
| Runtime type information. | |
| DeardorffDiffStress (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) | |
| Constructor from components. | |
| virtual | ~DeardorffDiffStress ()=default |
| Destructor. | |
| virtual bool | read () |
| Read model coefficients if they have changed. | |
| virtual void | correct () |
| Correct sub-grid stress, eddy-Viscosity and related properties. | |
| Public Member Functions inherited from ReynoldsStress< LESModel< BasicTurbulenceModel > > | |
| Foam::tmp< Foam::fvVectorMatrix > | DivDevRhoReff (const RhoFieldType &rho, volVectorField &U) const |
| ReynoldsStress (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 | ~ReynoldsStress ()=default |
| Destructor. | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity. | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy. | |
| virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual void | validate () |
| Validate the turbulence fields after construction. | |
Protected Member Functions | |
| virtual void | correctNut () |
| Update the eddy-viscosity. | |
| Protected Member Functions inherited from ReynoldsStress< LESModel< BasicTurbulenceModel > > | |
| void | boundNormalStress (volSymmTensorField &R) const |
| void | correctWallShearStress (volSymmTensorField &R) const |
| void | checkRealizabilityConditions (const volSymmTensorField &R) const |
Protected Attributes | |
| dimensionedScalar | Ck_ |
| dimensionedScalar | Cm_ |
| dimensionedScalar | Ce_ |
| dimensionedScalar | Cs_ |
| Protected Attributes inherited from ReynoldsStress< LESModel< BasicTurbulenceModel > > | |
| dimensionedScalar | couplingFactor_ |
| volSymmTensorField | R_ |
| volScalarField | nut_ |
Differential SGS Stress Equation Model for incompressible and compressible flows.
Reference:
Deardorff, J. W. (1973).
The use of subgrid transport equations in a three-dimensional model
of atmospheric turbulence.
Journal of Fluids Engineering, 95(3), 429-438.
This SGS model uses a full balance equation for the SGS stress tensor to simulate the behaviour of B.
This implementation is as described in the above paper except that the triple correlation model of Donaldson is replaced with the generalized gradient diffusion model of Daly and Harlow:
Daly, B. J., & Harlow, F. H. (1970).
Transport equations in turbulence.
Physics of Fluids (1958-1988), 13(11), 2634-2649.
with the default value for the coefficient Cs of 0.25 from
Launder, B. E., Reece, G. J., & Rodi, W. (1975).
Progress in the development of a Reynolds-stress turbulence closure.
Journal of fluid mechanics, 68(03), 537-566.
Definition at line 81 of file DeardorffDiffStress.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 120 of file DeardorffDiffStress.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 121 of file DeardorffDiffStress.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 122 of file DeardorffDiffStress.H.
| DeardorffDiffStress | ( | 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 ) |
Constructor from components.
Definition at line 49 of file DeardorffDiffStress.C.
References alpha, ReynoldsStress< LESModel< BasicTurbulenceModel > >::boundNormalStress(), Ce_, Ck_, Cm_, Cs_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, phi, ReynoldsStress< LESModel< BasicTurbulenceModel > >::R_, ReynoldsStress< LESModel< BasicTurbulenceModel > >::ReynoldsStress(), rho, and U.

|
virtualdefault |
Destructor.
|
protectedvirtual |
Update the eddy-viscosity.
Implements ReynoldsStress< LESModel< BasicTurbulenceModel > >.
Definition at line 36 of file DeardorffDiffStress.C.
References Ck_, optionList::correct(), delta, k, options::New(), ReynoldsStress< LESModel< BasicTurbulenceModel > >::nut_, and Foam::sqrt().
Referenced by correct().


| TypeName | ( | "DeardorffDiffStress< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Implements ReynoldsStress< LESModel< BasicTurbulenceModel > >.
Definition at line 121 of file DeardorffDiffStress.C.
References Ce_, Ck_, Cm_, Cs_, read(), and ReynoldsStress< LESModel< BasicTurbulenceModel > >::ReynoldsStress().
Referenced by read().


|
virtual |
Correct sub-grid stress, eddy-Viscosity and related properties.
Implements ReynoldsStress< LESModel< BasicTurbulenceModel > >.
Definition at line 138 of file DeardorffDiffStress.C.
References alpha, ReynoldsStress< LESModel< BasicTurbulenceModel > >::boundNormalStress(), Ce_, Cm_, correct(), correctNut(), Cs_, D, Foam::fvm::ddt(), delta, Foam::fvm::div(), epsilon, fvOptions, Foam::fvc::grad(), Foam::I, ReynoldsStress< LESModel< BasicTurbulenceModel > >::k(), k, Foam::fvm::laplacian(), options::New(), nu, ReynoldsStress< LESModel< BasicTurbulenceModel > >::R(), ReynoldsStress< LESModel< BasicTurbulenceModel > >::R_, tmp< T >::ref(), ReynoldsStress< LESModel< BasicTurbulenceModel > >::ReynoldsStress(), rho, Foam::fvm::Sp(), Foam::sqrt(), Foam::symm(), Foam::twoSymm(), and U.

|
protected |
Definition at line 104 of file DeardorffDiffStress.H.
Referenced by correctNut(), DeardorffDiffStress(), and read().
|
protected |
Definition at line 105 of file DeardorffDiffStress.H.
Referenced by correct(), DeardorffDiffStress(), and read().
|
protected |
Definition at line 106 of file DeardorffDiffStress.H.
Referenced by correct(), DeardorffDiffStress(), and read().
|
protected |
Definition at line 107 of file DeardorffDiffStress.H.
Referenced by correct(), DeardorffDiffStress(), and read().