Differential SGS Stress Equation Model for incompressible flows. More...
#include <incompressibleLESModels/DeardorffDiffStress.H>
Differential SGS Stress Equation Model for incompressible flows.
The DSEM uses a model version of the full balance equation for the SGS stress tensor to simulate the behaviour of B. Thus,
d/dt(B) + div(U*B) - div(nuSgs*grad(B)) = P - c1*epsilon/k*B - 0.667*(1 - c1)*epsilon*I - c2*(P - 0.333*trP*I) where k = 0.5*tr(B), epsilon = ce*k^3/2/delta, epsilon/k = ce*k^1/2/delta P = -(B'L + L'B) nuSgs = ck*sqrt(k)*delta nuEff = nuSgs + nu
Definition at line 71 of file DeardorffDiffStress.H.
Public Member Functions | |
TypeName ("DeardorffDiffStress") | |
Runtime type information.
| |
DeardorffDiffStress (const volVectorField &U, const surfaceScalarField &phi, transportModel &transport) | |
Construct from components.
| |
virtual | ~DeardorffDiffStress () |
Destructor.
| |
tmp< volScalarField > | DBEff () const |
Return the effective diffusivity for B.
| |
virtual void | correct (const tmp< volTensorField > &gradU) |
Correct Eddy-Viscosity and related properties.
| |
virtual bool | read () |
Read LESProperties dictionary.
|
DeardorffDiffStress | ( | const volVectorField & | U, |
const surfaceScalarField & | phi, | ||
transportModel & | transport | ||
) |
Construct from components.
Definition at line 55 of file DeardorffDiffStress.C.
References Foam::tr().
virtual ~DeardorffDiffStress | ( | ) | [inline, virtual]
|
Destructor.
Definition at line 108 of file DeardorffDiffStress.H.
TypeName | ( | "DeardorffDiffStress" | ) |
Runtime type information.
tmp<volScalarField> DBEff | ( | ) | const [inline]
|
Return the effective diffusivity for B.
Definition at line 115 of file DeardorffDiffStress.H.
References incompressible::turbulenceModel::nu(), and GenSGSStress::nuSgs_.
Referenced by DeardorffDiffStress::correct().
void correct | ( | const tmp< volTensorField > & | gradU ) | [virtual]
|
Correct Eddy-Viscosity and related properties.
Reimplemented from incompressible::LESModel.
Definition at line 91 of file DeardorffDiffStress.C.
References GenSGSStress::B_, Foam::bound(), GenSGSStress::ce_, Foam::component(), GeometricField< Type, PatchField, GeoMesh >::component(), incompressible::LESModel::correct(), DeardorffDiffStress::DBEff(), Foam::fvm::ddt(), incompressible::LESModel::delta(), Foam::fvm::div(), forAll, Foam::I, incompressible::LESModel::k0(), Foam::fvm::laplacian(), Foam::magSqr(), Foam::max(), incompressible::LESModel::nuEff(), incompressible::turbulenceModel::phi(), fvMatrix< Type >::relax(), fvMatrix< Type >::solve(), Foam::fvm::Sp(), Foam::sqrt(), Foam::symm(), Foam::tr(), Foam::twoSymm(), SymmTensor< Cmpt >::XX, SymmTensor< Cmpt >::YY, and SymmTensor< Cmpt >::ZZ.
bool read | ( | ) | [virtual]
|
Read LESProperties dictionary.
Reimplemented from GenSGSStress.
Definition at line 138 of file DeardorffDiffStress.C.
References incompressible::LESModel::coeffDict(), GenSGSStress::read(), and dimensioned< Type >::readIfPresent().