#include <src/turbulenceModels/incompressible/LES/SpalartAllmaras/SpalartAllmaras.H>
Definition at line 55 of file SpalartAllmaras.H.
Public Member Functions | |
TypeName ("SpalartAllmaras") | |
Runtime type information.
| |
SpalartAllmaras (const volVectorField &U, const surfaceScalarField &phi, transportModel &transport, const word &modelName=typeName) | |
Construct from components.
| |
virtual | ~SpalartAllmaras () |
Destructor.
| |
virtual tmp< volScalarField > | k () const |
Return SGS kinetic energy.
| |
virtual tmp< volScalarField > | epsilon () const |
Return sub-grid disipation rate.
| |
tmp< volScalarField > | nuTilda () const |
virtual tmp< volScalarField > | nuSgs () const |
Return 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 |
Return the deviatoric part of the divergence of Beff.
| |
virtual void | correct (const tmp< volTensorField > &gradU) |
Correct nuTilda and related properties.
| |
virtual bool | read () |
Read LESProperties dictionary.
| |
Protected Member Functions | |
virtual tmp< volScalarField > | fv1 () const |
virtual tmp< volScalarField > | fv2 () const |
virtual tmp< volScalarField > | fv3 () const |
virtual tmp< volScalarField > | S (const volTensorField &gradU) const |
virtual tmp< volScalarField > | STilda (const volScalarField &S, const volScalarField &dTilda) const |
virtual tmp< volScalarField > | r (const volScalarField &visc, const volScalarField &S, const volScalarField &dTilda) const |
virtual tmp< volScalarField > | fw (const volScalarField &S, const volScalarField &dTilda) const |
virtual tmp< volScalarField > | dTilda (const volScalarField &S) const |
Length scale.
| |
Protected Attributes | |
dimensionedScalar | sigmaNut_ |
dimensionedScalar | kappa_ |
dimensionedScalar | Cb1_ |
dimensionedScalar | Cb2_ |
dimensionedScalar | Cv1_ |
dimensionedScalar | Cv2_ |
dimensionedScalar | CDES_ |
dimensionedScalar | ck_ |
dimensionedScalar | Cw1_ |
dimensionedScalar | Cw2_ |
dimensionedScalar | Cw3_ |
wallDist | y_ |
volScalarField | nuTilda_ |
volScalarField | nuSgs_ |
SpalartAllmaras | ( | const volVectorField & | U, |
const surfaceScalarField & | phi, | ||
transportModel & | transport, | ||
const word & | modelName = typeName
|
||
) |
Construct from components.
Definition at line 148 of file SpalartAllmaras.C.
virtual ~SpalartAllmaras | ( | ) | [inline, virtual]
|
Destructor.
Definition at line 146 of file SpalartAllmaras.H.
tmp< volScalarField > fv1 | ( | ) | const [protected, virtual]
|
Definition at line 54 of file SpalartAllmaras.C.
References SpalartAllmaras::Cv1_, incompressible::turbulenceModel::nu(), SpalartAllmaras::nuTilda_, and Foam::pow3().
Referenced by SpalartAllmarasIDDES::dTilda(), and SpalartAllmaras::fv3().
tmp< volScalarField > fv2 | ( | ) | const [protected, virtual]
|
Definition at line 61 of file SpalartAllmaras.C.
References SpalartAllmaras::Cv2_, incompressible::turbulenceModel::nu(), SpalartAllmaras::nuTilda_, and Foam::pow3().
Referenced by SpalartAllmarasIDDES::dTilda().
tmp< volScalarField > fv3 | ( | ) | const [protected, virtual]
|
Definition at line 67 of file SpalartAllmaras.C.
References SpalartAllmaras::Cv2_, SpalartAllmaras::fv1(), incompressible::turbulenceModel::nu(), SpalartAllmaras::nuTilda_, Foam::pow3(), and Foam::sqr().
tmp< volScalarField > S | ( | const volTensorField & | gradU ) | const [protected, virtual]
|
Reimplemented in SpalartAllmarasDDES.
Definition at line 80 of file SpalartAllmaras.C.
References Foam::mag(), Foam::skew(), and Foam::sqrt().
Referenced by SpalartAllmaras::correct(), and SpalartAllmaras::k().
tmp< volScalarField > STilda | ( | const volScalarField & | S, |
const volScalarField & | dTilda | ||
) | const [protected, virtual]
|
Definition at line 87 of file SpalartAllmaras.C.
References Foam::sqr().
Referenced by SpalartAllmaras::correct().
tmp< volScalarField > r | ( | const volScalarField & | visc, |
const volScalarField & | S, | ||
const volScalarField & | dTilda | ||
) | const [protected, virtual]
|
Definition at line 97 of file SpalartAllmaras.C.
References DimensionedField< Type, GeoMesh >::dimensions(), Foam::max(), Foam::min(), and Foam::sqr().
tmp< volScalarField > fw | ( | const volScalarField & | S, |
const volScalarField & | dTilda | ||
) | const [protected, virtual]
|
Definition at line 126 of file SpalartAllmaras.C.
References Foam::pow(), and Foam::pow6().
Referenced by SpalartAllmaras::correct().
tmp< volScalarField > dTilda | ( | const volScalarField & | S ) | const [protected, virtual]
|
Length scale.
Reimplemented in SpalartAllmarasDDES, and SpalartAllmarasIDDES.
Definition at line 139 of file SpalartAllmaras.C.
References SpalartAllmaras::CDES_, incompressible::LESModel::delta(), Foam::min(), and SpalartAllmaras::y_.
Referenced by SpalartAllmaras::correct(), and SpalartAllmaras::k().
TypeName | ( | "SpalartAllmaras" | ) |
Runtime type information.
tmp< volScalarField > k | ( | ) | const [virtual]
|
Return SGS kinetic energy.
Implements incompressible::LESModel.
Definition at line 323 of file SpalartAllmaras.C.
References SpalartAllmaras::ck_, SpalartAllmaras::dTilda(), Foam::fvc::grad(), SpalartAllmaras::nuSgs(), SpalartAllmaras::S(), Foam::sqr(), and incompressible::turbulenceModel::U().
Referenced by SpalartAllmaras::B().
tmp< volScalarField > epsilon | ( | ) | const [virtual]
|
Return sub-grid disipation rate.
Implements incompressible::LESModel.
Definition at line 329 of file SpalartAllmaras.C.
References Foam::fvc::grad(), Foam::magSqr(), incompressible::LESModel::nuEff(), Foam::symm(), and incompressible::turbulenceModel::U().
tmp<volScalarField> nuTilda | ( | ) | const [inline]
|
Definition at line 158 of file SpalartAllmaras.H.
References SpalartAllmaras::nuTilda_.
virtual tmp<volScalarField> nuSgs | ( | ) | const [inline, virtual]
|
Return SGS viscosity.
Implements incompressible::LESModel.
Definition at line 164 of file SpalartAllmaras.H.
References SpalartAllmaras::nuSgs_.
Referenced by SpalartAllmaras::B(), and SpalartAllmaras::k().
tmp< volSymmTensorField > B | ( | ) | const [virtual]
|
Return the sub-grid stress tensor.
Implements incompressible::LESModel.
Definition at line 335 of file SpalartAllmaras.C.
References Foam::fvc::grad(), Foam::I, SpalartAllmaras::k(), SpalartAllmaras::nuSgs(), Foam::twoSymm(), and incompressible::turbulenceModel::U().
tmp< volSymmTensorField > devBeff | ( | ) | const [virtual]
|
Return the effective sub-grid turbulence stress tensor.
including the laminar stress
Implements incompressible::LESModel.
Definition at line 341 of file SpalartAllmaras.C.
References Foam::dev(), Foam::fvc::grad(), incompressible::LESModel::nuEff(), Foam::twoSymm(), and incompressible::turbulenceModel::U().
tmp< fvVectorMatrix > divDevBeff | ( | volVectorField & | U ) | const [virtual]
|
Return the deviatoric part of the divergence of Beff.
i.e. the additional term in the filtered NSE.
Implements incompressible::LESModel.
Definition at line 347 of file SpalartAllmaras.C.
References Foam::dev(), Foam::fvc::div(), Foam::fvc::grad(), Foam::fvm::laplacian(), incompressible::LESModel::nuEff(), and Foam::T().
void correct | ( | const tmp< volTensorField > & | gradU ) | [virtual]
|
Correct nuTilda and related properties.
Reimplemented from incompressible::LESModel.
Definition at line 283 of file SpalartAllmaras.C.
References Foam::bound(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), SpalartAllmaras::Cb1_, SpalartAllmaras::Cb2_, polyMesh::changing(), wallDist::correct(), incompressible::LESModel::correct(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), SpalartAllmaras::Cw1_, Foam::fvm::ddt(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::fvm::div(), SpalartAllmaras::dTilda(), SpalartAllmaras::fw(), Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), Foam::max(), incompressible::turbulenceModel::mesh_, incompressible::turbulenceModel::nu(), SpalartAllmaras::nuTilda_, incompressible::turbulenceModel::phi(), fvMatrix< Type >::relax(), SpalartAllmaras::S(), SpalartAllmaras::sigmaNut_, fvMatrix< Type >::solve(), Foam::fvm::Sp(), Foam::sqr(), SpalartAllmaras::STilda(), and SpalartAllmaras::y_.
bool read | ( | ) | [virtual]
|
Read LESProperties dictionary.
Implements incompressible::LESModel.
Reimplemented in SpalartAllmarasIDDES.
Definition at line 356 of file SpalartAllmaras.C.
References SpalartAllmaras::Cb1_, SpalartAllmaras::Cb2_, SpalartAllmaras::CDES_, SpalartAllmaras::ck_, incompressible::LESModel::coeffDict(), SpalartAllmaras::Cv1_, SpalartAllmaras::Cv2_, SpalartAllmaras::Cw1_, SpalartAllmaras::Cw2_, SpalartAllmaras::Cw3_, SpalartAllmaras::kappa_, dimensioned< Type >::readIfPresent(), SpalartAllmaras::sigmaNut_, and Foam::sqr().
Referenced by SpalartAllmarasIDDES::read().
dimensionedScalar sigmaNut_ [protected]
|
Definition at line 73 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::correct(), and SpalartAllmaras::read().
dimensionedScalar kappa_ [protected]
|
Definition at line 74 of file SpalartAllmaras.H.
Referenced by SpalartAllmarasIDDES::dTilda(), and SpalartAllmaras::read().
dimensionedScalar Cb1_ [protected]
|
Definition at line 79 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::correct(), SpalartAllmarasIDDES::dTilda(), and SpalartAllmaras::read().
dimensionedScalar Cb2_ [protected]
|
Definition at line 80 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::correct(), and SpalartAllmaras::read().
dimensionedScalar Cv1_ [protected]
|
Definition at line 81 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::fv1(), and SpalartAllmaras::read().
dimensionedScalar Cv2_ [protected]
|
Definition at line 82 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::fv2(), SpalartAllmaras::fv3(), and SpalartAllmaras::read().
dimensionedScalar CDES_ [protected]
|
Definition at line 83 of file SpalartAllmaras.H.
Referenced by SpalartAllmarasIDDES::dTilda(), SpalartAllmarasDDES::dTilda(), SpalartAllmaras::dTilda(), and SpalartAllmaras::read().
dimensionedScalar ck_ [protected]
|
Definition at line 84 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::k(), and SpalartAllmaras::read().
dimensionedScalar Cw1_ [protected]
|
Definition at line 85 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::correct(), SpalartAllmarasIDDES::dTilda(), and SpalartAllmaras::read().
dimensionedScalar Cw2_ [protected]
|
Definition at line 86 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::read().
dimensionedScalar Cw3_ [protected]
|
Definition at line 87 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::read().
Definition at line 92 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::correct(), SpalartAllmarasIDDES::dTilda(), SpalartAllmarasDDES::dTilda(), and SpalartAllmaras::dTilda().
volScalarField nuTilda_ [protected]
|
Definition at line 93 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::correct(), SpalartAllmaras::fv1(), SpalartAllmaras::fv2(), SpalartAllmaras::fv3(), and SpalartAllmaras::nuTilda().
volScalarField nuSgs_ [protected]
|
Definition at line 94 of file SpalartAllmaras.H.
Referenced by SpalartAllmaras::nuSgs().