From 48fd678e3d8e6eec833bf8bc732027ca056761c6 Mon Sep 17 00:00:00 2001 From: Sergio Ferraris Date: Wed, 24 Nov 2010 12:01:40 +0000 Subject: [PATCH] ENH: change S definition and add extra reference --- .../incompressible/LES/kOmegaSSTSAS/kOmegaSSTSAS.C | 8 ++++---- .../incompressible/LES/kOmegaSSTSAS/kOmegaSSTSAS.H | 7 ++++++- 2 files changed, 10 insertions(+), 5 deletions(-) diff --git a/src/turbulenceModels/incompressible/LES/kOmegaSSTSAS/kOmegaSSTSAS.C b/src/turbulenceModels/incompressible/LES/kOmegaSSTSAS/kOmegaSSTSAS.C index df57ff44..4541bdca 100644 --- a/src/turbulenceModels/incompressible/LES/kOmegaSSTSAS/kOmegaSSTSAS.C +++ b/src/turbulenceModels/incompressible/LES/kOmegaSSTSAS/kOmegaSSTSAS.C @@ -323,7 +323,7 @@ kOmegaSSTSAS::kOmegaSSTSAS mesh_ ) { - updateSubGridScaleFields(magSqr(symm(fvc::grad(U)))); + updateSubGridScaleFields(magSqr(2.0*symm(fvc::grad(U)))); printCoeffs(); } @@ -340,7 +340,7 @@ void kOmegaSSTSAS::correct(const tmp& gradU) y_.correct(); } - volScalarField S2 = magSqr(symm(gradU())); + volScalarField S2 = magSqr(2.0*symm(gradU())); gradU.clear(); volVectorField gradK = fvc::grad(k_); @@ -349,7 +349,7 @@ void kOmegaSSTSAS::correct(const tmp& gradU) volScalarField CDkOmega = (2.0*alphaOmega2_)*(gradK & gradOmega)/(omega_ + omegaSmall_); volScalarField F1 = this->F1(CDkOmega); - volScalarField G = nuSgs_*2.0*S2; + volScalarField G = nuSgs_*0.5*S2; // Turbulent kinetic energy equation { @@ -386,7 +386,7 @@ void kOmegaSSTSAS::correct(const tmp& gradU) - fvm::Sp(fvc::div(phi()), omega_) - fvm::laplacian(DomegaEff(F1), omega_) == - gamma(F1)*2.0*S2 + gamma(F1)*0.5*S2 - fvm::Sp(beta(F1)*omega_, omega_) - fvm::SuSp // cross diffusion term ( diff --git a/src/turbulenceModels/incompressible/LES/kOmegaSSTSAS/kOmegaSSTSAS.H b/src/turbulenceModels/incompressible/LES/kOmegaSSTSAS/kOmegaSSTSAS.H index daee30c6..d2bbdb3d 100644 --- a/src/turbulenceModels/incompressible/LES/kOmegaSSTSAS/kOmegaSSTSAS.H +++ b/src/turbulenceModels/incompressible/LES/kOmegaSSTSAS/kOmegaSSTSAS.H @@ -27,7 +27,12 @@ Class Description kOmegaSSTSAS LES turbulence model for incompressible flows - Reference: + References: + + A Scale-Adaptive Simulation Model using Two-Equation Models + AIAA 2005-1095 + F. R. Menter and Y. Egorov + DESider A European Effort on Hybrid RANS-LES Modelling: Results of the European-Union Funded Project, 2004 - 2007 (Notes on Numerical Fluid Mechanics and Multidisciplinary Design).