Skip to content
Snippets Groups Projects
Commit 276cef6a authored by Andrew Heather's avatar Andrew Heather
Browse files

ENH: kOmegaSST - added turbulence decay control terms

Based on the reference:
      Spalart, P. R. and Rumsey, C. L. (2007).
      Effective Inflow Conditions for Turbulence Models in Aerodynamic
      Calculations
      AIAA Journal, 45(10), 2544 - 2553.

The decay control default is off for backwards compatibility.  To enable
it, add the following to the coefficients dictionary

      // Optional decay control
      decayControl    yes;
      kInf            \<far-field k value\>;
      omegaInf        \<far-field omega value\>;
parent ef08bc56
No related merge requests found
......@@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2016 OpenCFD Ltd.
\\/ M anipulation | Copyright (C) 2016-2017 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
......@@ -384,15 +384,65 @@ kOmegaSSTBase<BasicEddyViscosityModel>::kOmegaSSTBase
IOobject::AUTO_WRITE
),
this->mesh_
),
decayControl_
(
Switch::lookupOrAddToDict
(
"decayControl",
this->coeffDict_,
false
)
),
kInf_
(
dimensioned<scalar>::lookupOrAddToDict
(
"kInf",
this->coeffDict_,
k_.dimensions(),
0
)
),
omegaInf_
(
dimensioned<scalar>::lookupOrAddToDict
(
"omegaInf",
this->coeffDict_,
omega_.dimensions(),
0
)
)
{
bound(k_, this->kMin_);
bound(omega_, this->omegaMin_);
setDecayControl(this->coeffDict_);
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class BasicEddyViscosityModel>
void kOmegaSSTBase<BasicEddyViscosityModel>::setDecayControl
(
const dictionary& dict
)
{
decayControl_.readIfPresent("decayControl", this->coeffDict());
if (decayControl_)
{
kInf_.read(this->coeffDict());
omegaInf_.read(this->coeffDict());
Info<< " Employing decay control with kInf:" << kInf_
<< " and omegaInf:" << omegaInf_ << endl;
}
}
template<class BasicEddyViscosityModel>
bool kOmegaSSTBase<BasicEddyViscosityModel>::read()
{
......@@ -412,6 +462,8 @@ bool kOmegaSSTBase<BasicEddyViscosityModel>::read()
c1_.readIfPresent(this->coeffDict());
F3_.readIfPresent("F3", this->coeffDict());
setDecayControl(this->coeffDict());
return true;
}
else
......@@ -476,6 +528,7 @@ void kOmegaSSTBase<BasicEddyViscosityModel>::correct()
alpha()*rho()*(F1() - scalar(1))*CDkOmega()/omega_(),
omega_
)
+ alpha()*rho()*beta*sqr(omegaInf_)
+ Qsas(S2(), gamma, beta)
+ omegaSource()
+ fvOptions(alpha, rho, omega_)
......@@ -499,6 +552,7 @@ void kOmegaSSTBase<BasicEddyViscosityModel>::correct()
alpha()*rho()*Pk(G)
- fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
- fvm::Sp(alpha()*rho()*epsilonByk(F1, tgradU()), k_)
+ alpha()*rho()*betaStar_*omegaInf_*kInf_
+ kSource()
+ fvOptions(alpha, rho, k_)
);
......
......@@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
\\/ M anipulation | Copyright (C) 2017 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
......@@ -52,13 +52,21 @@ Description
and the addition of the optional F3 term for rough walls from
\verbatim
Hellsten, A. (1998).
"Some Improvements in Menter's k-omega-SST turbulence model"
Some Improvements in Menter's k-omega-SST turbulence model
29th AIAA Fluid Dynamics Conference, AIAA-98-2554.
\endverbatim
and the optional decay control from:
\verbatim
Spalart, P. R. and Rumsey, C. L. (2007).
Effective Inflow Conditions for Turbulence Models in Aerodynamic
Calculations
AIAA Journal, 45(10), 2544 - 2553.
\endverbatim
Note that this implementation is written in terms of alpha diffusion
coefficients rather than the more traditional sigma (alpha = 1/sigma) so
that the blending can be applied to all coefficuients in a consistent
that the blending can be applied to all coefficients in a consistent
manner. The paper suggests that sigma is blended but this would not be
consistent with the blending of the k-epsilon and k-omega models.
......@@ -76,19 +84,24 @@ Description
\verbatim
kOmegaSSTBaseCoeffs
{
alphaK1 0.85;
alphaK2 1.0;
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.0;
c1 10.0;
F3 no;
alphaK1 0.85;
alphaK2 1.0;
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.0;
c1 10.0;
F3 no;
// Optional decay control
decayControl yes;
kInf \<far-field k value\>;
omegaInf \<far-field omega value\>;
}
\endverbatim
......@@ -145,6 +158,7 @@ protected:
dimensionedScalar b1_;
dimensionedScalar c1_;
//- Flag to include the F3 term
Switch F3_;
......@@ -159,8 +173,18 @@ protected:
volScalarField omega_;
// Decay control
//- Flag to include the decay control
Switch decayControl_;
dimensionedScalar omegaInf_;
dimensionedScalar kInf_;
// Protected Member Functions
void setDecayControl(const dictionary& dict);
virtual tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
virtual tmp<volScalarField> F2() const;
virtual tmp<volScalarField> F3() const;
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment