Skip to content
Snippets Groups Projects
Commit dfddcad5 authored by Henry's avatar Henry
Browse files

TurbulenceModels: Do not correct nut if the turbulence model is switched-off

parent b0eb5b8a
No related branches found
No related tags found
No related merge requests found
......@@ -166,6 +166,11 @@ tmp<fvScalarMatrix> kEqn<BasicTurbulenceModel>::kSource() const
template<class BasicTurbulenceModel>
void kEqn<BasicTurbulenceModel>::correct()
{
if (!this->turbulence_)
{
return;
}
// Local references
const alphaField& alpha = this->alpha_;
const rhoField& rho = this->rho_;
......@@ -174,12 +179,6 @@ void kEqn<BasicTurbulenceModel>::correct()
const volVectorField& U = this->U_;
volScalarField& nut = this->nut_;
if (!this->turbulence_)
{
correctNut();
return;
}
eddyViscosity<LESModel<BasicTurbulenceModel> >::correct();
volScalarField divU(fvc::div(fvc::absolute(phi/fvc::interpolate(rho), U)));
......
......@@ -215,6 +215,11 @@ tmp<fvScalarMatrix> kEpsilon<BasicTurbulenceModel>::epsilonSource() const
template<class BasicTurbulenceModel>
void kEpsilon<BasicTurbulenceModel>::correct()
{
if (!this->turbulence_)
{
return;
}
// Local references
const alphaField& alpha = this->alpha_;
const rhoField& rho = this->rho_;
......@@ -223,12 +228,6 @@ void kEpsilon<BasicTurbulenceModel>::correct()
const volVectorField& U = this->U_;
volScalarField& nut = this->nut_;
if (!this->turbulence_)
{
correctNut();
return;
}
eddyViscosity<RASModel<BasicTurbulenceModel> >::correct();
volScalarField divU(fvc::div(fvc::absolute(phi/fvc::interpolate(rho), U)));
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment