Langtry-Menter 4-equation transitional SST model based on the k-omega-SST RAS model. More...
Public Types | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
typedef BasicMomentumTransportModel::transportModel | transportModel |
![]() | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
typedef BasicMomentumTransportModel::transportModel | transportModel |
![]() | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
typedef BasicMomentumTransportModel::transportModel | transportModel |
![]() | |
typedef RASModel< BasicMomentumTransportModel > ::alphaField | alphaField |
typedef RASModel< BasicMomentumTransportModel > ::rhoField | rhoField |
typedef RASModel< BasicMomentumTransportModel > ::transportModel | transportModel |
![]() | |
typedef RASModel< BasicMomentumTransportModel > ::alphaField | alphaField |
typedef RASModel< BasicMomentumTransportModel > ::rhoField | rhoField |
typedef RASModel< BasicMomentumTransportModel > ::transportModel | transportModel |
![]() | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
typedef BasicMomentumTransportModel::transportModel | transportModel |
Public Member Functions | |
TypeName ("kOmegaSSTLM") | |
Runtime type information. More... | |
kOmegaSSTLM (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... | |
kOmegaSSTLM (const kOmegaSSTLM &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~kOmegaSSTLM () |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. More... | |
const volScalarField & | ReThetat () const |
Access function transition onset momentum-thickness Reynolds number. More... | |
const volScalarField & | gammaInt () const |
Access function to intermittency. More... | |
tmp< volScalarField > | DReThetatEff () const |
Return the effective diffusivity for transition onset. More... | |
tmp< volScalarField > | DgammaIntEff () const |
Return the effective diffusivity for intermittency. More... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
void | operator= (const kOmegaSSTLM &)=delete |
Disallow default bitwise assignment. More... | |
![]() | |
TypeName ("kOmegaSST") | |
Runtime type information. More... | |
kOmegaSST (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... | |
virtual | ~kOmegaSST () |
Destructor. More... | |
![]() | |
kOmegaSST (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... | |
kOmegaSST (const kOmegaSST &)=delete | |
Disallow default bitwise copy construction. More... | |
tmp< volScalarField > | DkEff (const volScalarField &F1) const |
Return the effective diffusivity for k. More... | |
tmp< volScalarField > | DomegaEff (const volScalarField &F1) const |
Return the effective diffusivity for omega. 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 kinetic energy dissipation rate. More... | |
void | operator= (const kOmegaSST &)=delete |
Disallow default bitwise assignment. More... | |
![]() | |
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< 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... | |
![]() | |
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< 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... | |
![]() | |
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 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... | |
void | operator= (const RASModel &)=delete |
Disallow default bitwise assignment. More... | |
Additional Inherited Members | |
![]() | |
static autoPtr< RASModel > | New (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... | |
Langtry-Menter 4-equation transitional SST model based on the k-omega-SST RAS model.
Langtry, R. B., & Menter, F. R. (2009). Correlation-based transition modeling for unstructured parallelized computational fluid dynamics codes. AIAA journal, 47(12), 2894-2906. Menter, F. R., Langtry, R., & Volker, S. (2006). Transition modelling for general purpose CFD codes. Flow, turbulence and combustion, 77(1-4), 277-303. Langtry, R. B. (2006). A correlation-based transition model using local variables for unstructured parallelized CFD codes. Phd. Thesis, Universität Stuttgart.
The model coefficients are
kOmegaSSTCoeffs { // Default SST coefficients alphaK1 0.85; alphaK2 1; alphaOmega1 0.5; alphaOmega2 0.856; beta1 0.075; beta2 0.0828; betaStar 0.09; gamma1 5/9; gamma2 0.44; a1 0.31; b1 1; c1 10; F3 no; // Default LM coefficients ca1 2; ca2 0.06; ce1 1; ce2 50; cThetat 0.03; sigmaThetat 2; lambdaErr 1e-6; maxLambdaIter 10; }
Definition at line 101 of file kOmegaSSTLM.H.
typedef BasicMomentumTransportModel::alphaField alphaField |
Definition at line 202 of file kOmegaSSTLM.H.
typedef BasicMomentumTransportModel::rhoField rhoField |
Definition at line 203 of file kOmegaSSTLM.H.
typedef BasicMomentumTransportModel::transportModel transportModel |
Definition at line 204 of file kOmegaSSTLM.H.
kOmegaSSTLM | ( | 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 338 of file kOmegaSSTLM.C.
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::Fonset().
|
delete |
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 230 of file kOmegaSSTLM.H.
References kOmegaSSTLM< BasicMomentumTransportModel >::read().
|
protectedvirtual |
Modified form of the k-omega SST F1 function.
Reimplemented from kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >.
Definition at line 39 of file kOmegaSSTLM.C.
References Foam::exp(), F3, Foam::max(), kOmegaSSTLM< BasicMomentumTransportModel >::Pk(), Foam::pow(), Foam::Ry(), and Foam::sqrt().
|
protectedvirtual |
Modified form of the k-omega SST k production rate.
Definition at line 52 of file kOmegaSSTLM.C.
References kOmegaSSTLM< BasicMomentumTransportModel >::epsilonByk(), and kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >::Pk().
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::F1().
|
protectedvirtual |
Modified form of the k-omega SST epsilon/k.
Definition at line 63 of file kOmegaSSTLM.C.
References kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >::epsilonByk(), kOmegaSSTLM< BasicMomentumTransportModel >::Fthetat(), Foam::max(), and Foam::min().
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::Pk().
|
protected |
Freestream blending-function.
Definition at line 76 of file kOmegaSSTLM.C.
References delta, Foam::exp(), IOobject::groupName(), Foam::max(), Foam::min(), DimensionedField< Type, GeoMesh >::New(), Foam::pow4(), Foam::sqr(), and y.
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::epsilonByk().
|
protected |
Empirical correlation for critical Reynolds number where the.
intermittency first starts to increase in the boundary layer
Definition at line 107 of file kOmegaSSTLM.C.
References Foam::dimless, kOmegaSSTLM< BasicMomentumTransportModel >::Flength(), forAll, IOobject::groupName(), DimensionedField< Type, GeoMesh >::New(), Foam::pow3(), Foam::pow4(), and Foam::sqr().
|
protected |
Empirical correlation that controls the length of the.
transition region
Definition at line 143 of file kOmegaSSTLM.C.
References Foam::dimless, Foam::e, Foam::exp(), forAll, IOobject::groupName(), DimensionedField< Type, GeoMesh >::New(), Foam::pow3(), kOmegaSSTLM< BasicMomentumTransportModel >::ReThetat0(), Foam::sqr(), and y.
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::ReThetac().
|
protected |
Transition onset location control function.
Definition at line 311 of file kOmegaSSTLM.C.
References IOobject::groupName(), kOmegaSSTLM< BasicMomentumTransportModel >::kOmegaSSTLM(), Foam::max(), Foam::min(), DimensionedField< Type, GeoMesh >::New(), Foam::pow3(), and Foam::pow4().
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::ReThetat0().
|
protected |
Return the transition onset momentum-thickness Reynolds number.
(based on freestream conditions)
Definition at line 202 of file kOmegaSSTLM.C.
References Foam::dimless, Foam::endl(), Foam::exp(), kOmegaSSTLM< BasicMomentumTransportModel >::Fonset(), forAll, IOobject::groupName(), k, lambda(), Foam::mag(), Foam::max(), Foam::min(), DimensionedField< Type, GeoMesh >::New(), Foam::pow(), Foam::pow3(), Foam::sqr(), Foam::sqrt(), and WarningInFunction.
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::Flength().
|
protected |
Solve the turbulence equations and correct the turbulence viscosity.
Definition at line 488 of file kOmegaSSTLM.C.
References alpha(), Foam::bound(), tmp< T >::clear(), optionList::constrain(), optionList::correct(), Foam::fvm::ddt(), Foam::fvm::div(), Foam::exp(), fvOptions, Foam::fvc::grad(), k, Foam::fvm::laplacian(), Foam::mag(), Foam::magSqr(), Foam::max(), Foam::min(), options::New(), Foam::pow4(), tmp< T >::ref(), rho, Foam::skew(), Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), Foam::symm(), U, and y.
TypeName | ( | "kOmegaSSTLM< BasicMomentumTransportModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented from kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >.
Definition at line 465 of file kOmegaSSTLM.C.
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::~kOmegaSSTLM().
|
inline |
Access function transition onset momentum-thickness Reynolds number.
Definition at line 240 of file kOmegaSSTLM.H.
References kOmegaSSTLM< BasicMomentumTransportModel >::ReThetat_.
|
inline |
Access function to intermittency.
Definition at line 246 of file kOmegaSSTLM.H.
References kOmegaSSTLM< BasicMomentumTransportModel >::DReThetatEff(), kOmegaSSTLM< BasicMomentumTransportModel >::gammaInt_, GeometricField< scalar, fvPatchField, volMesh >::New(), and eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_.
|
inline |
Return the effective diffusivity for transition onset.
momentum-thickness Reynolds number
Definition at line 253 of file kOmegaSSTLM.H.
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::gammaInt().
|
inline |
Return the effective diffusivity for intermittency.
Definition at line 263 of file kOmegaSSTLM.H.
References kOmegaSSTLM< BasicMomentumTransportModel >::correct(), GeometricField< scalar, fvPatchField, volMesh >::New(), eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_, and kOmegaSSTLM< BasicMomentumTransportModel >::operator=().
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Reimplemented from kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >.
Definition at line 586 of file kOmegaSSTLM.C.
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::DgammaIntEff().
|
delete |
Disallow default bitwise assignment.
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::DgammaIntEff().
|
protected |
Definition at line 111 of file kOmegaSSTLM.H.
|
protected |
Definition at line 112 of file kOmegaSSTLM.H.
|
protected |
Definition at line 114 of file kOmegaSSTLM.H.
|
protected |
Definition at line 115 of file kOmegaSSTLM.H.
|
protected |
Definition at line 117 of file kOmegaSSTLM.H.
|
protected |
Definition at line 118 of file kOmegaSSTLM.H.
|
protected |
Convergence criterion for the lambda/thetat loop.
Definition at line 121 of file kOmegaSSTLM.H.
|
protected |
Maximum number of iterations to converge the lambda/thetat loop.
Definition at line 124 of file kOmegaSSTLM.H.
|
protected |
Stabilization for division by the magnitude of the velocity.
Definition at line 127 of file kOmegaSSTLM.H.
|
protected |
Transition onset momentum-thickness Reynolds number.
Definition at line 133 of file kOmegaSSTLM.H.
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::ReThetat().
|
protected |
Intermittency.
Definition at line 136 of file kOmegaSSTLM.H.
Referenced by kOmegaSSTLM< BasicMomentumTransportModel >::gammaInt().
|
protected |
Effective intermittency.
Definition at line 139 of file kOmegaSSTLM.H.