This boundary condition provides a wall constraint on the turbulent kinematic viscosity, i.e. nut
, when using wall functions, based on turbulent kinetic energy, k
.
More...
Public Member Functions | |
TypeName ("nutkWallFunction") | |
Runtime type information. More... | |
nutkWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
nutkWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
nutkWallFunctionFvPatchScalarField (const nutkWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
nutkWallFunctionFvPatchScalarField (const nutkWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
virtual tmp< fvPatchScalarField > | clone () const |
Construct and return a clone. More... | |
nutkWallFunctionFvPatchScalarField (const nutkWallFunctionFvPatchScalarField &, 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... | |
virtual tmp< scalarField > | yPlus () const |
Calculate and return the yPlus at the boundary. More... | |
Public Member Functions inherited from nutWallFunctionFvPatchScalarField | |
TypeName ("nutWallFunction") | |
Runtime type information. More... | |
nutWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
nutWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
scalar | Cmu () const |
Return Cmu. More... | |
scalar | kappa () const |
Return kappa. More... | |
scalar | E () const |
Return E. More... | |
scalar | yPlusLam () const |
Return the y+ at the edge of the viscous sublayer. More... | |
virtual void | updateCoeffs () |
Update the coefficients associated with the patch field. More... | |
virtual void | write (Ostream &) const |
Write. More... | |
Protected Member Functions | |
virtual tmp< scalarField > | calcNut () const |
Calculate the turbulence viscosity. More... | |
Protected Member Functions inherited from nutWallFunctionFvPatchScalarField | |
virtual const volVectorField & | U (const turbulenceModel &turb) const |
virtual void | checkType () |
Check the type of the patch. More... | |
virtual void | writeLocalEntries (Ostream &) const |
Write local wall function variables. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from nutWallFunctionFvPatchScalarField | |
static const nutWallFunctionFvPatchScalarField & | nutw (const turbulenceModel &turbModel, const label patchi) |
Return the nut patchField for the given wall patch. More... | |
static scalar | yPlusLam (const scalar kappa, const scalar E) |
Calculate the y+ at the edge of the viscous sublayer. More... | |
Protected Attributes inherited from nutWallFunctionFvPatchScalarField | |
word | UName_ |
Name of velocity field. More... | |
scalar | Cmu_ |
Cmu coefficient. More... | |
scalar | kappa_ |
Von Karman constant. More... | |
scalar | E_ |
E coefficient. More... | |
scalar | yPlusLam_ |
Estimated y+ value at the edge of the viscous sublayer. More... | |
This boundary condition provides a wall constraint on the turbulent kinematic viscosity, i.e. nut
, when using wall functions, based on turbulent kinetic energy, k
.
Property | Description | Required | Default value |
---|---|---|---|
Cmu | Model coefficient | no | 0.09 |
kappa | Von Karman constant | no | 0.41 |
E | Model coefficient | no | 9.8 |
Example of the boundary condition specification:
<patchName> { // Mandatory entries type nutkWallFunction; // Optional entries }
Definition at line 98 of file nutkWallFunctionFvPatchScalarField.H.
nutkWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 83 of file nutkWallFunctionFvPatchScalarField.C.
Referenced by nutkWallFunctionFvPatchScalarField::clone().
nutkWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 105 of file nutkWallFunctionFvPatchScalarField.C.
nutkWallFunctionFvPatchScalarField | ( | const nutkWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given nutkWallFunctionFvPatchScalarField onto a new patch
Definition at line 93 of file nutkWallFunctionFvPatchScalarField.C.
nutkWallFunctionFvPatchScalarField | ( | const nutkWallFunctionFvPatchScalarField & | wfpsf | ) |
Construct as copy.
Definition at line 116 of file nutkWallFunctionFvPatchScalarField.C.
nutkWallFunctionFvPatchScalarField | ( | const nutkWallFunctionFvPatchScalarField & | wfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 125 of file nutkWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Calculate the turbulence viscosity.
Implements nutWallFunctionFvPatchScalarField.
Reimplemented in nutkRoughWallFunctionFvPatchScalarField, nutkAtmRoughWallFunctionFvPatchScalarField, and nutkFilmWallFunctionFvPatchScalarField.
Definition at line 40 of file nutkWallFunctionFvPatchScalarField.C.
References nutWallFunctionFvPatchScalarField::Cmu_, nutWallFunctionFvPatchScalarField::E_, forAll, Foam::constant::atomic::group, IOobject::groupName(), k, turbulenceModel::k(), nutWallFunctionFvPatchScalarField::kappa_, Foam::log(), turbulenceModel::nu(), nutWallFunctionFvPatchScalarField::nutw(), Foam::foamVersion::patch, Foam::pow025(), turbulenceModel::propertiesName, tmp< T >::ref(), Foam::sqrt(), y, turbulenceModel::y(), nutkWallFunctionFvPatchScalarField::yPlus(), nutWallFunctionFvPatchScalarField::yPlusLam_, and Foam::Zero.
TypeName | ( | "nutkWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Reimplemented in nutkRoughWallFunctionFvPatchScalarField, nutkAtmRoughWallFunctionFvPatchScalarField, and nutkFilmWallFunctionFvPatchScalarField.
Definition at line 151 of file nutkWallFunctionFvPatchScalarField.H.
References nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField().
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Reimplemented in nutkRoughWallFunctionFvPatchScalarField, nutkAtmRoughWallFunctionFvPatchScalarField, and nutkFilmWallFunctionFvPatchScalarField.
Definition at line 168 of file nutkWallFunctionFvPatchScalarField.H.
References nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField().
|
virtual |
Calculate and return the yPlus at the boundary.
Implements nutWallFunctionFvPatchScalarField.
Reimplemented in nutkFilmWallFunctionFvPatchScalarField.
Definition at line 137 of file nutkWallFunctionFvPatchScalarField.C.
References Foam::constant::atomic::group, IOobject::groupName(), k, turbulenceModel::k(), turbulenceModel::nu(), Foam::foamVersion::patch, Foam::pow025(), turbulenceModel::propertiesName, Foam::sqrt(), y, and turbulenceModel::y().
Referenced by nutkWallFunctionFvPatchScalarField::calcNut(), nutkAtmRoughWallFunctionFvPatchScalarField::calcNut(), and nutkRoughWallFunctionFvPatchScalarField::calcNut().