Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows. More...
Public Types | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Types inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > > | |
typedef RASModel< BasicMomentumTransportModel > ::alphaField | alphaField |
typedef RASModel< BasicMomentumTransportModel > ::rhoField | rhoField |
Public Types inherited from linearViscousStress< RASModel< BasicMomentumTransportModel > > | |
typedef RASModel< BasicMomentumTransportModel > ::alphaField | alphaField |
typedef RASModel< BasicMomentumTransportModel > ::rhoField | rhoField |
Public Types inherited from RASModel< BasicMomentumTransportModel > | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
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 viscosity &viscosity, 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< 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 turbulence specific 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 viscosity &viscosity) | |
Construct from components. More... | |
virtual | ~eddyViscosity () |
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 > | sigma () 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 viscosity &viscosity) | |
Construct from components. More... | |
virtual | ~linearViscousStress () |
Destructor. More... | |
virtual tmp< volSymmTensorField > | devTau () const |
Return the effective stress tensor. More... | |
virtual tmp< fvVectorMatrix > | divDevTau (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevTau (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 viscosity &viscosity),(alpha, rho, U, alphaRhoPhi, phi, viscosity)) | |
RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
Construct from components. More... | |
RASModel (const RASModel &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~RASModel () |
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 > | nu () const |
Return the laminar viscosity. More... | |
virtual tmp< scalarField > | nu (const label patchi) const |
Return the laminar viscosity on patchi. 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... | |
void | operator= (const RASModel &)=delete |
Disallow default bitwise assignment. More... | |
Protected Member Functions | |
tmp< volScalarField > | chi () const |
tmp< volScalarField > | fv1 (const volScalarField &chi) const |
tmp< volScalarField::Internal > | fv2 (const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const |
tmp< volScalarField::Internal > | Stilda (const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const |
tmp< volScalarField::Internal > | fw (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... | |
Additional Inherited Members | |
Static Public Member Functions inherited from RASModel< BasicMomentumTransportModel > | |
static autoPtr< RASModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) |
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 85 of file SpalartAllmaras.H.
typedef BasicMomentumTransportModel::alphaField alphaField |
Definition at line 151 of file SpalartAllmaras.H.
typedef BasicMomentumTransportModel::rhoField rhoField |
Definition at line 152 of file SpalartAllmaras.H.
SpalartAllmaras | ( | const alphaField & | alpha, |
const rhoField & | rho, | ||
const volVectorField & | U, | ||
const surfaceScalarField & | alphaRhoPhi, | ||
const surfaceScalarField & | phi, | ||
const viscosity & | viscosity, | ||
const word & | type = typeName |
||
) |
Construct from components.
Definition at line 160 of file SpalartAllmaras.C.
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::correctNut().
|
delete |
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 178 of file SpalartAllmaras.H.
References SpalartAllmaras< BasicMomentumTransportModel >::correct(), SpalartAllmaras< BasicMomentumTransportModel >::DnuTildaEff(), SpalartAllmaras< BasicMomentumTransportModel >::epsilon(), SpalartAllmaras< BasicMomentumTransportModel >::k(), SpalartAllmaras< BasicMomentumTransportModel >::omega(), SpalartAllmaras< BasicMomentumTransportModel >::operator=(), and SpalartAllmaras< BasicMomentumTransportModel >::read().
|
protected |
Definition at line 42 of file SpalartAllmaras.C.
References SpalartAllmaras< BasicMomentumTransportModel >::fv1(), and GeometricField< scalar, fvPatchField, volMesh >::New().
|
protected |
Definition at line 50 of file SpalartAllmaras.C.
References SpalartAllmaras< BasicMomentumTransportModel >::fv2(), GeometricField< scalar, fvPatchField, volMesh >::New(), and Foam::pow3().
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::chi().
|
protected |
Definition at line 61 of file SpalartAllmaras.C.
References DimensionedField< Type, GeoMesh >::New(), and SpalartAllmaras< BasicMomentumTransportModel >::Stilda().
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::fv1().
|
protected |
Definition at line 77 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().
|
protected |
Definition at line 105 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().
|
protected |
Definition at line 139 of file SpalartAllmaras.C.
References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), and dictionary::New().
|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 150 of file SpalartAllmaras.C.
References SpalartAllmaras< BasicMomentumTransportModel >::SpalartAllmaras().
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::fw().
TypeName | ( | "SpalartAllmaras< BasicMomentumTransportModel >" | ) |
Runtime type information.
|
virtual |
Read RASProperties dictionary.
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 281 of file SpalartAllmaras.C.
References Foam::read(), dimensioned< Type >::readIfPresent(), and Foam::sqr().
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().
tmp< volScalarField > DnuTildaEff | ( | ) | const |
Return the effective diffusivity for nuTilda.
Definition at line 307 of file SpalartAllmaras.C.
References GeometricField< scalar, fvPatchField, volMesh >::New().
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().
|
virtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 318 of file SpalartAllmaras.C.
References GeometricField< scalar, fvPatchField, volMesh >::New().
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().
|
virtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 331 of file SpalartAllmaras.C.
References Foam::endl(), GeometricField< scalar, fvPatchField, volMesh >::New(), and WarningInFunction.
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().
|
virtual |
Return the turbulence specific dissipation rate.
Definition at line 349 of file SpalartAllmaras.C.
References Foam::dimless, Foam::dimTime, Foam::endl(), GeometricField< scalar, fvPatchField, volMesh >::New(), and WarningInFunction.
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 366 of file SpalartAllmaras.C.
References alpha(), Foam::bound(), fvConstraints::constrain(), correct, Foam::fvm::ddt(), Foam::fvm::div(), fvConstraints, fvModels, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), dictionary::New(), tmp< T >::ref(), rho, Foam::solve(), Foam::fvm::Sp(), and Foam::sqr().
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().
|
delete |
Disallow default bitwise assignment.
Referenced by SpalartAllmaras< BasicMomentumTransportModel >::~SpalartAllmaras().
|
protected |
Definition at line 100 of file SpalartAllmaras.H.
|
protected |
Definition at line 101 of file SpalartAllmaras.H.
|
protected |
Definition at line 103 of file SpalartAllmaras.H.
|
protected |
Definition at line 104 of file SpalartAllmaras.H.
|
protected |
Definition at line 105 of file SpalartAllmaras.H.
|
protected |
Definition at line 106 of file SpalartAllmaras.H.
|
protected |
Definition at line 107 of file SpalartAllmaras.H.
|
protected |
Definition at line 108 of file SpalartAllmaras.H.
|
protected |
Definition at line 109 of file SpalartAllmaras.H.
|
protected |
Definition at line 114 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 119 of file SpalartAllmaras.H.