This boundary condition provides a wall constraint on the turbulent kinetic energy, i.e. k
, for low- and high-Reynolds number turbulence models.
More...
Public Member Functions | |
TypeName ("kLowReWallFunction") | |
kLowReWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
kLowReWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
kLowReWallFunctionFvPatchScalarField (const kLowReWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
kLowReWallFunctionFvPatchScalarField (const kLowReWallFunctionFvPatchScalarField &) | |
virtual tmp< fvPatchScalarField > | clone () const |
kLowReWallFunctionFvPatchScalarField (const kLowReWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
virtual tmp< fvPatchScalarField > | clone (const DimensionedField< scalar, volMesh > &iF) const |
virtual void | updateCoeffs () |
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 void | updateWeightedCoeffs (const scalarField &weights) |
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 Attributes | |
scalar | Ceps2_ |
scalar | Ck_ |
scalar | Bk_ |
scalar | C_ |
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, i.e. k
, for low- and high-Reynolds number turbulence models.
<patchName> { // Mandatory entries (unmodifiable) type kLowReWallFunction; // Optional entries (unmodifiable) Ceps2 1.9; Ck -0.416; Bk 8.366; C 11.0; // Optional (inherited) entries ... }
where the entries mean:
Property | Description | Type | Req'd | Dflt |
---|---|---|---|---|
type | Type name: kLowReWallFunction | word | yes | - |
Ceps2 | Model coefficient | scalar | no | 1.9 |
Ck | Model coefficient | scalar | no | -0.416 |
Bk | Model coefficient | scalar | no | 8.366 |
C | Model coefficient | scalar | no | 11.0 |
The inherited entries are elaborated in:
Viscous and inertial sublayer predictions for k
are blended in a stepwise manner:
where
![]() | = | k prediction in the viscous sublayer |
![]() | = | k prediction in the inertial sublayer |
![]() | = | estimated wall-normal height of the cell centre in wall units |
![]() | = | estimated intersection of the viscous and inertial sublayers |
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.Definition at line 161 of file kLowReWallFunctionFvPatchScalarField.H.
kLowReWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Definition at line 31 of file kLowReWallFunctionFvPatchScalarField.C.
Referenced by kLowReWallFunctionFvPatchScalarField::clone().
kLowReWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Definition at line 61 of file kLowReWallFunctionFvPatchScalarField.C.
kLowReWallFunctionFvPatchScalarField | ( | const kLowReWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Definition at line 45 of file kLowReWallFunctionFvPatchScalarField.C.
kLowReWallFunctionFvPatchScalarField | ( | const kLowReWallFunctionFvPatchScalarField & | kwfpsf | ) |
Definition at line 84 of file kLowReWallFunctionFvPatchScalarField.C.
kLowReWallFunctionFvPatchScalarField | ( | const kLowReWallFunctionFvPatchScalarField & | kwfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Definition at line 97 of file kLowReWallFunctionFvPatchScalarField.C.
TypeName | ( | "kLowReWallFunction" | ) |
|
inlinevirtual |
Definition at line 222 of file kLowReWallFunctionFvPatchScalarField.H.
References kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField().
|
inlinevirtual |
Reimplemented from fixedValueFvPatchField< scalar >.
Definition at line 239 of file kLowReWallFunctionFvPatchScalarField.H.
References kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField().
|
virtual |
Reimplemented from fvPatchField< scalar >.
Definition at line 112 of file kLowReWallFunctionFvPatchScalarField.C.
References kLowReWallFunctionFvPatchScalarField::Bk_, kLowReWallFunctionFvPatchScalarField::C_, kLowReWallFunctionFvPatchScalarField::Ceps2_, kLowReWallFunctionFvPatchScalarField::Ck_, fvPatchField< scalar >::db(), fvPatch::faceCells(), forAll, Foam::constant::atomic::group, IOobject::groupName(), fvPatch::index(), fvPatchField< scalar >::internalField(), k, turbulenceModel::k(), Foam::log(), objectRegistry::lookupObject(), Foam::max(), turbulenceModel::nu(), nutWallFunctionFvPatchScalarField::nutw(), fvPatchField< scalar >::patch(), Foam::pow025(), Foam::pow3(), turbulenceModel::propertiesName, Foam::sqr(), Foam::sqrt(), fvPatchField< Type >::updateCoeffs(), fvPatchField< scalar >::updated(), uTau, y, turbulenceModel::y(), and yPlus.
|
virtual |
Reimplemented from fixedValueFvPatchField< scalar >.
Definition at line 176 of file kLowReWallFunctionFvPatchScalarField.C.
References os(), fixedValueFvPatchField< Type >::write(), and Ostream::writeEntry().
|
protected |
Definition at line 170 of file kLowReWallFunctionFvPatchScalarField.H.
Referenced by kLowReWallFunctionFvPatchScalarField::updateCoeffs().
|
protected |
Definition at line 173 of file kLowReWallFunctionFvPatchScalarField.H.
Referenced by kLowReWallFunctionFvPatchScalarField::updateCoeffs().
|
protected |
Definition at line 176 of file kLowReWallFunctionFvPatchScalarField.H.
Referenced by kLowReWallFunctionFvPatchScalarField::updateCoeffs().
|
protected |
Definition at line 179 of file kLowReWallFunctionFvPatchScalarField.H.
Referenced by kLowReWallFunctionFvPatchScalarField::updateCoeffs().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.