General base class for all incompressible models that directly solve for the SGS stress tensor B. More...
#include <incompressibleLESModels/GenSGSStress.H>
General base class for all incompressible models that directly solve for the SGS stress tensor B.
Contains tensor fields B (the SGS stress tensor) as well as scalar fields for k (SGS turbulent energy) gamma (SGS viscosity) and epsilon (SGS dissipation).
Definition at line 58 of file GenSGSStress.H.
Public Member Functions | |
GenSGSStress (const volVectorField &U, const surfaceScalarField &phi, transportModel &transport) | |
Construct from components.
| |
virtual | ~GenSGSStress () |
Destructor.
| |
virtual tmp< volScalarField > | k () const |
Return the SGS turbulent kinetic energy.
| |
virtual tmp< volScalarField > | epsilon () const |
Return the SGS turbulent dissipation.
| |
virtual tmp< volScalarField > | nuSgs () const |
Return the SGS viscosity.
| |
virtual tmp< volSymmTensorField > | B () const |
Return the sub-grid stress tensor.
| |
virtual tmp< volSymmTensorField > | devBeff () const |
Return the effective sub-grid turbulence stress tensor.
| |
virtual tmp< fvVectorMatrix > | divDevBeff (volVectorField &U) const |
Returns div(B).
| |
virtual bool | read () |
Read LESProperties dictionary.
| |
Protected Attributes | |
dimensionedScalar | ce_ |
dimensionedScalar | couplingFactor_ |
volSymmTensorField | B_ |
volScalarField | nuSgs_ |
GenSGSStress | ( | const volVectorField & | U, |
const surfaceScalarField & | phi, | ||
transportModel & | transport | ||
) |
Construct from components.
Definition at line 40 of file GenSGSStress.C.
References Foam::exit(), Foam::FatalError, FatalErrorIn, and Foam::nl.
virtual ~GenSGSStress | ( | ) | [inline, virtual]
|
Destructor.
Definition at line 93 of file GenSGSStress.H.
virtual tmp<volScalarField> k | ( | ) | const [inline, virtual]
|
Return the SGS turbulent kinetic energy.
Implements incompressible::LESModel.
Definition at line 100 of file GenSGSStress.H.
References GenSGSStress::B_, and Foam::tr().
Referenced by GenSGSStress::epsilon().
virtual tmp<volScalarField> epsilon | ( | ) | const [inline, virtual]
|
Return the SGS turbulent dissipation.
Implements incompressible::LESModel.
Definition at line 106 of file GenSGSStress.H.
References GenSGSStress::ce_, incompressible::LESModel::delta(), GenSGSStress::k(), K, and Foam::sqrt().
virtual tmp<volScalarField> nuSgs | ( | ) | const [inline, virtual]
|
Return the SGS viscosity.
Implements incompressible::LESModel.
Definition at line 113 of file GenSGSStress.H.
References GenSGSStress::nuSgs_.
virtual tmp<volSymmTensorField> B | ( | ) | const [inline, virtual]
|
Return the sub-grid stress tensor.
Implements incompressible::LESModel.
Definition at line 119 of file GenSGSStress.H.
References GenSGSStress::B_.
tmp< volSymmTensorField > devBeff | ( | ) | const [virtual]
|
Return the effective sub-grid turbulence stress tensor.
including the laminar stress
Implements incompressible::LESModel.
Definition at line 111 of file GenSGSStress.C.
References GenSGSStress::B_, Foam::dev(), Foam::fvc::grad(), incompressible::turbulenceModel::mesh_, IOobject::NO_READ, IOobject::NO_WRITE, incompressible::turbulenceModel::nu(), incompressible::turbulenceModel::runTime_, Time::timeName(), Foam::twoSymm(), and incompressible::turbulenceModel::U().
tmp< fvVectorMatrix > divDevBeff | ( | volVectorField & | U ) | const [virtual]
|
Returns div(B).
This is the additional term due to the filtering of the NSE.
Implements incompressible::LESModel.
Definition at line 132 of file GenSGSStress.C.
References Foam::fvc::div(), Foam::fvc::grad(), and Foam::resError::laplacian().
bool read | ( | ) | [virtual]
|
Read LESProperties dictionary.
Implements incompressible::LESModel.
Reimplemented in DeardorffDiffStress, and LRRDiffStress.
Definition at line 160 of file GenSGSStress.C.
References GenSGSStress::ce_, incompressible::LESModel::coeffDict(), GenSGSStress::couplingFactor_, Foam::exit(), Foam::FatalError, FatalErrorIn, dimensioned< Type >::readIfPresent(), and dimensioned< Type >::value().
Referenced by LRRDiffStress::read(), and DeardorffDiffStress::read().
dimensionedScalar ce_ [protected]
|
Definition at line 71 of file GenSGSStress.H.
Referenced by DeardorffDiffStress::correct(), GenSGSStress::epsilon(), and GenSGSStress::read().
dimensionedScalar couplingFactor_ [protected]
|
Definition at line 73 of file GenSGSStress.H.
Referenced by GenSGSStress::read().
volSymmTensorField B_ [protected]
|
Definition at line 75 of file GenSGSStress.H.
Referenced by GenSGSStress::B(), LRRDiffStress::correct(), DeardorffDiffStress::correct(), GenSGSStress::devBeff(), and GenSGSStress::k().
volScalarField nuSgs_ [protected]
|
Definition at line 76 of file GenSGSStress.H.
Referenced by LRRDiffStress::DBEff(), DeardorffDiffStress::DBEff(), and GenSGSStress::nuSgs().