This boundary condition provides a wall constraint on the specific dissipation rate, i.e. omega
, and the turbulent kinetic energy production contribution, i.e. G
, for low- and high-Reynolds number turbulence models.
More...
Public Member Functions | |
TypeName ("omegaWallFunction") | |
omegaWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
omegaWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &) | |
virtual tmp< fvPatchScalarField > | clone () const |
omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
virtual tmp< fvPatchScalarField > | clone (const DimensionedField< scalar, volMesh > &iF) const |
virtual | ~omegaWallFunctionFvPatchScalarField ()=default |
scalarField & | G (bool init=false) |
scalarField & | omega (bool init=false) |
virtual void | updateCoeffs () |
virtual void | updateWeightedCoeffs (const scalarField &weights) |
virtual void | manipulateMatrix (fvMatrix< scalar > &matrix) |
virtual void | manipulateMatrix (fvMatrix< scalar > &matrix, const scalarField &weights) |
virtual void | write (Ostream &) const |
![]() | |
TypeName ("fixedValue") | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const scalar &value) | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &, const bool valueRequired=true) | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &) | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &, const DimensionedField< scalar, volMesh > &) | |
virtual tmp< fvPatchField< scalar > > | clone () const |
virtual bool | fixesValue () const |
virtual bool | assignable () const |
virtual tmp< Field< scalar > > | valueInternalCoeffs (const tmp< scalarField > &) const |
virtual tmp< Field< scalar > > | valueBoundaryCoeffs (const tmp< scalarField > &) const |
virtual tmp< Field< scalar > > | gradientInternalCoeffs () const |
virtual tmp< Field< scalar > > | gradientBoundaryCoeffs () const |
virtual void | operator= (const UList< scalar > &) |
virtual void | operator= (const fvPatchField< scalar > &) |
virtual void | operator= (const scalar &) |
virtual void | operator+= (const fvPatchField< scalar > &) |
virtual void | operator+= (const Field< scalar > &) |
virtual void | operator+= (const scalar &) |
virtual void | operator-= (const fvPatchField< scalar > &) |
virtual void | operator-= (const Field< scalar > &) |
virtual void | operator-= (const scalar &) |
virtual void | operator*= (const fvPatchField< scalar > &) |
virtual void | operator*= (const Field< scalar > &) |
virtual void | operator*= (const scalar) |
virtual void | operator/= (const fvPatchField< scalar > &) |
virtual void | operator/= (const Field< scalar > &) |
virtual void | operator/= (const scalar) |
![]() | |
TypeName ("fvPatchField") | |
declareRunTimeSelectionTable (tmp, fvPatchField, patch,(const fvPatch &p, const DimensionedField< Type, volMesh > &iF),(p, iF)) | |
declareRunTimeSelectionTable (tmp, fvPatchField, patchMapper,(const fvPatchField< Type > &ptf, const fvPatch &p, const DimensionedField< Type, volMesh > &iF, const fvPatchFieldMapper &m),(dynamic_cast< const fvPatchFieldType & >(ptf), p, iF, m)) | |
declareRunTimeSelectionTable (tmp, fvPatchField, dictionary,(const fvPatch &p, const DimensionedField< Type, volMesh > &iF, const dictionary &dict),(p, iF, dict)) | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &) | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const Type &value) | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const word &patchType) | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const Field< Type > &) | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &, const bool valueRequired=true) | |
fvPatchField (const fvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &) | |
fvPatchField (const fvPatchField< Type > &) | |
fvPatchField (const fvPatchField< Type > &, const DimensionedField< Type, volMesh > &) | |
virtual tmp< fvPatchField< Type > > | clone (const DimensionedField< Type, volMesh > &iF) const |
Foam::tmp< Foam::fvPatchField< Type > > | NewCalculatedType (const fvPatch &p) |
Foam::tmp< Foam::fvPatchField< Type > > | NewCalculatedType (const fvPatchField< Type2 > &pf) |
virtual | ~fvPatchField ()=default |
bool | useImplicit () const noexcept |
bool | useImplicit (bool on) noexcept |
virtual bool | coupled () const |
const objectRegistry & | db () const |
const fvPatch & | patch () const |
const DimensionedField< Type, volMesh > & | internalField () const |
const Field< Type > & | primitiveField () const |
const word & | patchType () const |
word & | patchType () |
bool | updated () const |
bool | manipulatedMatrix () const |
virtual void | autoMap (const fvPatchFieldMapper &) |
virtual void | rmap (const fvPatchField< Type > &, const labelList &) |
virtual tmp< Field< Type > > | snGrad () const |
virtual tmp< Field< Type > > | snGrad (const scalarField &deltaCoeffs) const |
virtual tmp< Field< Type > > | patchInternalField () const |
virtual void | patchInternalField (Field< Type > &) const |
virtual tmp< Field< Type > > | patchNeighbourField () const |
virtual void | initEvaluate (const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) |
virtual void | evaluate (const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) |
virtual tmp< Field< Type > > | valueInternalCoeffs (const tmp< Field< scalar >> &) const |
virtual tmp< Field< Type > > | valueBoundaryCoeffs (const tmp< Field< scalar >> &) const |
virtual tmp< Field< Type > > | gradientInternalCoeffs (const scalarField &deltaCoeffs) const |
virtual tmp< Field< Type > > | gradientBoundaryCoeffs (const scalarField &deltaCoeffs) const |
virtual void | manipulateMatrix (fvMatrix< Type > &matrix) |
virtual void | manipulateMatrix (fvMatrix< Type > &matrix, const scalarField &weights) |
virtual void | manipulateMatrix (fvMatrix< Type > &matrix, const label iMatrix, const direction cmp) |
void | check (const fvPatchField< Type > &) const |
virtual void | operator= (const UList< Type > &) |
virtual void | operator= (const fvPatchField< Type > &) |
virtual void | operator= (const Type &) |
virtual void | operator+= (const fvPatchField< Type > &) |
virtual void | operator+= (const Field< Type > &) |
virtual void | operator+= (const Type &) |
virtual void | operator-= (const fvPatchField< Type > &) |
virtual void | operator-= (const Field< Type > &) |
virtual void | operator-= (const Type &) |
virtual void | operator== (const fvPatchField< Type > &) |
virtual void | operator== (const Field< Type > &) |
virtual void | operator== (const Type &) |
Protected Member Functions | |
virtual void | setMaster () |
virtual void | createAveragingWeights () |
virtual omegaWallFunctionFvPatchScalarField & | omegaPatch (const label patchi) |
virtual void | calculateTurbulenceFields (const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0) |
virtual void | calculate (const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega) |
virtual label & | master () |
Protected Attributes | |
bool | blended_ |
bool | initialised_ |
label | master_ |
scalar | beta1_ |
scalarField | G_ |
scalarField | omega_ |
List< List< scalar > > | cornerWeights_ |
Static Protected Attributes | |
static scalar | tolerance_ = 1e-5 |
Additional Inherited Members | |
![]() | |
typedef fvPatch | Patch |
typedef calculatedFvPatchField< Type > | Calculated |
![]() | |
static tmp< fvPatchField< Type > > | New (const word &, const fvPatch &, const DimensionedField< Type, volMesh > &) |
static tmp< fvPatchField< Type > > | New (const word &, const word &actualPatchType, const fvPatch &, const DimensionedField< Type, volMesh > &) |
static tmp< fvPatchField< Type > > | New (const fvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &) |
static tmp< fvPatchField< Type > > | New (const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &) |
static tmp< fvPatchField< Type > > | NewCalculatedType (const fvPatch &) |
static tmp< fvPatchField< Type > > | NewCalculatedType (const fvPatchField< Type2 > &) |
static const word & | calculatedType () |
![]() | |
static int | disallowGenericFvPatchField |
This boundary condition provides a wall constraint on the specific dissipation rate, i.e. omega
, and the turbulent kinetic energy production contribution, i.e. G
, for low- and high-Reynolds number turbulence models.
Binomial blending of the viscous and inertial sublayers (tag:ME): Menter, F., & Esch, T. (2001). Elements of industrial heat transfer prediction. In Proceedings of the 16th Brazilian Congress of Mechanical Engineering (COBEM), November 2001. vol. 20, p. 117-127. Exponential/Max blending of the viscous and inertial sublayers (tag:PH): Popovac, M., & Hanjalić, K. (2007). Compound wall treatment for RANS computation of complex turbulent flows and heat transfer. Flow, turbulence and combustion, 78(2), 177-202. DOI:10.1007/s10494-006-9067-x Tanh blending of the viscous and inertial sublayers (tag:KAS): Knopp, T., Alrutz, T., & Schwamborn, D. (2006). A grid and flow adaptive wall-function method for RANS turbulence modelling. Journal of Computational Physics, 220(1), 19-40. DOI:10.1016/j.jcp.2006.05.003
<patchName> { // Mandatory entries (unmodifiable) type omegaWallFunction; // Optional entries (unmodifiable) beta1 0.075; blending binomial2; n 2.0; // Optional (inherited) entries ... }
Property | Description | Type | Req'd | Dflt |
---|---|---|---|---|
type | Type name: omegaWallFunction | word | yes | - |
beta1 | Model coefficient | scalar | no | 0.075 |
blending | Viscous/inertial sublayer blending method | word | no | binomial2 |
n | Binomial blending exponent | scalar | no | 2.0 |
The inherited entries are elaborated in:
Options for the blending
entry:
stepwise | Stepwise switch (discontinuous) max | Maximum value switch (discontinuous) binomial2 | Binomial blending (smooth) n = 2 binomial | Binomial blending (smooth) exponential | Exponential blending (smooth) tanh | Tanh blending (smooth)
wherein omega
predictions for the viscous and inertial sublayers are blended according to the following expressions:
stepwise:
where
![]() | = | ![]() ![]() |
![]() | = | ![]() |
![]() | = | ![]() |
![]() | = | estimated wall-normal height of the cell centre in wall units |
![]() | = | estimated intersection of the viscous and inertial sublayers |
max
(PH:Eq. 27):
binomial2
(ME:Eq. 15) (default):
binomial:
where
![]() | = | Binomial blending exponent |
exponential
(PH:Eq. 32):
where (PH:Eq. 31)
![]() | = | Blending expression |
![]() | = | ![]() |
tanh
(KAS:Eqs. 33-34):
where
![]() | = | ![]() |
![]() | = | ![]() |
![]() | = | ![]() |
G
predictions for the viscous and inertial sublayers are blended in a stepwise manner, and G
below (i.e. in the viscous sublayer) is presumed to be zero.
Cmu
, kappa
, and E
are obtained from the specified nutWallFunction
in order to ensure that each patch possesses the same set of values for these coefficients.binomial2
and binomial
blending methods exist at the same time is to ensure the bitwise regression with the previous versions since binomial2
and binomial
with n=2
will yield slightly different output due to the miniscule differences in the implementation of the basic functions (i.e. pow
, sqrt
, sqr
).Definition at line 287 of file omegaWallFunctionFvPatchScalarField.H.
omegaWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Definition at line 308 of file omegaWallFunctionFvPatchScalarField.C.
Referenced by omegaWallFunctionFvPatchScalarField::clone().
omegaWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Definition at line 346 of file omegaWallFunctionFvPatchScalarField.C.
References dict, Foam::endl(), dictionary::found(), dictionary::get(), IOWarningInFunction, Foam::nl, and Foam::operator==().
omegaWallFunctionFvPatchScalarField | ( | const omegaWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Definition at line 326 of file omegaWallFunctionFvPatchScalarField.C.
omegaWallFunctionFvPatchScalarField | ( | const omegaWallFunctionFvPatchScalarField & | owfpsf | ) |
Definition at line 407 of file omegaWallFunctionFvPatchScalarField.C.
omegaWallFunctionFvPatchScalarField | ( | const omegaWallFunctionFvPatchScalarField & | owfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Definition at line 424 of file omegaWallFunctionFvPatchScalarField.C.
|
virtualdefault |
|
protectedvirtual |
Definition at line 48 of file omegaWallFunctionFvPatchScalarField.C.
References forAll, fvPatchField< scalar >::internalField(), omegaWallFunctionFvPatchScalarField::master(), omegaWallFunctionFvPatchScalarField::master_, omegaWallFunctionFvPatchScalarField::omega(), and omegaWallFunctionFvPatchScalarField::omegaPatch().
|
protectedvirtual |
Definition at line 78 of file omegaWallFunctionFvPatchScalarField.C.
References DynamicList::append(), GeometricField::boundaryField(), polyMesh::changing(), Foam::dimless, forAll, mesh, IOobject::NO_READ, IOobject::NO_WRITE, fvPatchField::patchInternalField(), fvMesh::time(), Time::timeName(), and Foam::Zero.
|
protectedvirtual |
Definition at line 138 of file omegaWallFunctionFvPatchScalarField.C.
References GeometricField::boundaryField().
Referenced by omegaWallFunctionFvPatchScalarField::setMaster().
|
protectedvirtual |
Definition at line 155 of file omegaWallFunctionFvPatchScalarField.C.
References omegaWallFunctionFvPatchScalarField::calculate(), fvPatch::faceCells(), forAll, Foam::constant::electromagnetic::G0, and fvPatchField::patch().
|
protectedvirtual |
Reimplemented in atmOmegaWallFunctionFvPatchScalarField.
Definition at line 188 of file omegaWallFunctionFvPatchScalarField.C.
References GeometricField::boundaryField(), nutWallFunctionFvPatchScalarField::Cmu(), Foam::exp(), forAll, Foam::constant::electromagnetic::G0, k, turbulenceModel::k(), nutWallFunctionFvPatchScalarField::kappa(), Foam::mag(), Foam::max(), turbulenceModel::nu(), nutWallFunctionFvPatchScalarField::nutw(), Foam::foamVersion::patch, Foam::pow(), Foam::pow025(), Foam::pow4(), fvPatchField::snGrad(), Foam::sqr(), Foam::sqrt(), Foam::tanh(), turbulenceModel::U(), y, turbulenceModel::y(), yPlus, and nutWallFunctionFvPatchScalarField::yPlusLam().
Referenced by omegaWallFunctionFvPatchScalarField::calculateTurbulenceFields().
|
inlineprotectedvirtual |
Definition at line 383 of file omegaWallFunctionFvPatchScalarField.H.
References omegaWallFunctionFvPatchScalarField::master_.
Referenced by omegaWallFunctionFvPatchScalarField::setMaster().
TypeName | ( | "omegaWallFunction" | ) |
|
inlinevirtual |
Reimplemented in atmOmegaWallFunctionFvPatchScalarField.
Definition at line 430 of file omegaWallFunctionFvPatchScalarField.H.
References omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField().
|
inlinevirtual |
Reimplemented from fixedValueFvPatchField< scalar >.
Reimplemented in atmOmegaWallFunctionFvPatchScalarField.
Definition at line 447 of file omegaWallFunctionFvPatchScalarField.H.
References omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField().
Foam::scalarField & G | ( | bool | init = false | ) |
Definition at line 444 of file omegaWallFunctionFvPatchScalarField.C.
References init(), and Foam::foamVersion::patch.
Foam::scalarField & omega | ( | bool | init = false | ) |
Definition at line 463 of file omegaWallFunctionFvPatchScalarField.C.
References init(), and Foam::foamVersion::patch.
Referenced by omegaWallFunctionFvPatchScalarField::setMaster().
|
virtual |
Reimplemented from fvPatchField< scalar >.
Definition at line 481 of file omegaWallFunctionFvPatchScalarField.C.
References forAll, Foam::constant::universal::G, Foam::constant::electromagnetic::G0, turbulenceModel::GName(), Foam::constant::atomic::group, IOobject::groupName(), Foam::foamVersion::patch, turbulenceModel::propertiesName, and fvPatchField::updateCoeffs().
|
virtual |
Reimplemented from fvPatchField< scalar >.
Definition at line 527 of file omegaWallFunctionFvPatchScalarField.C.
References forAll, Foam::constant::universal::G, Foam::constant::electromagnetic::G0, turbulenceModel::GName(), Foam::constant::atomic::group, IOobject::groupName(), Foam::foamVersion::patch, turbulenceModel::propertiesName, and fvPatchField::updateCoeffs().
|
virtual |
Definition at line 584 of file omegaWallFunctionFvPatchScalarField.C.
References fvPatchField::manipulateMatrix(), Foam::foamVersion::patch, and fvMatrix::setValues().
|
virtual |
Definition at line 600 of file omegaWallFunctionFvPatchScalarField.C.
References DynamicList::append(), Foam::expressions::patchExpr::debug, Foam::endl(), fld, forAll, fvPatchField::manipulateMatrix(), Foam::foamVersion::patch, Foam::Pout, and fvMatrix::setValues().
|
virtual |
Reimplemented from fixedValueFvPatchField< scalar >.
Reimplemented in atmOmegaWallFunctionFvPatchScalarField.
Definition at line 643 of file omegaWallFunctionFvPatchScalarField.C.
References os(), fixedValueFvPatchField< Type >::write(), and Ostream::writeEntry().
|
staticprotected |
Definition at line 323 of file omegaWallFunctionFvPatchScalarField.H.
|
protected |
Definition at line 327 of file omegaWallFunctionFvPatchScalarField.H.
|
protected |
Definition at line 330 of file omegaWallFunctionFvPatchScalarField.H.
|
protected |
Definition at line 333 of file omegaWallFunctionFvPatchScalarField.H.
Referenced by omegaWallFunctionFvPatchScalarField::master(), and omegaWallFunctionFvPatchScalarField::setMaster().
|
protected |
Definition at line 336 of file omegaWallFunctionFvPatchScalarField.H.
|
protected |
Definition at line 339 of file omegaWallFunctionFvPatchScalarField.H.
|
protected |
Definition at line 342 of file omegaWallFunctionFvPatchScalarField.H.
Definition at line 345 of file omegaWallFunctionFvPatchScalarField.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.