SpalartAllmarasDES DES turbulence model for incompressible and compressible flows. More...
Public Types | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef LESModel< BasicTurbulenceModel > ::alphaField | alphaField |
typedef LESModel< BasicTurbulenceModel > ::rhoField | rhoField |
typedef LESModel< BasicTurbulenceModel > ::transportModel | transportModel |
![]() | |
typedef LESModel< BasicTurbulenceModel > ::alphaField | alphaField |
typedef LESModel< BasicTurbulenceModel > ::rhoField | rhoField |
typedef LESModel< BasicTurbulenceModel > ::transportModel | transportModel |
![]() | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::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 &propertiesName=turbulenceModel::propertiesName, const word &type=typeName) | |
Construct from components. More... | |
virtual | ~SpalartAllmarasDES () |
Destructor. More... | |
virtual bool | read () |
Re-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 tmp< volScalarField > | LESRegion () const |
Return the LES field indicator. More... | |
virtual void | correct () |
Correct nuTilda and related properties. More... | |
![]() | |
DESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) | |
Construct from components. More... | |
virtual | ~DESModel () |
Destructor. More... | |
![]() | |
DESModelBase () | |
Constructor. More... | |
virtual | ~DESModelBase () |
Destructor. More... | |
ClassName ("DESModelBase") | |
![]() | |
LESeddyViscosity (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) | |
Construct from components. More... | |
virtual | ~LESeddyViscosity () |
Destructor. More... | |
virtual tmp< volScalarField > | epsilon () const |
Return sub-grid disipation rate. More... | |
![]() | |
eddyViscosity (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. More... | |
virtual | ~eddyViscosity () |
Destructor. More... | |
volScalarField & | evNut () |
Return non-const access to the turbulence viscosity. 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< volScalarField > | k () const=0 |
Return the turbulence kinetic energy. More... | |
virtual tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
![]() | |
linearViscousStress (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. More... | |
virtual | ~linearViscousStress () |
Destructor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
![]() | |
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, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)) | |
LESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
Construct from components. 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... | |
Private Member Functions | |
SpalartAllmarasDES (const SpalartAllmarasDES &) | |
SpalartAllmarasDES & | operator= (const SpalartAllmarasDES &) |
Additional Inherited Members | |
![]() | |
static autoPtr< LESModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) |
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.
Including the low Reynolds number correction documented in
Spalart, P. R., Deck, S., Shur, M.L., Squires, K.D., Strelets, M.Kh, Travin, A. (2006). A new version of detached-eddy simulation, resistant to ambiguous grid densities. Theor. Comput. Fluid Dyn., 20, 181-195.
Definition at line 72 of file SpalartAllmarasDES.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 178 of file SpalartAllmarasDES.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 179 of file SpalartAllmarasDES.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 180 of file SpalartAllmarasDES.H.
|
private |
SpalartAllmarasDES | ( | 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 |
||
) |
Construct from components.
Definition at line 236 of file SpalartAllmarasDES.C.
References Foam::type().
|
inlinevirtual |
Destructor.
Definition at line 204 of file SpalartAllmarasDES.H.
|
private |
|
protected |
Definition at line 38 of file SpalartAllmarasDES.C.
References nu.
|
protected |
Definition at line 46 of file SpalartAllmarasDES.C.
References Foam::pow3().
|
protected |
Definition at line 57 of file SpalartAllmarasDES.C.
|
protected |
Definition at line 68 of file SpalartAllmarasDES.C.
References Foam::exp(), and Foam::sqr().
|
protected |
Definition at line 78 of file SpalartAllmarasDES.C.
References Foam::mag(), Foam::skew(), and Foam::sqrt().
|
protected |
Definition at line 88 of file SpalartAllmarasDES.C.
References Foam::max(), and Foam::sqr().
|
protected |
Definition at line 109 of file SpalartAllmarasDES.C.
References Foam::max(), Foam::min(), Foam::sqr(), and Foam::tr().
|
protected |
Definition at line 139 of file SpalartAllmarasDES.C.
References g, Foam::pow(), and Foam::pow6().
|
protected |
Definition at line 153 of file SpalartAllmarasDES.C.
References Foam::dimless, Foam::e, Foam::max(), mesh, Foam::min(), IOobject::NO_READ, IOobject::NO_WRITE, psi, Foam::sqr(), Foam::sqrt(), timeName, and Foam::type().
|
protectedvirtual |
Length scale.
Reimplemented in SpalartAllmarasIDDES< BasicTurbulenceModel >, and SpalartAllmarasDDES< BasicTurbulenceModel >.
Definition at line 200 of file SpalartAllmarasDES.C.
References delta, dimensionedInternalField(), Foam::min(), and psi.
|
protected |
Definition at line 214 of file SpalartAllmarasDES.C.
References GeometricField::correctBoundaryConditions().
|
protectedvirtual |
Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.
Definition at line 226 of file SpalartAllmarasDES.C.
TypeName | ( | "SpalartAllmarasDES< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented from LESeddyViscosity< BasicTurbulenceModel >.
Reimplemented in SpalartAllmarasIDDES< BasicTurbulenceModel >, and SpalartAllmarasDDES< BasicTurbulenceModel >.
Definition at line 412 of file SpalartAllmarasDES.C.
References dimensioned::readIfPresent(), and Foam::sqr().
tmp< volScalarField > DnuTildaEff |
Return the effective diffusivity for nuTilda.
Definition at line 446 of file SpalartAllmarasDES.C.
References nu.
|
virtual |
Return SGS kinetic energy.
Definition at line 456 of file SpalartAllmarasDES.C.
References Foam::fvc::grad(), nut, and Foam::sqr().
|
inline |
Definition at line 219 of file SpalartAllmarasDES.H.
References SpalartAllmarasDES< BasicTurbulenceModel >::nuTilda_.
|
virtual |
Return the LES field indicator.
Implements DESModel< BasicTurbulenceModel >.
Definition at line 465 of file SpalartAllmarasDES.C.
References Foam::fvc::grad(), Foam::neg(), IOobject::NO_READ, and IOobject::NO_WRITE.
|
virtual |
Correct nuTilda and related properties.
Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.
Definition at line 491 of file SpalartAllmarasDES.C.
References Foam::constant::atomic::alpha, Foam::bound(), eddyViscosity< LESModel< BasicTurbulenceModel > >::correct(), Foam::fvm::ddt(), Foam::fvm::div(), Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), rho, Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), and U.
|
protected |
Definition at line 89 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 90 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 92 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 93 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 94 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 95 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 96 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 97 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 98 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 99 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 100 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 105 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 106 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 107 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 108 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 113 of file SpalartAllmarasDES.H.
Referenced by SpalartAllmarasDES< BasicTurbulenceModel >::nuTilda().
|
protected |
Wall distance.
Note: different to wall distance in parent RASModel which is for near-wall cells only
Definition at line 118 of file SpalartAllmarasDES.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.