This boundary condition provides a wall constraint on specific turbulent kinetic energy dissipation rate, i.e. omega
, for low- and high-Reynolds number turbulence models.
More...
Public Member Functions | |
TypeName ("omegaWallFunction") | |
Runtime type information. More... | |
omegaWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
omegaWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
virtual tmp< fvPatchScalarField > | clone () const |
Construct and return a clone. More... | |
omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
virtual tmp< fvPatchScalarField > | clone (const DimensionedField< scalar, volMesh > &iF) const |
Construct and return a clone setting internal field reference. More... | |
scalarField & | G (bool init=false) |
Return non-const access to the master's G field. More... | |
scalarField & | omega (bool init=false) |
Return non-const access to the master's omega field. More... | |
virtual void | updateCoeffs () |
Update the coefficients associated with the patch field. More... | |
virtual void | updateWeightedCoeffs (const scalarField &weights) |
Update the coefficients associated with the patch field. More... | |
virtual void | manipulateMatrix (fvMatrix< scalar > &matrix) |
Manipulate matrix. More... | |
virtual void | manipulateMatrix (fvMatrix< scalar > &matrix, const scalarField &weights) |
Manipulate matrix with given weights. More... | |
virtual void | write (Ostream &) const |
Write. More... | |
Public Member Functions inherited from fixedValueFvPatchField< scalar > | |
TypeName ("fixedValue") | |
Runtime type information. More... | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const scalar &value) | |
Construct from patch, internal field and value. More... | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &, const bool valueRequired=true) | |
Construct from patch, internal field and dictionary. More... | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
Construct by mapping the given fixedValueFvPatchField<Type> More... | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &) | |
Construct as copy. More... | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &, const DimensionedField< scalar, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
virtual bool | fixesValue () const |
Return true if this patch field fixes a value. More... | |
virtual bool | assignable () const |
Return false: this patch field is not altered by assignment. More... | |
virtual tmp< Field< scalar > > | valueInternalCoeffs (const tmp< scalarField > &) const |
Return the matrix diagonal coefficients corresponding to the. More... | |
virtual tmp< Field< scalar > > | valueBoundaryCoeffs (const tmp< scalarField > &) const |
Return the matrix source coefficients corresponding to the. More... | |
virtual tmp< Field< scalar > > | gradientInternalCoeffs () const |
Return the matrix diagonal coefficients corresponding to the. More... | |
virtual tmp< Field< scalar > > | gradientBoundaryCoeffs () const |
Return the matrix source coefficients corresponding to the. More... | |
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) |
Public Member Functions inherited from fvPatchField< scalar > | |
TypeName ("fvPatchField") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (tmp, fvPatchField, patch,(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF),(p, iF)) | |
declareRunTimeSelectionTable (tmp, fvPatchField, patchMapper,(const fvPatchField< scalar > &ptf, const fvPatch &p, const DimensionedField< scalar, volMesh > &iF, const fvPatchFieldMapper &m),(dynamic_cast< const fvPatchFieldType & >(ptf), p, iF, m)) | |
declareRunTimeSelectionTable (tmp, fvPatchField, dictionary,(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF, const dictionary &dict),(p, iF, dict)) | |
fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const scalar &value) | |
Construct from patch, internal field and value. More... | |
fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const word &patchType) | |
Construct from patch and internal field and patch type. More... | |
fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const Field< scalar > &) | |
Construct from patch and internal field and patch field. More... | |
fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &, const bool valueRequired=true) | |
Construct from patch, internal field and dictionary. More... | |
fvPatchField (const fvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
Construct by mapping the given fvPatchField onto a new patch. More... | |
fvPatchField (const fvPatchField< scalar > &) | |
Construct as copy. More... | |
fvPatchField (const fvPatchField< scalar > &, const DimensionedField< scalar, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
Foam::tmp< Foam::fvPatchField< scalar > > | NewCalculatedType (const fvPatchField< Type2 > &pf) |
virtual | ~fvPatchField ()=default |
Destructor. More... | |
virtual bool | coupled () const |
Return true if this patch field is coupled. More... | |
const objectRegistry & | db () const |
Return local objectRegistry. More... | |
const fvPatch & | patch () const |
Return patch. More... | |
const DimensionedField< scalar, volMesh > & | internalField () const |
Return dimensioned internal field reference. More... | |
const Field< scalar > & | primitiveField () const |
Return internal field reference. More... | |
const word & | patchType () const |
Optional patch type. More... | |
word & | patchType () |
Optional patch type. More... | |
bool | updated () const |
Return true if the boundary condition has already been updated. More... | |
bool | manipulatedMatrix () const |
Return true if the matrix has already been manipulated. More... | |
virtual void | autoMap (const fvPatchFieldMapper &) |
Map (and resize as needed) from self given a mapping object. More... | |
virtual void | rmap (const fvPatchField< scalar > &, const labelList &) |
Reverse map the given fvPatchField onto this fvPatchField. More... | |
virtual tmp< Field< scalar > > | snGrad () const |
Return patch-normal gradient. More... | |
virtual tmp< Field< scalar > > | snGrad (const scalarField &deltaCoeffs) const |
Return patch-normal gradient for coupled-patches. More... | |
virtual tmp< Field< scalar > > | patchInternalField () const |
Return internal field next to patch as patch field. More... | |
virtual void | patchInternalField (Field< scalar > &) const |
Return internal field next to patch as patch field. More... | |
virtual tmp< Field< scalar > > | patchNeighbourField () const |
Return patchField on the opposite patch of a coupled patch. More... | |
virtual void | initEvaluate (const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) |
Initialise the evaluation of the patch field. More... | |
virtual void | evaluate (const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) |
Evaluate the patch field, sets Updated to false. More... | |
virtual tmp< Field< scalar > > | valueInternalCoeffs (const tmp< Field< scalar >> &) const |
Return the matrix diagonal coefficients corresponding to the. More... | |
virtual tmp< Field< scalar > > | valueBoundaryCoeffs (const tmp< Field< scalar >> &) const |
Return the matrix source coefficients corresponding to the. More... | |
virtual tmp< Field< scalar > > | gradientInternalCoeffs (const scalarField &deltaCoeffs) const |
Return the matrix diagonal coefficients corresponding to the. More... | |
virtual tmp< Field< scalar > > | gradientBoundaryCoeffs (const scalarField &deltaCoeffs) const |
Return the matrix source coefficients corresponding to the. More... | |
void | check (const fvPatchField< scalar > &) const |
Check fvPatchField<Type> against given fvPatchField<Type> More... | |
virtual void | operator== (const fvPatchField< scalar > &) |
virtual void | operator== (const Field< scalar > &) |
virtual void | operator== (const scalar &) |
Protected Member Functions | |
virtual void | setMaster () |
virtual void | createAveragingWeights () |
virtual omegaWallFunctionFvPatchScalarField & | omegaPatch (const label patchi) |
Helper function to return non-const access to an omega patch. More... | |
virtual void | calculateTurbulenceFields (const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0) |
Main driver to calculate the turbulence fields. More... | |
virtual void | calculate (const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega) |
Calculate the omega and G. More... | |
virtual label & | master () |
Return non-const access to the master patch ID. More... | |
Protected Attributes | |
bool | blended_ |
Blending switch. More... | |
bool | initialised_ |
Initialised flag. More... | |
label | master_ |
Master patch ID. More... | |
scalar | beta1_ |
beta1 coefficient More... | |
scalarField | G_ |
Local copy of turbulence G field. More... | |
scalarField | omega_ |
Local copy of turbulence omega field. More... | |
List< List< scalar > > | cornerWeights_ |
List of averaging corner weights. More... | |
Static Protected Attributes | |
static scalar | tolerance_ = 1e-5 |
Tolerance used in weighted calculations. More... | |
Additional Inherited Members | |
Public Types inherited from fvPatchField< scalar > | |
typedef fvPatch | Patch |
typedef calculatedFvPatchField< scalar > | Calculated |
Static Public Member Functions inherited from fvPatchField< scalar > | |
static tmp< fvPatchField< scalar > > | New (const word &, const fvPatch &, const DimensionedField< scalar, volMesh > &) |
Return a pointer to a new patchField created on freestore given. More... | |
static tmp< fvPatchField< scalar > > | New (const word &, const word &actualPatchType, const fvPatch &, const DimensionedField< scalar, volMesh > &) |
Return a pointer to a new patchField created on freestore given. More... | |
static tmp< fvPatchField< scalar > > | New (const fvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) |
Return a pointer to a new patchField created on freestore from. More... | |
static tmp< fvPatchField< scalar > > | New (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) |
Return a pointer to a new patchField created on freestore. More... | |
static tmp< fvPatchField< scalar > > | NewCalculatedType (const fvPatch &) |
Return a pointer to a new calculatedFvPatchField created on. More... | |
static tmp< fvPatchField< scalar > > | NewCalculatedType (const fvPatchField< Type2 > &) |
Return a pointer to a new calculatedFvPatchField created on. More... | |
static const word & | calculatedType () |
Return the type of the calculated for of fvPatchField. More... | |
Static Public Attributes inherited from fvPatchField< scalar > | |
static int | disallowGenericFvPatchField |
Debug switch to disallow the use of genericFvPatchField. More... | |
This boundary condition provides a wall constraint on specific turbulent kinetic energy dissipation rate, i.e. omega
, for low- and high-Reynolds number turbulence models.
The near-wall omega may be either blended between the viscous region and logarithmic region values using:
\[ \omega = sqrt(\omega_{vis}^2 + \omega_{log}^2) \]
where
\( \omega_{vis} \) | = | omega in viscous region |
\( \omega_{log} \) | = | omega in logarithmic region |
see eq.(15) of:
Menter, F., Esch, T. "Elements of Industrial Heat Transfer Prediction" 16th Brazilian Congress of Mechanical Engineering (COBEM), Nov. 2001
or switched between these values based on the viscous-to-turbulent y+
value derived from kappa
and E
.
Property | Description | Required | Default value |
---|---|---|---|
beta1 | Model coefficient | no | 0.075 |
blended | Blending switch | no | false |
Example of the boundary condition specification:
<patchName> { // Mandatory entries type omegaWallFunction; // Optional entries }
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.Some tests have shown that the stepwise switching method provides more accurate predictions for 10 < y+ < 30 when used with high Reynolds number wall-functions.
In addition, the stepwise switching method provides accurate predictions when used with continuous wall-functions.
Definition at line 140 of file omegaWallFunctionFvPatchScalarField.H.
omegaWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 267 of file omegaWallFunctionFvPatchScalarField.C.
Referenced by omegaWallFunctionFvPatchScalarField::clone().
omegaWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 303 of file omegaWallFunctionFvPatchScalarField.C.
References Foam::operator==().
omegaWallFunctionFvPatchScalarField | ( | const omegaWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given omegaWallFunctionFvPatchScalarField onto a new patch
Definition at line 284 of file omegaWallFunctionFvPatchScalarField.C.
omegaWallFunctionFvPatchScalarField | ( | const omegaWallFunctionFvPatchScalarField & | owfpsf | ) |
Construct as copy.
Definition at line 324 of file omegaWallFunctionFvPatchScalarField.C.
omegaWallFunctionFvPatchScalarField | ( | const omegaWallFunctionFvPatchScalarField & | owfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 340 of file omegaWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Set the master patch - master is responsible for updating all wall function patches
Definition at line 42 of file omegaWallFunctionFvPatchScalarField.C.
References forAll, fvPatchField< scalar >::internalField(), omegaWallFunctionFvPatchScalarField::master(), omegaWallFunctionFvPatchScalarField::master_, omegaWallFunctionFvPatchScalarField::omega(), and omegaWallFunctionFvPatchScalarField::omegaPatch().
|
protectedvirtual |
Create the averaging weights for cells which are bounded by multiple wall function faces
Definition at line 72 of file omegaWallFunctionFvPatchScalarField.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), polyMesh::changing(), Foam::dimless, forAll, mesh, IOobject::NO_READ, IOobject::NO_WRITE, fvPatchField< Type >::patchInternalField(), fvMesh::time(), Time::timeName(), and Foam::Zero.
|
protectedvirtual |
Helper function to return non-const access to an omega patch.
Definition at line 134 of file omegaWallFunctionFvPatchScalarField.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField().
Referenced by omegaWallFunctionFvPatchScalarField::setMaster().
|
protectedvirtual |
Main driver to calculate the turbulence fields.
Definition at line 151 of file omegaWallFunctionFvPatchScalarField.C.
References omegaWallFunctionFvPatchScalarField::calculate(), fvPatch::faceCells(), forAll, Foam::constant::electromagnetic::G0, and fvPatchField< Type >::patch().
|
protectedvirtual |
Calculate the omega and G.
Definition at line 184 of file omegaWallFunctionFvPatchScalarField.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), nutWallFunctionFvPatchScalarField::Cmu(), forAll, Foam::constant::electromagnetic::G0, k, turbulenceModel::k(), nutWallFunctionFvPatchScalarField::kappa(), Foam::mag(), turbulenceModel::nu(), nutWallFunctionFvPatchScalarField::nutw(), Foam::foamVersion::patch, Foam::pow025(), fvPatchField< Type >::snGrad(), Foam::sqr(), Foam::sqrt(), turbulenceModel::U(), y, turbulenceModel::y(), yPlus, and nutWallFunctionFvPatchScalarField::yPlusLam().
Referenced by omegaWallFunctionFvPatchScalarField::calculateTurbulenceFields().
|
inlineprotectedvirtual |
Return non-const access to the master patch ID.
Definition at line 208 of file omegaWallFunctionFvPatchScalarField.H.
References omegaWallFunctionFvPatchScalarField::master_.
Referenced by omegaWallFunctionFvPatchScalarField::setMaster().
TypeName | ( | "omegaWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Reimplemented from fixedValueFvPatchField< scalar >.
Definition at line 255 of file omegaWallFunctionFvPatchScalarField.H.
References omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField().
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Reimplemented from fixedValueFvPatchField< scalar >.
Definition at line 272 of file omegaWallFunctionFvPatchScalarField.H.
References omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField().
Foam::scalarField & G | ( | bool | init = false | ) |
Return non-const access to the master's G field.
Definition at line 359 of file omegaWallFunctionFvPatchScalarField.C.
References Foam::foamVersion::patch.
Foam::scalarField & omega | ( | bool | init = false | ) |
Return non-const access to the master's omega field.
Definition at line 378 of file omegaWallFunctionFvPatchScalarField.C.
References Foam::foamVersion::patch.
Referenced by omegaWallFunctionFvPatchScalarField::setMaster().
|
virtual |
Update the coefficients associated with the patch field.
Reimplemented from fvPatchField< scalar >.
Definition at line 396 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< Type >::updateCoeffs().
|
virtual |
Update the coefficients associated with the patch field.
Reimplemented from fvPatchField< scalar >.
Definition at line 442 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< Type >::updateCoeffs().
|
virtual |
Manipulate matrix.
Reimplemented from fvPatchField< scalar >.
Definition at line 499 of file omegaWallFunctionFvPatchScalarField.C.
References fvPatchField< Type >::manipulateMatrix(), Foam::foamVersion::patch, and fvMatrix< Type >::setValues().
|
virtual |
Manipulate matrix with given weights.
Reimplemented from fvPatchField< scalar >.
Definition at line 515 of file omegaWallFunctionFvPatchScalarField.C.
References Foam::expressions::patchExpr::debug, Foam::endl(), fld, forAll, fvPatchField< Type >::manipulateMatrix(), Foam::foamVersion::patch, Foam::Pout, and fvMatrix< Type >::setValues().
|
virtual |
Write.
Reimplemented from fixedValueFvPatchField< scalar >.
Definition at line 558 of file omegaWallFunctionFvPatchScalarField.C.
References fixedValueFvPatchField< Type >::write(), and Ostream::writeEntry().
|
staticprotected |
Tolerance used in weighted calculations.
Definition at line 149 of file omegaWallFunctionFvPatchScalarField.H.
|
protected |
Blending switch.
Definition at line 152 of file omegaWallFunctionFvPatchScalarField.H.
|
protected |
Initialised flag.
Definition at line 155 of file omegaWallFunctionFvPatchScalarField.H.
|
protected |
Master patch ID.
Definition at line 158 of file omegaWallFunctionFvPatchScalarField.H.
Referenced by omegaWallFunctionFvPatchScalarField::master(), and omegaWallFunctionFvPatchScalarField::setMaster().
|
protected |
beta1 coefficient
Definition at line 161 of file omegaWallFunctionFvPatchScalarField.H.
|
protected |
Local copy of turbulence G field.
Definition at line 164 of file omegaWallFunctionFvPatchScalarField.H.
|
protected |
Local copy of turbulence omega field.
Definition at line 167 of file omegaWallFunctionFvPatchScalarField.H.
List of averaging corner weights.
Definition at line 170 of file omegaWallFunctionFvPatchScalarField.H.