Spalart-Allmaras one-transport-equation linear-eddy-viscosity turbulence closure model for incompressible and compressible external flows. More...
#include <SpalartAllmaras.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > > | |
| typedef eddyViscosity< RASModel< BasicTurbulenceModel > >::alphaField | alphaField |
| typedef eddyViscosity< RASModel< BasicTurbulenceModel > >::rhoField | rhoField |
| typedef eddyViscosity< RASModel< BasicTurbulenceModel > >::transportModel | transportModel |
| Public Types inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from linearViscousStress< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("SpalartAllmaras") | |
| Runtime type information. | |
| 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. | |
| virtual | ~SpalartAllmaras ()=default |
| Destructor. | |
| Public Member Functions inherited from SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > > | |
| virtual | ~SpalartAllmarasBase ()=default |
| Destructor. | |
| virtual bool | read () |
| Re-read model coefficients if they have changed. | |
| tmp< volScalarField > | DnuTildaEff () const |
| Return the effective diffusivity for nuTilda. | |
| virtual tmp< volScalarField > | k () const |
| Return the (estimated) turbulent kinetic energy. | |
| virtual tmp< volScalarField > | epsilon () const |
| Return the (estimated) turbulent kinetic energy dissipation rate. | |
| virtual tmp< volScalarField > | omega () const |
| Return the (estimated) specific dissipation rate. | |
| tmp< volScalarField > | nuTilda () const |
| Return the modified kinematic viscosity. | |
| virtual void | correct () |
| Correct nuTilda and related properties. | |
| 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. | |
| virtual | ~eddyViscosity ()=default |
| Destructor. | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity. | |
| virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. | |
| virtual void | validate () |
| Validate the turbulence fields after construction. | |
| Public Member Functions inherited from linearViscousStress< 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. | |
| virtual | ~linearViscousStress ()=default |
| Destructor. | |
| virtual tmp< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
| Return the effective stress tensor based on a given velocity field. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. | |
Additional Inherited Members | |
| Protected Attributes inherited from SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > > | |
| dimensionedScalar | sigmaNut_ |
| dimensionedScalar | kappa_ |
| dimensionedScalar | Cb1_ |
| dimensionedScalar | Cb2_ |
| dimensionedScalar | Cw1_ |
| dimensionedScalar | Cw2_ |
| dimensionedScalar | Cw3_ |
| dimensionedScalar | Cv1_ |
| dimensionedScalar | Cs_ |
| dimensionedScalar | ck_ |
| Switch | ft2_ |
| dimensionedScalar | Ct3_ |
| dimensionedScalar | Ct4_ |
| volScalarField | nuTilda_ |
| Modified kinematic viscosity [m^2/s]. | |
| const volScalarField & | y_ |
| Wall distance. | |
| Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| volScalarField | nut_ |
Spalart-Allmaras one-transport-equation linear-eddy-viscosity turbulence closure model for incompressible and compressible external flows.
Required fields
nuTilda | Modified kinematic viscosity [m2/s]
References:
Standard model:
Spalart, P.R., & Allmaras, S.R. (1994).
A one-equation turbulence model for aerodynamic flows.
La Recherche Aerospatiale, 1, 5-21.
Standard model without trip and ft2 terms (tag:R):
Rumsey, C. (2020).
The Spalart-Allmaras Turbulence Model.
Spalart-Allmaras One-Equation Model without ft2 Term (SA-noft2).
https://turbmodels.larc.nasa.gov/spalart.html#sanoft2
(Retrieved:12-01-2021).
constant/turbulenceProperties: RAS
{
// Mandatory entries (unmodifiable)
RASModel SpalartAllmaras;
// Optional entries (runtime modifiable)
turbulence on;
printCoeffs on;
SpalartAllmarasCoeffs
{
sigmaNut 0.66666;
kappa 0.41;
Cb1 0.1355;
Cb2 0.622;
Cw2 0.3;
Cw3 2.0;
Cv1 7.1;
Cs 0.3;
}
}
Stilda generation term should never be allowed to be zero or negative to avoid potential numerical issues and unphysical results for complex flows. To this end, a limiter proposed by Spalart (R:Note-1(b)) is applied onto Stilda where Stilda is clipped at Cs*Omega with the default value of Cs=0.3.k, epsilon or omega. Nevertheless, these quantities can be estimated by using an approximate expressions for turbulent kinetic energy and dissipation rate reported in (B:Eq. 4.50).Definition at line 115 of file SpalartAllmaras.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 151 of file SpalartAllmaras.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 152 of file SpalartAllmaras.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 153 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 60 of file SpalartAllmaras.C.
References alpha, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, phi, rho, and U.

|
virtualdefault |
Destructor.
|
protectedvirtual |
Return the length scale.
Implements SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >.
Definition at line 34 of file SpalartAllmaras.C.
References SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::chi(), SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::fv1(), and SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::y_.

|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 46 of file SpalartAllmaras.C.
References SpalartAllmarasBase< BasicEddyViscosityModel >::correctNut(), and eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity().

| TypeName | ( | "SpalartAllmaras< BasicTurbulenceModel >" | ) |
Runtime type information.