SpalartAllmarasDES DES turbulence model for incompressible and compressible flows. More...


Public Types | |
| typedef BasicMomentumTransportModel::alphaField | alphaField |
| typedef BasicMomentumTransportModel::rhoField | rhoField |
| typedef BasicMomentumTransportModel::transportModel | transportModel |
Public Types inherited from LESeddyViscosity< BasicMomentumTransportModel > | |
| typedef BasicMomentumTransportModel::alphaField | alphaField |
| typedef BasicMomentumTransportModel::rhoField | rhoField |
| typedef BasicMomentumTransportModel::transportModel | transportModel |
Public Types inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > > | |
| typedef LESModel< BasicMomentumTransportModel > ::alphaField | alphaField |
| typedef LESModel< BasicMomentumTransportModel > ::rhoField | rhoField |
| typedef LESModel< BasicMomentumTransportModel > ::transportModel | transportModel |
Public Types inherited from linearViscousStress< LESModel< BasicMomentumTransportModel > > | |
| typedef LESModel< BasicMomentumTransportModel > ::alphaField | alphaField |
| typedef LESModel< BasicMomentumTransportModel > ::rhoField | rhoField |
| typedef LESModel< BasicMomentumTransportModel > ::transportModel | transportModel |
Public Types inherited from LESModel< BasicMomentumTransportModel > | |
| typedef BasicMomentumTransportModel::alphaField | alphaField |
| typedef BasicMomentumTransportModel::rhoField | rhoField |
| typedef BasicMomentumTransportModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("SpalartAllmarasDES") | |
| Runtime type information. More... | |
| SpalartAllmarasDES (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName) | |
| Construct from components. More... | |
| SpalartAllmarasDES (const SpalartAllmarasDES &)=delete | |
| Disallow default bitwise copy construction. More... | |
| virtual | ~SpalartAllmarasDES () |
| Destructor. More... | |
| virtual bool | read () |
| Read model coefficients if they have changed. More... | |
| tmp< volScalarField > | DnuTildaEff () const |
| Return the effective diffusivity for nuTilda. More... | |
| virtual tmp< volScalarField > | k () const |
| Return SGS kinetic energy. More... | |
| tmp< volScalarField > | nuTilda () const |
| virtual void | correct () |
| Correct nuTilda and related properties. More... | |
| void | operator= (const SpalartAllmarasDES &)=delete |
| Disallow default bitwise assignment. More... | |
Public Member Functions inherited from LESeddyViscosity< BasicMomentumTransportModel > | |
| LESeddyViscosity (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) | |
| Construct from components. More... | |
| LESeddyViscosity (const LESeddyViscosity &)=delete | |
| Disallow default bitwise copy construction. More... | |
| virtual | ~LESeddyViscosity () |
| Destructor. More... | |
| virtual tmp< volScalarField > | epsilon () const |
| Return sub-grid disipation rate. More... | |
| void | operator= (const LESeddyViscosity &)=delete |
| Disallow default bitwise assignment. More... | |
Public Member Functions inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > > | |
| eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) | |
| Construct from components. More... | |
| virtual | ~eddyViscosity () |
| Destructor. 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 > | sigma () const |
| Return the Reynolds stress tensor [m^2/s^2]. More... | |
| virtual void | validate () |
| Validate the turbulence fields after construction. More... | |
Public Member Functions inherited from linearViscousStress< LESModel< BasicMomentumTransportModel > > | |
| linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) | |
| Construct from components. More... | |
| virtual | ~linearViscousStress () |
| Destructor. More... | |
| virtual tmp< volSymmTensorField > | devTau () const |
| Return the effective stress tensor. More... | |
| virtual tmp< fvVectorMatrix > | divDevTau (volVectorField &U) const |
| Return the source term for the momentum equation. More... | |
| virtual tmp< fvVectorMatrix > | divDevTau (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. More... | |
Public Member Functions inherited from LESModel< BasicMomentumTransportModel > | |
| 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),(alpha, rho, U, alphaRhoPhi, phi, transport)) | |
| LESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) | |
| Construct from components. More... | |
| LESModel (const LESModel &)=delete | |
| Disallow default bitwise copy construction. 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... | |
| void | operator= (const LESModel &)=delete |
| Disallow default bitwise assignment. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from LESModel< BasicMomentumTransportModel > | |
| static autoPtr< LESModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) |
| Return a reference to the selected LES model. More... | |
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Spalart, P. R., Jou, W. H., Strelets, M., & Allmaras, S. R. (1997).
Comments on the feasibility of LES for wings, and on a hybrid
RANS/LES approach.
Advances in DNS/LES, 1, 4-8.
Definition at line 60 of file SpalartAllmarasDES.H.
| typedef BasicMomentumTransportModel::alphaField alphaField |
Definition at line 153 of file SpalartAllmarasDES.H.
| typedef BasicMomentumTransportModel::rhoField rhoField |
Definition at line 154 of file SpalartAllmarasDES.H.
| typedef BasicMomentumTransportModel::transportModel transportModel |
Definition at line 155 of file SpalartAllmarasDES.H.
| SpalartAllmarasDES | ( | const alphaField & | alpha, |
| const rhoField & | rho, | ||
| const volVectorField & | U, | ||
| const surfaceScalarField & | alphaRhoPhi, | ||
| const surfaceScalarField & | phi, | ||
| const transportModel & | transport, | ||
| const word & | type = typeName |
||
| ) |
Construct from components.
Definition at line 223 of file SpalartAllmarasDES.C.
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::correctNut().

|
delete |
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 181 of file SpalartAllmarasDES.H.
References SpalartAllmarasDES< BasicMomentumTransportModel >::DnuTildaEff(), SpalartAllmarasDES< BasicMomentumTransportModel >::k(), and SpalartAllmarasDES< BasicMomentumTransportModel >::read().

|
protected |
Definition at line 40 of file SpalartAllmarasDES.C.
References SpalartAllmarasDES< BasicMomentumTransportModel >::fv1(), and GeometricField< scalar, fvPatchField, volMesh >::New().

|
protected |
Definition at line 52 of file SpalartAllmarasDES.C.
References SpalartAllmarasDES< BasicMomentumTransportModel >::fv2(), GeometricField< scalar, fvPatchField, volMesh >::New(), and Foam::pow3().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::chi().


|
protected |
Definition at line 68 of file SpalartAllmarasDES.C.
References DimensionedField< Type, GeoMesh >::New(), and SpalartAllmarasDES< BasicMomentumTransportModel >::Omega().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::fv1().


|
protected |
Definition at line 84 of file SpalartAllmarasDES.C.
References Foam::mag(), DimensionedField< Type, GeoMesh >::New(), Foam::skew(), Foam::sqrt(), and SpalartAllmarasDES< BasicMomentumTransportModel >::Stilda().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::fv2().


|
protected |
Definition at line 99 of file SpalartAllmarasDES.C.
References Foam::max(), DimensionedField< Type, GeoMesh >::New(), SpalartAllmarasDES< BasicMomentumTransportModel >::r(), and Foam::sqr().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::Omega().


|
protected |
Definition at line 121 of file SpalartAllmarasDES.C.
References SpalartAllmarasDES< BasicMomentumTransportModel >::fw(), Foam::max(), Foam::min(), DimensionedField< Type, GeoMesh >::New(), and Foam::sqr().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::Stilda().


|
protected |
Definition at line 150 of file SpalartAllmarasDES.C.
References SpalartAllmarasDES< BasicMomentumTransportModel >::dTilda(), g, DimensionedField< Type, GeoMesh >::New(), Foam::pow(), and Foam::pow6().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::r().


|
protectedvirtual |
Length scale.
Reimplemented in SpalartAllmarasIDDES< BasicMomentumTransportModel >, and SpalartAllmarasDDES< BasicMomentumTransportModel >.
Definition at line 169 of file SpalartAllmarasDES.C.
References SpalartAllmarasDES< BasicMomentumTransportModel >::cacheLESRegion(), delta, Foam::min(), and DimensionedField< Type, GeoMesh >::New().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::fw().


|
protectedvirtual |
Cache the LES region indicator field.
Definition at line 185 of file SpalartAllmarasDES.C.
References SpalartAllmarasDES< BasicMomentumTransportModel >::correctNut(), Foam::neg(), and DimensionedField< Type, GeoMesh >::New().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::dTilda().


|
protected |
Definition at line 202 of file SpalartAllmarasDES.C.
References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), and dictionary::New().

|
protectedvirtual |
Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.
Definition at line 213 of file SpalartAllmarasDES.C.
References SpalartAllmarasDES< BasicMomentumTransportModel >::SpalartAllmarasDES().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::cacheLESRegion().


| TypeName | ( | "SpalartAllmarasDES< BasicMomentumTransportModel >" | ) |
Runtime type information.
|
virtual |
Read model coefficients if they have changed.
Reimplemented from LESeddyViscosity< BasicMomentumTransportModel >.
Reimplemented in SpalartAllmarasIDDES< BasicMomentumTransportModel >.
Definition at line 361 of file SpalartAllmarasDES.C.
References SpalartAllmarasDES< BasicMomentumTransportModel >::DnuTildaEff(), dimensioned< Type >::readIfPresent(), and Foam::sqr().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::~SpalartAllmarasDES().


| tmp< volScalarField > DnuTildaEff | ( | ) | const |
Return the effective diffusivity for nuTilda.
Definition at line 390 of file SpalartAllmarasDES.C.
References GeometricField< scalar, fvPatchField, volMesh >::New().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::read(), and SpalartAllmarasDES< BasicMomentumTransportModel >::~SpalartAllmarasDES().


|
virtual |
Return SGS kinetic energy.
Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.
Definition at line 401 of file SpalartAllmarasDES.C.
References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dimLength, forAll, Foam::fvc::grad(), GeometricField< scalar, fvPatchField, volMesh >::New(), nut, patches, patchi, tmp< T >::ref(), GeometricField< Type, PatchField, GeoMesh >::ref(), and Foam::sqr().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::~SpalartAllmarasDES().


|
inline |
Definition at line 196 of file SpalartAllmarasDES.H.
References SpalartAllmarasDES< BasicMomentumTransportModel >::correct(), SpalartAllmarasDES< BasicMomentumTransportModel >::nuTilda_, and SpalartAllmarasDES< BasicMomentumTransportModel >::operator=().

|
virtual |
Correct nuTilda and related properties.
Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.
Definition at line 446 of file SpalartAllmarasDES.C.
References alpha(), Foam::bound(), tmp< T >::clear(), fvConstraints::constrain(), eddyViscosity< LESModel< BasicMomentumTransportModel > >::correct(), Foam::fvm::ddt(), Foam::fvm::div(), fvConstraints, fvModels, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), dictionary::New(), tmp< T >::ref(), rho, Foam::solve(), Foam::fvm::Sp(), and Foam::sqr().
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::nuTilda().


|
delete |
Disallow default bitwise assignment.
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::nuTilda().

|
protected |
Definition at line 70 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 71 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 73 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 74 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 75 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 76 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 77 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 78 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 79 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 80 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 81 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 85 of file SpalartAllmarasDES.H.
Referenced by SpalartAllmarasDES< BasicMomentumTransportModel >::nuTilda().
|
protected |
Wall distance.
Note: different to wall distance in parent RASModel which is for near-wall cells only
Definition at line 90 of file SpalartAllmarasDES.H.