SpalartAllmaras< BasicMomentumTransportModel > Class Template Reference

Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows. More...

Inheritance diagram for SpalartAllmaras< BasicMomentumTransportModel >:
Collaboration diagram for SpalartAllmaras< BasicMomentumTransportModel >:

Public Types

typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::transportModel transportModel
 
- Public Types inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > >
typedef RASModel< BasicMomentumTransportModel > ::alphaField alphaField
 
typedef RASModel< BasicMomentumTransportModel > ::rhoField rhoField
 
typedef RASModel< BasicMomentumTransportModel > ::transportModel transportModel
 
- Public Types inherited from linearViscousStress< RASModel< BasicMomentumTransportModel > >
typedef RASModel< BasicMomentumTransportModel > ::alphaField alphaField
 
typedef RASModel< BasicMomentumTransportModel > ::rhoField rhoField
 
typedef RASModel< BasicMomentumTransportModel > ::transportModel transportModel
 
- Public Types inherited from RASModel< BasicMomentumTransportModel >
typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::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 &type=typeName)
 Construct from components. More...
 
 SpalartAllmaras (const SpalartAllmaras &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~SpalartAllmaras ()
 Destructor. More...
 
virtual bool read ()
 Read RASProperties dictionary. More...
 
tmp< volScalarFieldDnuTildaEff () const
 Return the effective diffusivity for nuTilda. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
void operator= (const SpalartAllmaras &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from eddyViscosity< RASModel< 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< volScalarFieldnut () const
 Return the turbulence viscosity. More...
 
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch. More...
 
virtual tmp< volSymmTensorFieldsigma () 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< RASModel< 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< volSymmTensorFielddevTau () const
 Return the effective stress tensor. More...
 
virtual tmp< fvVectorMatrixdivDevTau (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevTau (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
- Public Member Functions inherited from RASModel< BasicMomentumTransportModel >
 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),(alpha, rho, U, alphaRhoPhi, phi, transport))
 
 RASModel (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...
 
 RASModel (const RASModel &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~RASModel ()
 Destructor. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: small) More...
 
const dimensionedScalarepsilonMin () const
 Return the lower allowable limit for epsilon (default: small) More...
 
const dimensionedScalaromegaMin () const
 Return the lower allowable limit for omega (default: small) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
dimensionedScalarepsilonMin ()
 Allow epsilonMin to be changed. More...
 
dimensionedScalaromegaMin ()
 Allow omegaMin to be changed. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 
void operator= (const RASModel &)=delete
 Disallow default bitwise assignment. 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 volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
 
tmp< volScalarField::Internalfw (const volScalarField::Internal &Stilda) const
 
void correctNut (const volScalarField &fv1)
 
virtual void correctNut ()
 
- Protected Member Functions inherited from RASModel< BasicMomentumTransportModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. 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::Internaly_
 Wall distance. More...
 
- Protected Attributes inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > >
volScalarField nut_
 
- Protected Attributes inherited from RASModel< BasicMomentumTransportModel >
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< BasicMomentumTransportModel >
static autoPtr< RASModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Return a reference to the selected RAS model. More...
 

Detailed Description

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

Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows.

Reference:

    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;
    }
Source files

Definition at line 85 of file SpalartAllmaras.H.

Member Typedef Documentation

◆ alphaField

typedef BasicMomentumTransportModel::alphaField alphaField

Definition at line 151 of file SpalartAllmaras.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 152 of file SpalartAllmaras.H.

◆ transportModel

typedef BasicMomentumTransportModel::transportModel transportModel

Definition at line 153 of file SpalartAllmaras.H.

Constructor & Destructor Documentation

◆ SpalartAllmaras() [1/2]

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

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::correctNut().

Here is the caller graph for this function:

◆ SpalartAllmaras() [2/2]

SpalartAllmaras ( const SpalartAllmaras< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~SpalartAllmaras()

Member Function Documentation

◆ chi()

tmp< volScalarField > chi ( ) const
protected

Definition at line 41 of file SpalartAllmaras.C.

References SpalartAllmaras< BasicMomentumTransportModel >::fv1(), and GeometricField< scalar, fvPatchField, volMesh >::New().

Here is the call graph for this function:

◆ fv1()

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

Definition at line 49 of file SpalartAllmaras.C.

References SpalartAllmaras< BasicMomentumTransportModel >::fv2(), GeometricField< scalar, fvPatchField, volMesh >::New(), and Foam::pow3().

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::chi().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ fv2()

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

Definition at line 60 of file SpalartAllmaras.C.

References DimensionedField< Type, GeoMesh >::New(), and SpalartAllmaras< BasicMomentumTransportModel >::Stilda().

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::fv1().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Stilda()

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

Definition at line 76 of file SpalartAllmaras.C.

References SpalartAllmaras< BasicMomentumTransportModel >::fw(), Foam::fvc::grad(), Foam::mag(), Foam::max(), DimensionedField< Type, GeoMesh >::New(), Foam::skew(), Foam::sqr(), and Foam::sqrt().

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::fv2().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ fw()

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

Definition at line 104 of file SpalartAllmaras.C.

References SpalartAllmaras< BasicMomentumTransportModel >::correctNut(), g, Foam::max(), Foam::min(), DimensionedField< Type, GeoMesh >::New(), Foam::pow(), Foam::pow6(), and Foam::sqr().

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::Stilda().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ correctNut() [1/2]

void correctNut ( const volScalarField fv1)
protected

Definition at line 138 of file SpalartAllmaras.C.

References optionList::correct(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), and options::New().

Here is the call graph for this function:

◆ correctNut() [2/2]

void correctNut ( )
protectedvirtual

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 149 of file SpalartAllmaras.C.

References SpalartAllmaras< BasicMomentumTransportModel >::SpalartAllmaras().

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::fw().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ TypeName()

TypeName ( "SpalartAllmaras< BasicMomentumTransportModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Read RASProperties dictionary.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 280 of file SpalartAllmaras.C.

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

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ DnuTildaEff()

tmp< volScalarField > DnuTildaEff ( ) const

Return the effective diffusivity for nuTilda.

Definition at line 306 of file SpalartAllmaras.C.

References GeometricField< scalar, fvPatchField, volMesh >::New().

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ k()

tmp< volScalarField > k ( ) const
virtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 317 of file SpalartAllmaras.C.

References GeometricField< scalar, fvPatchField, volMesh >::New().

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ epsilon()

tmp< volScalarField > epsilon ( ) const
virtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 330 of file SpalartAllmaras.C.

References Foam::endl(), GeometricField< scalar, fvPatchField, volMesh >::New(), and WarningInFunction.

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ correct()

void correct ( )
virtual

Solve the turbulence equations and correct the turbulence viscosity.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 347 of file SpalartAllmaras.C.

References alpha(), Foam::bound(), optionList::constrain(), correct, optionList::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(), and Foam::sqr().

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

void operator= ( const SpalartAllmaras< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise assignment.

Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().

Here is the caller graph for this function:

Member Data Documentation

◆ sigmaNut_

dimensionedScalar sigmaNut_
protected

Definition at line 100 of file SpalartAllmaras.H.

◆ kappa_

dimensionedScalar kappa_
protected

Definition at line 101 of file SpalartAllmaras.H.

◆ Cb1_

dimensionedScalar Cb1_
protected

Definition at line 103 of file SpalartAllmaras.H.

◆ Cb2_

dimensionedScalar Cb2_
protected

Definition at line 104 of file SpalartAllmaras.H.

◆ Cw1_

dimensionedScalar Cw1_
protected

Definition at line 105 of file SpalartAllmaras.H.

◆ Cw2_

dimensionedScalar Cw2_
protected

Definition at line 106 of file SpalartAllmaras.H.

◆ Cw3_

dimensionedScalar Cw3_
protected

Definition at line 107 of file SpalartAllmaras.H.

◆ Cv1_

dimensionedScalar Cv1_
protected

Definition at line 108 of file SpalartAllmaras.H.

◆ Cs_

dimensionedScalar Cs_
protected

Definition at line 109 of file SpalartAllmaras.H.

◆ nuTilda_

volScalarField nuTilda_
protected

Definition at line 114 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 119 of file SpalartAllmaras.H.


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