diff --git a/src/TurbulenceModels/turbulenceModels/Base/kOmegaSST/kOmegaSSTBase.C b/src/TurbulenceModels/turbulenceModels/Base/kOmegaSST/kOmegaSSTBase.C index 226ac013c7..1c3febac0e 100644 --- a/src/TurbulenceModels/turbulenceModels/Base/kOmegaSST/kOmegaSSTBase.C +++ b/src/TurbulenceModels/turbulenceModels/Base/kOmegaSST/kOmegaSSTBase.C @@ -132,6 +132,17 @@ void kOmegaSSTBase::correctNut() } +template +Foam::tmp kOmegaSSTBase::S2 +( + const volScalarField& F1, + const volTensorField& gradU +) const +{ + return 2*magSqr(symm(gradU)); +} + + template tmp kOmegaSSTBase::Pk ( @@ -143,11 +154,10 @@ tmp kOmegaSSTBase::Pk template -tmp -kOmegaSSTBase::epsilonByk +tmp kOmegaSSTBase::epsilonByk ( - const volScalarField& F1, - const volTensorField& gradU + const volScalarField& /* F1 not used */, + const volTensorField& /* gradU not used */ ) const { return betaStar_*omega_(); @@ -165,8 +175,7 @@ tmp kOmegaSSTBase::GbyNu return min ( GbyNu0, - (c1_/a1_)*betaStar_*omega_() - *max(a1_*omega_(), b1_*F2*sqrt(S2)) + (c1_/a1_)*betaStar_*omega_()*max(a1_*omega_(), b1_*F2*sqrt(S2)) ); } @@ -174,13 +183,10 @@ tmp kOmegaSSTBase::GbyNu template tmp kOmegaSSTBase::kSource() const { - return tmp + return tmp::New ( - new fvScalarMatrix - ( - k_, - dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime - ) + k_, + dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime ); } @@ -188,13 +194,10 @@ tmp kOmegaSSTBase::kSource() const template tmp kOmegaSSTBase::omegaSource() const { - return tmp + return tmp::New ( - new fvScalarMatrix - ( - omega_, - dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime - ) + omega_, + dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime ); } @@ -207,13 +210,10 @@ tmp kOmegaSSTBase::Qsas const volScalarField::Internal& beta ) const { - return tmp + return tmp::New ( - new fvScalarMatrix - ( - omega_, - dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime - ) + omega_, + dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime ); } @@ -500,18 +500,6 @@ void kOmegaSSTBase::correct() volScalarField::Internal divU(fvc::div(fvc::absolute(this->phi(), U))); - tmp tgradU = fvc::grad(U); - volScalarField S2(2*magSqr(symm(tgradU()))); - volScalarField::Internal GbyNu0 - ( - this->type() + ":GbyNu", - (tgradU() && dev(twoSymm(tgradU()))) - ); - volScalarField::Internal G(this->GName(), nut*GbyNu0); - - // Update omega and G at the wall - omega_.boundaryFieldRef().updateCoeffs(); - volScalarField CDkOmega ( (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_ @@ -520,6 +508,19 @@ void kOmegaSSTBase::correct() volScalarField F1(this->F1(CDkOmega)); volScalarField F23(this->F23()); + tmp tgradU = fvc::grad(U); + volScalarField S2(this->S2(F1, tgradU())); + volScalarField::Internal GbyNu0 + ( + this->type() + ":GbyNu", + (tgradU() && dev(twoSymm(tgradU()))) + ); + volScalarField::Internal G(this->GName(), nut*GbyNu0); + tgradU.clear(); + + // Update omega and G at the wall + omega_.boundaryFieldRef().updateCoeffs(); + { volScalarField::Internal gamma(this->gamma(F1)); volScalarField::Internal beta(this->beta(F1)); diff --git a/src/TurbulenceModels/turbulenceModels/Base/kOmegaSST/kOmegaSSTBase.H b/src/TurbulenceModels/turbulenceModels/Base/kOmegaSST/kOmegaSSTBase.H index ef2528fdcb..cfdf33af82 100644 --- a/src/TurbulenceModels/turbulenceModels/Base/kOmegaSST/kOmegaSSTBase.H +++ b/src/TurbulenceModels/turbulenceModels/Base/kOmegaSST/kOmegaSSTBase.H @@ -253,6 +253,13 @@ protected: virtual void correctNut(); + //- Return square of strain rate + virtual tmp S2 + ( + const volScalarField& F1, + const volTensorField& gradU + ) const; + //- Return k production rate virtual tmp Pk ( @@ -321,22 +328,20 @@ public: //- Return the effective diffusivity for k tmp DkEff(const volScalarField& F1) const { - return tmp + return tmp::New ( - new volScalarField("DkEff", alphaK(F1)*this->nut_ + this->nu()) + "DkEff", + alphaK(F1)*this->nut_ + this->nu() ); } //- Return the effective diffusivity for omega tmp DomegaEff(const volScalarField& F1) const { - return tmp + return tmp::New ( - new volScalarField - ( - "DomegaEff", - alphaOmega(F1)*this->nut_ + this->nu() - ) + "DomegaEff", + alphaOmega(F1)*this->nut_ + this->nu() ); } diff --git a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDDES/kOmegaSSTDDES.C b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDDES/kOmegaSSTDDES.C index 136421423e..7ab2b4d474 100644 --- a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDDES/kOmegaSSTDDES.C +++ b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDDES/kOmegaSSTDDES.C @@ -37,41 +37,13 @@ namespace LESModels // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * // -template -tmp kOmegaSSTDDES::rd -( - const volScalarField& magGradU -) const -{ - tmp tr - ( - min - ( - this->nuEff() - /( - max - ( - magGradU, - dimensionedScalar("SMALL", magGradU.dimensions(), SMALL) - ) - *sqr(this->kappa_*this->y_) - ), - scalar(10) - ) - ); - tr.ref().boundaryFieldRef() == 0.0; - - return tr; -} - - template tmp kOmegaSSTDDES::fd ( const volScalarField& magGradU ) const { - return 1 - tanh(pow(Cd1_*rd(magGradU), Cd2_)); + return 1 - tanh(pow(Cd1_*this->r(this->nuEff(), magGradU), Cd2_)); } diff --git a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDDES/kOmegaSSTDDES.H b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDDES/kOmegaSSTDDES.H index cd85bb5259..c84a0f49df 100644 --- a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDDES/kOmegaSSTDDES.H +++ b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDDES/kOmegaSSTDDES.H @@ -71,8 +71,6 @@ class kOmegaSSTDDES tmp fd(const volScalarField& magGradU) const; - tmp rd(const volScalarField& magGradU) const; - //- No copy construct kOmegaSSTDDES(const kOmegaSSTDDES&) = delete; diff --git a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDES/kOmegaSSTDES.C b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDES/kOmegaSSTDES.C index 08f7b76bf9..0d4939975a 100644 --- a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDES/kOmegaSSTDES.C +++ b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDES/kOmegaSSTDES.C @@ -55,6 +55,25 @@ void kOmegaSSTDES::correctNut() } +template +tmp kOmegaSSTDES::r +( + const volScalarField& nur, + const volScalarField& magGradU +) const +{ + const dimensionedScalar eps("SMALL", magGradU.dimensions(), SMALL); + + tmp tr = + min(nur/(max(magGradU, eps)*sqr(this->kappa_*this->y_)), scalar(10)); + + tr.ref().boundaryFieldRef() == 0; + + return tr; + +} + + template tmp kOmegaSSTDES::dTilda ( @@ -187,31 +206,19 @@ tmp kOmegaSSTDES::LESRegion() const const volScalarField F1(this->F1(CDkOmega)); - tmp tLESRegion + return tmp::New ( - new volScalarField + IOobject::scopedName("DES", "LESRegion"), + neg ( - IOobject + dTilda ( - "DES::LESRegion", - this->mesh_.time().timeName(), - this->mesh_, - IOobject::NO_READ, - IOobject::NO_WRITE - ), - neg - ( - dTilda - ( - mag(fvc::grad(U)), - F1*CDESkom_ + (1 - F1)*CDESkeps_ - ) - - sqrt(k)/(this->betaStar_*omega) + mag(fvc::grad(U)), + F1*CDESkom_ + (1 - F1)*CDESkeps_ ) + - sqrt(k)/(this->betaStar_*omega) ) ); - - return tLESRegion; } diff --git a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDES/kOmegaSSTDES.H b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDES/kOmegaSSTDES.H index 40e4428afa..08eecbf236 100644 --- a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDES/kOmegaSSTDES.H +++ b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTDES/kOmegaSSTDES.H @@ -104,6 +104,12 @@ protected: virtual void correctNut(const volScalarField& S2); virtual void correctNut(); + tmp r + ( + const volScalarField& nur, + const volScalarField& magGradU + ) const; + //- Length scale virtual tmp dTilda ( diff --git a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTIDDES/kOmegaSSTIDDES.C b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTIDDES/kOmegaSSTIDDES.C index 4714eed05c..268935f45e 100644 --- a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTIDDES/kOmegaSSTIDDES.C +++ b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTIDDES/kOmegaSSTIDDES.C @@ -64,7 +64,7 @@ tmp kOmegaSSTIDDES::ft const volScalarField& magGradU ) const { - return tanh(pow3(sqr(Ct_)*rd(this->nut_, magGradU))); + return tanh(pow3(sqr(Ct_)*this->r(this->nut_, magGradU))); } @@ -74,36 +74,7 @@ tmp kOmegaSSTIDDES::fl const volScalarField& magGradU ) const { - return tanh(pow(sqr(Cl_)*rd(this->nu(), magGradU), 10)); -} - - -template -tmp kOmegaSSTIDDES::rd -( - const volScalarField& nur, - const volScalarField& magGradU -) const -{ - tmp tr - ( - min - ( - nur - /( - max - ( - magGradU, - dimensionedScalar("SMALL", magGradU.dimensions(), SMALL) - ) - *sqr(this->kappa_*this->y_) - ), - scalar(10) - ) - ); - tr.ref().boundaryFieldRef() == 0.0; - - return tr; + return tanh(pow(sqr(Cl_)*this->r(this->nu(), magGradU), 10)); } @@ -113,7 +84,7 @@ tmp kOmegaSSTIDDES::fdt const volScalarField& magGradU ) const { - return 1 - tanh(pow(Cdt1_*rd(this->nut_, magGradU), Cdt2_)); + return 1 - tanh(pow(Cdt1_*this->r(this->nut_, magGradU), Cdt2_)); } diff --git a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTIDDES/kOmegaSSTIDDES.H b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTIDDES/kOmegaSSTIDDES.H index 32ccb217b9..2da9cd3485 100644 --- a/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTIDDES/kOmegaSSTIDDES.H +++ b/src/TurbulenceModels/turbulenceModels/DES/kOmegaSSTIDDES/kOmegaSSTIDDES.H @@ -77,12 +77,6 @@ class kOmegaSSTIDDES tmp ft(const volScalarField& magGradU) const; tmp fl(const volScalarField& magGradU) const; - tmp rd - ( - const volScalarField& nur, - const volScalarField& magGradU - ) const; - //- Delay function tmp fdt(const volScalarField& magGradU) const; @@ -104,9 +98,10 @@ protected: dimensionedScalar Cl_; dimensionedScalar Ct_; - // Fields - const IDDESDelta& IDDESDelta_; + //- IDDES delta + const IDDESDelta& IDDESDelta_; + // Protected Member Functions