mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
RAS/kOmegaSST: Added F3 coefficient for rough-walls from
Hellsten, A.
"Some Improvements in Menter’s k-omega-SST turbulence model"
29th AIAA Fluid Dynamics Conference,
AIAA-98-2554,
June 1998.
This commit is contained in:
@ -85,6 +85,18 @@ tmp<volScalarField> kOmegaSST::F2() const
|
||||
}
|
||||
|
||||
|
||||
tmp<volScalarField> kOmegaSST::F3() const
|
||||
{
|
||||
tmp<volScalarField> arg3 = min
|
||||
(
|
||||
150*(mu()/rho_)/(omega_*sqr(y_)),
|
||||
scalar(10)
|
||||
);
|
||||
|
||||
return 1 - tanh(pow4(arg3));
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
kOmegaSST::kOmegaSST
|
||||
@ -198,6 +210,15 @@ kOmegaSST::kOmegaSST
|
||||
0.31
|
||||
)
|
||||
),
|
||||
b1_
|
||||
(
|
||||
dimensioned<scalar>::lookupOrAddToDict
|
||||
(
|
||||
"b1",
|
||||
coeffDict_,
|
||||
1.0
|
||||
)
|
||||
),
|
||||
c1_
|
||||
(
|
||||
dimensioned<scalar>::lookupOrAddToDict
|
||||
@ -268,7 +289,7 @@ kOmegaSST::kOmegaSST
|
||||
/ max
|
||||
(
|
||||
a1_*omega_,
|
||||
F2()*sqrt(2.0)*mag(symm(fvc::grad(U_)))
|
||||
b1_*F2()*F3()*sqrt(2.0)*mag(symm(fvc::grad(U_)))
|
||||
)
|
||||
);
|
||||
mut_.correctBoundaryConditions();
|
||||
@ -347,6 +368,7 @@ bool kOmegaSST::read()
|
||||
beta2_.readIfPresent(coeffDict());
|
||||
betaStar_.readIfPresent(coeffDict());
|
||||
a1_.readIfPresent(coeffDict());
|
||||
b1_.readIfPresent(coeffDict());
|
||||
c1_.readIfPresent(coeffDict());
|
||||
|
||||
return true;
|
||||
@ -448,7 +470,7 @@ void kOmegaSST::correct()
|
||||
|
||||
|
||||
// Re-calculate viscosity
|
||||
mut_ = a1_*rho_*k_/max(a1_*omega_, F2()*sqrt(S2));
|
||||
mut_ = a1_*rho_*k_/max(a1_*omega_, b1_*F2()*F3()*sqrt(S2));
|
||||
mut_.correctBoundaryConditions();
|
||||
|
||||
// Re-calculate thermal diffusivity
|
||||
|
||||
@ -38,6 +38,15 @@ Description
|
||||
Nov. 2001
|
||||
\endverbatim
|
||||
|
||||
with the addition of the F3 term for rough walls from
|
||||
\verbatim
|
||||
Hellsten, A.
|
||||
"Some Improvements in Menter’s k-omega-SST turbulence model"
|
||||
29th AIAA Fluid Dynamics Conference,
|
||||
AIAA-98-2554,
|
||||
June 1998.
|
||||
\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
|
||||
@ -69,6 +78,7 @@ Description
|
||||
gamma1 0.5532;
|
||||
gamma2 0.4403;
|
||||
a1 0.31;
|
||||
b1 1.0;
|
||||
c1 10.0;
|
||||
}
|
||||
\endverbatim
|
||||
@ -125,6 +135,7 @@ protected:
|
||||
dimensionedScalar betaStar_;
|
||||
|
||||
dimensionedScalar a1_;
|
||||
dimensionedScalar b1_;
|
||||
dimensionedScalar c1_;
|
||||
|
||||
|
||||
@ -144,6 +155,7 @@ protected:
|
||||
|
||||
tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
|
||||
tmp<volScalarField> F2() const;
|
||||
tmp<volScalarField> F3() const;
|
||||
|
||||
tmp<volScalarField> blend
|
||||
(
|
||||
|
||||
@ -69,6 +69,7 @@ tmp<volScalarField> kOmegaSST::F1(const volScalarField& CDkOmega) const
|
||||
return tanh(pow4(arg1));
|
||||
}
|
||||
|
||||
|
||||
tmp<volScalarField> kOmegaSST::F2() const
|
||||
{
|
||||
tmp<volScalarField> arg2 = min
|
||||
@ -85,6 +86,18 @@ tmp<volScalarField> kOmegaSST::F2() const
|
||||
}
|
||||
|
||||
|
||||
tmp<volScalarField> kOmegaSST::F3() const
|
||||
{
|
||||
tmp<volScalarField> arg3 = min
|
||||
(
|
||||
150*nu()/(omega_*sqr(y_)),
|
||||
scalar(10)
|
||||
);
|
||||
|
||||
return 1 - tanh(pow4(arg3));
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
kOmegaSST::kOmegaSST
|
||||
@ -188,6 +201,15 @@ kOmegaSST::kOmegaSST
|
||||
0.31
|
||||
)
|
||||
),
|
||||
b1_
|
||||
(
|
||||
dimensioned<scalar>::lookupOrAddToDict
|
||||
(
|
||||
"b1",
|
||||
coeffDict_,
|
||||
1.0
|
||||
)
|
||||
),
|
||||
c1_
|
||||
(
|
||||
dimensioned<scalar>::lookupOrAddToDict
|
||||
@ -246,7 +268,7 @@ kOmegaSST::kOmegaSST
|
||||
/ max
|
||||
(
|
||||
a1_*omega_,
|
||||
F2()*sqrt(2.0)*mag(symm(fvc::grad(U_)))
|
||||
b1_*F2()*F3()*sqrt(2.0)*mag(symm(fvc::grad(U_)))
|
||||
)
|
||||
);
|
||||
nut_.correctBoundaryConditions();
|
||||
@ -338,6 +360,7 @@ bool kOmegaSST::read()
|
||||
beta2_.readIfPresent(coeffDict());
|
||||
betaStar_.readIfPresent(coeffDict());
|
||||
a1_.readIfPresent(coeffDict());
|
||||
b1_.readIfPresent(coeffDict());
|
||||
c1_.readIfPresent(coeffDict());
|
||||
|
||||
return true;
|
||||
@ -416,7 +439,7 @@ void kOmegaSST::correct()
|
||||
|
||||
|
||||
// Re-calculate viscosity
|
||||
nut_ = a1_*k_/max(a1_*omega_, F2()*sqrt(S2));
|
||||
nut_ = a1_*k_/max(a1_*omega_, b1_*F2()*F3()*sqrt(S2));
|
||||
nut_.correctBoundaryConditions();
|
||||
}
|
||||
|
||||
|
||||
@ -33,10 +33,19 @@ Description
|
||||
|
||||
Turbulence model described in:
|
||||
\verbatim
|
||||
Menter, F., Esch, T.
|
||||
"Elements of Industrial Heat Transfer Prediction"
|
||||
Menter, F., Esch, T.,
|
||||
"Elements of Industrial Heat Transfer Prediction",
|
||||
16th Brazilian Congress of Mechanical Engineering (COBEM),
|
||||
Nov. 2001
|
||||
Nov. 2001.
|
||||
\endverbatim
|
||||
|
||||
with the addition of the F3 term for rough walls from
|
||||
\verbatim
|
||||
Hellsten, A.
|
||||
"Some Improvements in Menter’s k-omega-SST turbulence model"
|
||||
29th AIAA Fluid Dynamics Conference,
|
||||
AIAA-98-2554,
|
||||
June 1998.
|
||||
\endverbatim
|
||||
|
||||
Note that this implementation is written in terms of alpha diffusion
|
||||
@ -69,6 +78,7 @@ Description
|
||||
gamma1 0.5532;
|
||||
gamma2 0.4403;
|
||||
a1 0.31;
|
||||
b1 1.0;
|
||||
c1 10.0;
|
||||
}
|
||||
\endverbatim
|
||||
@ -122,6 +132,7 @@ protected:
|
||||
dimensionedScalar betaStar_;
|
||||
|
||||
dimensionedScalar a1_;
|
||||
dimensionedScalar b1_;
|
||||
dimensionedScalar c1_;
|
||||
|
||||
//- Wall distance field
|
||||
@ -139,6 +150,7 @@ protected:
|
||||
|
||||
tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
|
||||
tmp<volScalarField> F2() const;
|
||||
tmp<volScalarField> F3() const;
|
||||
|
||||
tmp<volScalarField> blend
|
||||
(
|
||||
|
||||
Reference in New Issue
Block a user