This boundary condition provides a wall constraint on the turbulent kinetic energy dissipation rate, i.e. epsilon
, and the turbulent kinetic energy production contribution, i.e. G
, for low- and high-Reynolds number turbulence models.
More...
Public Member Functions | |
TypeName ("epsilonWallFunction") | |
epsilonWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
epsilonWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
epsilonWallFunctionFvPatchScalarField (const epsilonWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
epsilonWallFunctionFvPatchScalarField (const epsilonWallFunctionFvPatchScalarField &) | |
virtual tmp< fvPatchScalarField > | clone () const |
epsilonWallFunctionFvPatchScalarField (const epsilonWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
virtual tmp< fvPatchScalarField > | clone (const DimensionedField< scalar, volMesh > &iF) const |
virtual | ~epsilonWallFunctionFvPatchScalarField ()=default |
scalarField & | G (bool init=false) |
scalarField & | epsilon (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 epsilonWallFunctionFvPatchScalarField & | epsilonPatch (const label patchi) |
virtual void | calculateTurbulenceFields (const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0) |
virtual void | calculate (const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon) |
virtual label & | master () |
Protected Attributes | |
const bool | lowReCorrection_ |
bool | initialised_ |
label | master_ |
scalarField | G_ |
scalarField | epsilon_ |
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 turbulent kinetic energy dissipation rate, i.e. epsilon
, 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
<patchName> { // Mandatory entries (unmodifiable) type epsilonWallFunction; // Optional entries (unmodifiable) lowReCorrection false; blending stepwise; n 2.0; // Optional (inherited) entries ... }
where the entries mean:
Property | Description | Type | Req'd | Dflt |
---|---|---|---|---|
type | Type name: epsilonWallFunction | word | yes | - |
lowReCorrection | Flag: apply low-Re correction | bool | no | false |
blending | Viscous/inertial sublayer blending method | word | no | stepwise |
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) binomial | Binomial blending (smooth) exponential | Exponential blending (smooth)
wherein epsilon
predictions for the viscous and inertial sublayers are blended according to the following expressions:
stepwise
(default):
where
![]() | = | ![]() ![]() |
![]() | = | ![]() |
![]() | = | ![]() |
![]() | = | estimated wall-normal height of the cell centre in wall units |
![]() | = | estimated intersection of the viscous and inertial sublayers |
max
(PH:Eq. 27):
binomial
(ME:Eqs. 15-16):
where
![]() | = | Binomial blending exponent |
exponential
(PH:Eq. 32):
where (PH:p. 193)
![]() | = | ![]() |
![]() | = | ![]() |
![]() | = | Blending expression for ![]() |
![]() | = | Blending expression for ![]() |
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.lowReCorrection
operates with only stepwise
blending treatment to ensure the backward compatibility.lowReCorrection
is on
, stepwise
blending treatment is fully active.lowReCorrection
is off
, only the inertial sublayer prediction is used in the wall function, hence high-Re mode operation.Definition at line 255 of file epsilonWallFunctionFvPatchScalarField.H.
epsilonWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Definition at line 298 of file epsilonWallFunctionFvPatchScalarField.C.
Referenced by epsilonWallFunctionFvPatchScalarField::clone().
epsilonWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Definition at line 338 of file epsilonWallFunctionFvPatchScalarField.C.
References Foam::operator==().
epsilonWallFunctionFvPatchScalarField | ( | const epsilonWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Definition at line 317 of file epsilonWallFunctionFvPatchScalarField.C.
epsilonWallFunctionFvPatchScalarField | ( | const epsilonWallFunctionFvPatchScalarField & | ewfpsf | ) |
Definition at line 377 of file epsilonWallFunctionFvPatchScalarField.C.
epsilonWallFunctionFvPatchScalarField | ( | const epsilonWallFunctionFvPatchScalarField & | ewfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Definition at line 395 of file epsilonWallFunctionFvPatchScalarField.C.
|
virtualdefault |
|
protectedvirtual |
Definition at line 47 of file epsilonWallFunctionFvPatchScalarField.C.
References epsilonWallFunctionFvPatchScalarField::epsilon(), epsilonWallFunctionFvPatchScalarField::epsilonPatch(), forAll, fvPatchField< scalar >::internalField(), epsilonWallFunctionFvPatchScalarField::master(), and epsilonWallFunctionFvPatchScalarField::master_.
|
protectedvirtual |
Definition at line 77 of file epsilonWallFunctionFvPatchScalarField.C.
References DynamicList::append(), GeometricField::boundaryField(), polyMesh::changing(), Foam::dimless, epsilon, forAll, mesh, IOobject::NO_READ, IOobject::NO_WRITE, fvPatchField::patchInternalField(), fvMesh::time(), Time::timeName(), and Foam::Zero.
|
protectedvirtual |
Definition at line 138 of file epsilonWallFunctionFvPatchScalarField.C.
References epsilon.
Referenced by epsilonWallFunctionFvPatchScalarField::setMaster().
|
protectedvirtual |
Definition at line 155 of file epsilonWallFunctionFvPatchScalarField.C.
References epsilonWallFunctionFvPatchScalarField::calculate(), Foam::constant::electromagnetic::epsilon0, fvPatch::faceCells(), forAll, Foam::constant::electromagnetic::G0, fvPatchField::patch(), and turbulence.
|
protectedvirtual |
Reimplemented in atmEpsilonWallFunctionFvPatchScalarField.
Definition at line 188 of file epsilonWallFunctionFvPatchScalarField.C.
References GeometricField::boundaryField(), nutWallFunctionFvPatchScalarField::Cmu(), Foam::constant::electromagnetic::epsilon0, 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(), turbulenceModel::U(), y, turbulenceModel::y(), yPlus, and nutWallFunctionFvPatchScalarField::yPlusLam().
Referenced by epsilonWallFunctionFvPatchScalarField::calculateTurbulenceFields().
|
inlineprotectedvirtual |
Definition at line 345 of file epsilonWallFunctionFvPatchScalarField.H.
References epsilonWallFunctionFvPatchScalarField::master_.
Referenced by epsilonWallFunctionFvPatchScalarField::setMaster().
TypeName | ( | "epsilonWallFunction" | ) |
|
inlinevirtual |
Reimplemented in atmEpsilonWallFunctionFvPatchScalarField.
Definition at line 392 of file epsilonWallFunctionFvPatchScalarField.H.
References epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField().
|
inlinevirtual |
Reimplemented from fixedValueFvPatchField< scalar >.
Reimplemented in atmEpsilonWallFunctionFvPatchScalarField.
Definition at line 409 of file epsilonWallFunctionFvPatchScalarField.H.
References epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField().
Foam::scalarField & G | ( | bool | init = false | ) |
Definition at line 415 of file epsilonWallFunctionFvPatchScalarField.C.
References init(), and Foam::foamVersion::patch.
Foam::scalarField & epsilon | ( | bool | init = false | ) |
Definition at line 434 of file epsilonWallFunctionFvPatchScalarField.C.
References init(), and Foam::foamVersion::patch.
Referenced by epsilonWallFunctionFvPatchScalarField::setMaster().
|
virtual |
Reimplemented from fvPatchField< scalar >.
Definition at line 452 of file epsilonWallFunctionFvPatchScalarField.C.
References epsilon, Foam::constant::electromagnetic::epsilon0, 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 498 of file epsilonWallFunctionFvPatchScalarField.C.
References epsilon, Foam::constant::electromagnetic::epsilon0, 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 555 of file epsilonWallFunctionFvPatchScalarField.C.
References fvPatchField::manipulateMatrix(), Foam::foamVersion::patch, and fvMatrix::setValues().
|
virtual |
Definition at line 571 of file epsilonWallFunctionFvPatchScalarField.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 atmEpsilonWallFunctionFvPatchScalarField.
Definition at line 614 of file epsilonWallFunctionFvPatchScalarField.C.
References os(), fixedValueFvPatchField< Type >::write(), and Ostream::writeEntry().
|
staticprotected |
Definition at line 289 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
Definition at line 292 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
Definition at line 295 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
Definition at line 298 of file epsilonWallFunctionFvPatchScalarField.H.
Referenced by epsilonWallFunctionFvPatchScalarField::master(), and epsilonWallFunctionFvPatchScalarField::setMaster().
|
protected |
Definition at line 301 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
Definition at line 304 of file epsilonWallFunctionFvPatchScalarField.H.
Definition at line 307 of file epsilonWallFunctionFvPatchScalarField.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.