Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows. More...
Public Types | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
Public Types inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
typedef RASModel< BasicTurbulenceModel > ::alphaField | alphaField |
typedef RASModel< BasicTurbulenceModel > ::rhoField | rhoField |
typedef RASModel< BasicTurbulenceModel > ::transportModel | transportModel |
Public Types inherited from linearViscousStress< RASModel< BasicTurbulenceModel > > | |
typedef RASModel< BasicTurbulenceModel > ::alphaField | alphaField |
typedef RASModel< BasicTurbulenceModel > ::rhoField | rhoField |
typedef RASModel< BasicTurbulenceModel > ::transportModel | transportModel |
Public Types inherited from RASModel< BasicTurbulenceModel > | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
TypeName ("SpalartAllmaras") | |
Runtime type information. More... | |
SpalartAllmaras (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 | ~SpalartAllmaras ()=default |
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 the turbulence kinetic energy. More... | |
virtual tmp< volScalarField > | epsilon () const |
Return the turbulence kinetic energy dissipation rate. More... | |
virtual tmp< volScalarField > | omega () const |
Return the (estimated) specific dissipation rate. More... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
Public Member Functions inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
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 ()=default |
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 > | R () const |
Return the Reynolds stress tensor. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
Public Member Functions inherited from linearViscousStress< RASModel< BasicTurbulenceModel > > | |
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 ()=default |
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... | |
Public Member Functions inherited from RASModel< BasicTurbulenceModel > | |
TypeName ("RAS") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, RASModel, 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)) | |
RASModel (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 | ~RASModel ()=default |
Destructor. More... | |
const dimensionedScalar & | kMin () const |
Return the lower allowable limit for k (default: SMALL) More... | |
const dimensionedScalar & | epsilonMin () const |
Return the lower allowable limit for epsilon (default: SMALL) More... | |
const dimensionedScalar & | omegaMin () const |
Return the lower allowable limit for omega (default: SMALL) More... | |
dimensionedScalar & | kMin () |
Allow kMin to be changed. More... | |
dimensionedScalar & | epsilonMin () |
Allow epsilonMin to be changed. More... | |
dimensionedScalar & | omegaMin () |
Allow omegaMin to be changed. More... | |
virtual const dictionary & | coeffDict () const |
Const access to the coefficients dictionary. 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... | |
Protected Member Functions | |
tmp< volScalarField > | chi () const |
tmp< volScalarField > | fv1 (const volScalarField &chi) const |
tmp< volScalarField > | fv2 (const volScalarField &chi, const volScalarField &fv1) const |
tmp< volScalarField > | Stilda (const volScalarField &chi, const volScalarField &fv1) const |
tmp< volScalarField > | fw (const volScalarField &Stilda) const |
void | correctNut (const volScalarField &fv1) |
virtual void | correctNut () |
Protected Member Functions inherited from RASModel< BasicTurbulenceModel > | |
virtual void | printCoeffs (const word &type) |
Print model coefficients. More... | |
RASModel (const RASModel &)=delete | |
No copy construct. More... | |
void | operator= (const RASModel &)=delete |
No copy assignment. More... | |
Protected Attributes | |
dimensionedScalar | sigmaNut_ |
dimensionedScalar | kappa_ |
dimensionedScalar | Cb1_ |
dimensionedScalar | Cb2_ |
dimensionedScalar | Cw1_ |
dimensionedScalar | Cw2_ |
dimensionedScalar | Cw3_ |
dimensionedScalar | Cv1_ |
dimensionedScalar | Cs_ |
volScalarField | nuTilda_ |
const volScalarField & | y_ |
Wall distance. More... | |
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
volScalarField | nut_ |
Protected Attributes inherited from RASModel< BasicTurbulenceModel > | |
dictionary | RASDict_ |
RAS coefficients dictionary. More... | |
Switch | turbulence_ |
Turbulence on/off flag. More... | |
Switch | printCoeffs_ |
Flag to print the model coeffs at run-time. More... | |
dictionary | coeffDict_ |
Model coefficients dictionary. More... | |
dimensionedScalar | kMin_ |
Lower limit of k. More... | |
dimensionedScalar | epsilonMin_ |
Lower limit of epsilon. More... | |
dimensionedScalar | omegaMin_ |
Lower limit for omega. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from RASModel< BasicTurbulenceModel > | |
static autoPtr< RASModel > | 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 RAS model. More... | |
Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows.
Spalart, P.R., & Allmaras, S.R. (1994). A one-equation turbulence model for aerodynamic flows. La Recherche Aerospatiale, 1, 5-21.
The model is implemented without the trip-term and hence the ft2 term is not needed.
It is necessary to limit the Stilda generation term as the model generates unphysical results if this term becomes negative which occurs for complex flow. Several approaches have been proposed to limit Stilda but it is not clear which is the most appropriate. Here the limiter proposed by Spalart is implemented in which Stilda is clipped at Cs*Omega with the default value of Cs = 0.3.
The default model coefficients are
SpalartAllmarasCoeffs { Cb1 0.1355; Cb2 0.622; Cw2 0.3; Cw3 2.0; Cv1 7.1; Cs 0.3; sigmaNut 0.66666; kappa 0.41; }
Definition at line 91 of file SpalartAllmaras.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 158 of file SpalartAllmaras.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 159 of file SpalartAllmaras.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 160 of file SpalartAllmaras.H.
SpalartAllmaras | ( | 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 148 of file SpalartAllmaras.C.
References Foam::type().
|
virtualdefault |
Destructor.
|
protected |
Definition at line 44 of file SpalartAllmaras.C.
References nu.
|
protected |
Definition at line 52 of file SpalartAllmaras.C.
References Foam::pow3().
|
protected |
Definition at line 63 of file SpalartAllmaras.C.
|
protected |
Definition at line 74 of file SpalartAllmaras.C.
References Foam::fvc::grad(), Foam::mag(), Foam::max(), Foam::skew(), Foam::sqr(), and Foam::sqrt().
|
protected |
Definition at line 95 of file SpalartAllmaras.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), g, Foam::max(), Foam::min(), Foam::pow(), Foam::pow6(), and Foam::sqr().
|
protected |
Definition at line 125 of file SpalartAllmaras.C.
References optionList::correct(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), and options::New().
|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 138 of file SpalartAllmaras.C.
TypeName | ( | "SpalartAllmaras< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 271 of file SpalartAllmaras.C.
References Foam::read(), dimensioned< Type >::readIfPresent(), and Foam::sqr().
tmp< volScalarField > DnuTildaEff | ( | ) | const |
Return the effective diffusivity for nuTilda.
Definition at line 294 of file SpalartAllmaras.C.
References nu.
|
virtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 304 of file SpalartAllmaras.C.
References Foam::dimLength, Foam::dimTime, Foam::endl(), tmp< T >::New(), Foam::sqr(), WarningInFunction, and Foam::Zero.
|
virtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 327 of file SpalartAllmaras.C.
References Foam::dimLength, Foam::dimTime, Foam::endl(), tmp< T >::New(), Foam::pow3(), Foam::sqr(), WarningInFunction, and Foam::Zero.
|
virtual |
Return the (estimated) specific dissipation rate.
Definition at line 350 of file SpalartAllmaras.C.
References Foam::dimless, Foam::dimTime, Foam::endl(), IOobject::groupName(), tmp< T >::New(), WarningInFunction, and Foam::Zero.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 372 of file SpalartAllmaras.C.
References Foam::constant::atomic::alpha, Foam::bound(), correct(), Foam::fvm::ddt(), Foam::fvm::div(), fvOptions, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), options::New(), tmp< T >::ref(), rho, Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), and Foam::Zero.
|
protected |
Definition at line 110 of file SpalartAllmaras.H.
|
protected |
Definition at line 111 of file SpalartAllmaras.H.
|
protected |
Definition at line 113 of file SpalartAllmaras.H.
|
protected |
Definition at line 114 of file SpalartAllmaras.H.
|
protected |
Definition at line 115 of file SpalartAllmaras.H.
|
protected |
Definition at line 116 of file SpalartAllmaras.H.
|
protected |
Definition at line 117 of file SpalartAllmaras.H.
|
protected |
Definition at line 118 of file SpalartAllmaras.H.
|
protected |
Definition at line 119 of file SpalartAllmaras.H.
|
protected |
Definition at line 124 of file SpalartAllmaras.H.
|
protected |
Wall distance.
Note: different to wall distance in parent RASModel which is for near-wall cells only
Definition at line 129 of file SpalartAllmaras.H.