mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
omegaWallFunction: Improved low-Reynolds number behavior and consistency with the epsilonWallFunction
Changed default mode of operation to use standard y+ based switching
rather than the previous ad hoc blending and added consistent handling
of the near-wall generation term.
This boundary condition provides a wall constraint on turbulnce specific
dissipation, omega for both low and high Reynolds number turbulence models.
The near-wall omega may be either blended between the viscous region and
logarithmic region values using:
\f[
\omega = sqrt(\omega_{vis}^2 + \omega_{log}^2)
\f]
where
\vartable
\omega_{vis} | omega in viscous region
\omega_{log} | omega in logarithmic region
\endvartable
see eq.(15) of:
\verbatim
Menter, F., Esch, T.
"Elements of Industrial Heat Transfer Prediction"
16th Brazilian Congress of Mechanical Engineering (COBEM),
Nov. 2001
\endverbatim
or switched between these values based on the laminar-to-turbulent y+ value
derived from kappa and E. Recent tests have shown that the standard
switching method provides more accurate results for 10 < y+ < 30 when used
with high Reynolds number wall-functions and both methods provide accurate
results when used with continuous wall-functions. Based on this the
standard switching method is used by default.
This commit is contained in:
@ -24,12 +24,12 @@ License
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "omegaWallFunctionFvPatchScalarField.H"
|
||||
#include "nutWallFunctionFvPatchScalarField.H"
|
||||
#include "turbulenceModel.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "fvMatrix.H"
|
||||
#include "volFields.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "nutkWallFunctionFvPatchScalarField.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
@ -63,6 +63,7 @@ void omegaWallFunctionFvPatchScalarField::writeLocalEntries(Ostream& os) const
|
||||
os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
|
||||
os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
|
||||
os.writeKeyword("beta1") << beta1_ << token::END_STATEMENT << nl;
|
||||
os.writeKeyword("blended") << blended_ << token::END_STATEMENT << nl;
|
||||
}
|
||||
|
||||
|
||||
@ -209,8 +210,8 @@ void omegaWallFunctionFvPatchScalarField::calculate
|
||||
const turbulenceModel& turbModel,
|
||||
const List<scalar>& cornerWeights,
|
||||
const fvPatch& patch,
|
||||
scalarField& G,
|
||||
scalarField& omega
|
||||
scalarField& G0,
|
||||
scalarField& omega0
|
||||
)
|
||||
{
|
||||
const label patchi = patch.index();
|
||||
@ -232,25 +233,56 @@ void omegaWallFunctionFvPatchScalarField::calculate
|
||||
|
||||
const scalarField magGradUw(mag(Uw.snGrad()));
|
||||
|
||||
typedef DimensionedField<scalar, volMesh> FieldType;
|
||||
const FieldType& G = db().lookupObject<FieldType>(turbModel.GName());
|
||||
|
||||
// Set omega and G
|
||||
forAll(nutw, facei)
|
||||
{
|
||||
label celli = patch.faceCells()[facei];
|
||||
const label celli = patch.faceCells()[facei];
|
||||
|
||||
scalar w = cornerWeights[facei];
|
||||
const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
|
||||
|
||||
scalar omegaVis = 6.0*nuw[facei]/(beta1_*sqr(y[facei]));
|
||||
const scalar w = cornerWeights[facei];
|
||||
|
||||
scalar omegaLog = sqrt(k[celli])/(Cmu25*kappa_*y[facei]);
|
||||
const scalar omegaVis = 6*nuw[facei]/(beta1_*sqr(y[facei]));
|
||||
const scalar omegaLog = sqrt(k[celli])/(Cmu25*kappa_*y[facei]);
|
||||
|
||||
omega[celli] += w*sqrt(sqr(omegaVis) + sqr(omegaLog));
|
||||
// Switching between the laminar sub-layer and the log-region rather
|
||||
// than blending has been found to provide more accurate results over a
|
||||
// range of near-wall y+.
|
||||
//
|
||||
// For backward-compatibility the blending method is provided as an
|
||||
// option
|
||||
|
||||
G[celli] +=
|
||||
w
|
||||
*(nutw[facei] + nuw[facei])
|
||||
*magGradUw[facei]
|
||||
*Cmu25*sqrt(k[celli])
|
||||
/(kappa_*y[facei]);
|
||||
if (blended_)
|
||||
{
|
||||
omega0[celli] += w*sqrt(sqr(omegaVis) + sqr(omegaLog));
|
||||
}
|
||||
|
||||
if (yPlus > yPlusLam_)
|
||||
{
|
||||
if (!blended_)
|
||||
{
|
||||
omega0[celli] += w*omegaLog;
|
||||
}
|
||||
|
||||
G0[celli] +=
|
||||
w
|
||||
*(nutw[facei] + nuw[facei])
|
||||
*magGradUw[facei]
|
||||
*Cmu25*sqrt(k[celli])
|
||||
/(kappa_*y[facei]);
|
||||
}
|
||||
else
|
||||
{
|
||||
if (!blended_)
|
||||
{
|
||||
omega0[celli] += w*omegaVis;
|
||||
}
|
||||
|
||||
G0[celli] += w*G[celli];
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
@ -268,7 +300,8 @@ omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
kappa_(0.41),
|
||||
E_(9.8),
|
||||
beta1_(0.075),
|
||||
yPlusLam_(nutkWallFunctionFvPatchScalarField::yPlusLam(kappa_, E_)),
|
||||
blended_(false),
|
||||
yPlusLam_(nutWallFunctionFvPatchScalarField::yPlusLam(kappa_, E_)),
|
||||
G_(),
|
||||
omega_(),
|
||||
initialised_(false),
|
||||
@ -292,6 +325,7 @@ omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
kappa_(ptf.kappa_),
|
||||
E_(ptf.E_),
|
||||
beta1_(ptf.beta1_),
|
||||
blended_(ptf.blended_),
|
||||
yPlusLam_(ptf.yPlusLam_),
|
||||
G_(),
|
||||
omega_(),
|
||||
@ -315,7 +349,8 @@ omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
|
||||
E_(dict.lookupOrDefault<scalar>("E", 9.8)),
|
||||
beta1_(dict.lookupOrDefault<scalar>("beta1", 0.075)),
|
||||
yPlusLam_(nutkWallFunctionFvPatchScalarField::yPlusLam(kappa_, E_)),
|
||||
blended_(dict.lookupOrDefault<Switch>("blended", false)),
|
||||
yPlusLam_(nutWallFunctionFvPatchScalarField::yPlusLam(kappa_, E_)),
|
||||
G_(),
|
||||
omega_(),
|
||||
initialised_(false),
|
||||
@ -339,6 +374,7 @@ omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
kappa_(owfpsf.kappa_),
|
||||
E_(owfpsf.E_),
|
||||
beta1_(owfpsf.beta1_),
|
||||
blended_(owfpsf.blended_),
|
||||
yPlusLam_(owfpsf.yPlusLam_),
|
||||
G_(),
|
||||
omega_(),
|
||||
@ -361,6 +397,7 @@ omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
kappa_(owfpsf.kappa_),
|
||||
E_(owfpsf.E_),
|
||||
beta1_(owfpsf.beta1_),
|
||||
blended_(owfpsf.blended_),
|
||||
yPlusLam_(owfpsf.yPlusLam_),
|
||||
G_(),
|
||||
omega_(),
|
||||
|
||||
@ -28,8 +28,11 @@ Group
|
||||
grpWallFunctions
|
||||
|
||||
Description
|
||||
This boundary condition provides a wall function constraint on turbulnce
|
||||
specific dissipation, omega. The values are computed using:
|
||||
This boundary condition provides a wall constraint on turbulnce specific
|
||||
dissipation, omega for both low and high Reynolds number turbulence models.
|
||||
|
||||
The near-wall omega may be either blended between the viscous region and
|
||||
logarithmic region values using:
|
||||
|
||||
\f[
|
||||
\omega = sqrt(\omega_{vis}^2 + \omega_{log}^2)
|
||||
@ -42,7 +45,7 @@ Description
|
||||
\omega_{log} | omega in logarithmic region
|
||||
\endvartable
|
||||
|
||||
Model described by Eq.(15) of:
|
||||
see eq.(15) of:
|
||||
\verbatim
|
||||
Menter, F., Esch, T.
|
||||
"Elements of Industrial Heat Transfer Prediction"
|
||||
@ -50,13 +53,21 @@ Description
|
||||
Nov. 2001
|
||||
\endverbatim
|
||||
|
||||
or switched between these values based on the laminar-to-turbulent y+ value
|
||||
derived from kappa and E. Recent tests have shown that the standard
|
||||
switching method provides more accurate results for 10 < y+ < 30 when used
|
||||
with high Reynolds number wall-functions and both methods provide accurate
|
||||
results when used with continuous wall-functions. Based on this the
|
||||
standard switching method is used by default.
|
||||
|
||||
Usage
|
||||
\table
|
||||
Property | Description | Required | Default value
|
||||
Cmu | model coefficient | no | 0.09
|
||||
kappa | Von Karman constant | no | 0.41
|
||||
E | model coefficient | no | 9.8
|
||||
beta1 | model coefficient | no | 0.075
|
||||
Cmu | Model coefficient | no | 0.09
|
||||
kappa | von Karman constant | no | 0.41
|
||||
E | Model coefficient | no | 9.8
|
||||
beta1 | Model coefficient | no | 0.075
|
||||
blended | Blending switch | no | false
|
||||
\endtable
|
||||
|
||||
Example of the boundary condition specification:
|
||||
@ -67,6 +78,10 @@ Usage
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
See also
|
||||
Foam::fixedInternalValueFvPatchField
|
||||
Foam::epsilonWallFunctionFvPatchScalarField
|
||||
|
||||
SourceFiles
|
||||
omegaWallFunctionFvPatchScalarField.C
|
||||
|
||||
@ -111,7 +126,10 @@ protected:
|
||||
//- beta1 coefficient
|
||||
scalar beta1_;
|
||||
|
||||
//- Y+ at the edge of the laminar sublayer
|
||||
//- beta1 coefficient
|
||||
Switch blended_;
|
||||
|
||||
//- y+ at the edge of the laminar sublayer
|
||||
scalar yPlusLam_;
|
||||
|
||||
//- Local copy of turbulence G field
|
||||
|
||||
Reference in New Issue
Block a user