SpalartAllmaras< BasicTurbulenceModel > Class Template Reference

Spalart-Allmaras one-transport-equation linear-eddy-viscosity turbulence closure model for incompressible and compressible external flows. More...

Inheritance diagram for SpalartAllmaras< BasicTurbulenceModel >:
[legend]
Collaboration diagram for SpalartAllmaras< BasicTurbulenceModel >:
[legend]

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef 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. 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< volScalarFieldDnuTildaEff () const
 Return the effective diffusivity for nuTilda. More...
 
virtual tmp< volScalarFieldk () const
 Return the (estimated) turbulent kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the (estimated) turbulent kinetic energy dissipation rate. More...
 
virtual tmp< volScalarFieldomega () const
 Return the (estimated) specific dissipation rate. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulent 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 bool read ()=0
 Re-read model coefficients if they have changed. More...
 
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity. More...
 
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch. More...
 
virtual tmp< volScalarFieldk () const=0
 Return the turbulence kinetic energy. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
- 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. More...
 
virtual ~linearViscousStress ()=default
 Destructor. More...
 
virtual bool read ()=0
 Re-read model coefficients if they have changed. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity. More...
 

Protected Member Functions

tmp< volScalarFieldchi () const
 
tmp< volScalarFieldfv1 (const volScalarField &chi) const
 
tmp< volScalarField::Internalfv2 (const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
 
tmp< volScalarField::InternalStilda () const
 
tmp< volScalarField::Internalfw (const volScalarField::Internal &Stilda) const
 
virtual void correctNut ()
 Update nut with the latest available nuTilda. More...
 
virtual void correctNut ()=0
 

Protected Attributes

dimensionedScalar sigmaNut_
 
dimensionedScalar kappa_
 
dimensionedScalar Cb1_
 
dimensionedScalar Cb2_
 
dimensionedScalar Cw1_
 
dimensionedScalar Cw2_
 
dimensionedScalar Cw3_
 
dimensionedScalar Cv1_
 
dimensionedScalar Cs_
 
volScalarField nuTilda_
 Modified kinematic viscosity [m2/s]. More...
 
const volScalarField::Internaly_
 Wall distance. More...
 
- Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_
 

Detailed Description

template<class BasicTurbulenceModel>
class Foam::RASModels::SpalartAllmaras< BasicTurbulenceModel >

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).

    Estimation expression for k and epsilon (tag:B), Eq. 4.50:
        Bourgoin, A. (2019).
        Bathymetry induced turbulence modelling the
        Alderney Race site: regional approach with TELEMAC-LES.
        Normandie Université.

    Estimation expressions for omega (tag:P):
        Pope, S. B. (2000).
        Turbulent flows.
        Cambridge, UK: Cambridge Univ. Press
        DOI:10.1017/CBO9780511840531
Usage
Example by using 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;
    }
}
Note
  • The model is implemented without the trip-term since the model has almost always been used in fully turbulent applications rather than those where laminar-turbulent transition occurs.
  • It has been argued that the ft2 term is not needed in the absence of the trip-term, hence ft2 term is also not implementated.
  • The 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.
  • The model does not produce 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).
Source files

Definition at line 132 of file SpalartAllmaras.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 199 of file SpalartAllmaras.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 200 of file SpalartAllmaras.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 201 of file SpalartAllmaras.H.

Constructor & Destructor Documentation

◆ SpalartAllmaras()

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 145 of file SpalartAllmaras.C.

References Foam::type().

Here is the call graph for this function:

◆ ~SpalartAllmaras()

virtual ~SpalartAllmaras ( )
virtualdefault

Destructor.

Member Function Documentation

◆ chi()

tmp< volScalarField > chi
protected

Definition at line 44 of file SpalartAllmaras.C.

References nu.

◆ fv1()

tmp< volScalarField > fv1 ( const volScalarField chi) const
protected

Definition at line 51 of file SpalartAllmaras.C.

References Foam::pow3().

Here is the call graph for this function:

◆ fv2()

tmp< volScalarField::Internal > fv2 ( const volScalarField::Internal chi,
const volScalarField::Internal fv1 
) const
protected

Definition at line 63 of file SpalartAllmaras.C.

◆ Stilda()

tmp< volScalarField::Internal > Stilda
protected

Definition at line 74 of file SpalartAllmaras.C.

References Foam::fvc::grad(), Foam::mag(), Foam::max(), Foam::skew(), Foam::sqr(), and Foam::sqrt().

Here is the call graph for this function:

◆ fw()

tmp< volScalarField::Internal > fw ( const volScalarField::Internal Stilda) const
protected

Definition at line 98 of file SpalartAllmaras.C.

References DimensionedField< Type, GeoMesh >::dimensions(), g, Foam::max(), Foam::min(), Foam::pow(), Foam::pow6(), and Foam::sqr().

Here is the call graph for this function:

◆ correctNut()

void correctNut
protectedvirtual

Update nut with the latest available nuTilda.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 132 of file SpalartAllmaras.C.

References Time::New().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "SpalartAllmaras< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 269 of file SpalartAllmaras.C.

References Foam::read(), dimensioned< Type >::readIfPresent(), and Foam::sqr().

Here is the call graph for this function:

◆ DnuTildaEff()

tmp< volScalarField > DnuTildaEff

Return the effective diffusivity for nuTilda.

Definition at line 292 of file SpalartAllmaras.C.

References Time::New(), and nu.

Here is the call graph for this function:

◆ k()

tmp< volScalarField > k
virtual

Return the (estimated) turbulent kinetic energy.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 303 of file SpalartAllmaras.C.

References Foam::cbrt(), Foam::fvc::grad(), IOobject::groupName(), Foam::mag(), Time::New(), Foam::sqrt(), and Foam::symm().

Here is the call graph for this function:

◆ epsilon()

tmp< volScalarField > epsilon
virtual

Return the (estimated) turbulent kinetic energy dissipation rate.

Definition at line 326 of file SpalartAllmaras.C.

References Foam::dimLength, Foam::dimTime, IOobject::groupName(), k, Time::New(), Foam::pow(), Foam::sqr(), and Foam::sqrt().

Here is the call graph for this function:

◆ omega()

tmp< volScalarField > omega
virtual

Return the (estimated) specific dissipation rate.

Definition at line 349 of file SpalartAllmaras.C.

References Foam::dimLength, Foam::dimTime, IOobject::groupName(), Time::New(), and Foam::sqr().

Here is the call graph for this function:

◆ correct()

void correct
virtual

Solve the turbulence equations and correct the turbulent viscosity.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 370 of file SpalartAllmaras.C.

References alpha, Foam::bound(), correct(), Foam::fvm::ddt(), Foam::fvm::div(), fvOptions, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), Time::New(), tmp< T >::ref(), rho, solve(), Foam::fvm::Sp(), Foam::sqr(), and Foam::Zero.

Here is the call graph for this function:

Member Data Documentation

◆ sigmaNut_

dimensionedScalar sigmaNut_
protected

Definition at line 151 of file SpalartAllmaras.H.

◆ kappa_

dimensionedScalar kappa_
protected

Definition at line 152 of file SpalartAllmaras.H.

◆ Cb1_

dimensionedScalar Cb1_
protected

Definition at line 153 of file SpalartAllmaras.H.

◆ Cb2_

dimensionedScalar Cb2_
protected

Definition at line 154 of file SpalartAllmaras.H.

◆ Cw1_

dimensionedScalar Cw1_
protected

Definition at line 155 of file SpalartAllmaras.H.

◆ Cw2_

dimensionedScalar Cw2_
protected

Definition at line 156 of file SpalartAllmaras.H.

◆ Cw3_

dimensionedScalar Cw3_
protected

Definition at line 157 of file SpalartAllmaras.H.

◆ Cv1_

dimensionedScalar Cv1_
protected

Definition at line 158 of file SpalartAllmaras.H.

◆ Cs_

dimensionedScalar Cs_
protected

Definition at line 159 of file SpalartAllmaras.H.

◆ nuTilda_

volScalarField nuTilda_
protected

Modified kinematic viscosity [m2/s].

Definition at line 165 of file SpalartAllmaras.H.

◆ y_

const volScalarField::Internal& y_
protected

Wall distance.

Note: different to wall distance in parent RASModel which is for near-wall cells only

Definition at line 170 of file SpalartAllmaras.H.


The documentation for this class was generated from the following files: