This boundary condition provides a wall constraint on the turbulent viscosity, i.e. nut
, based on the turbulent kinetic energy, i.e. k
, for for low- and high-Reynolds number turbulence models.
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... | |
![]() | |
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 estimated y+ at the two-sublayer intersection. More... | |
scalar | blend (const scalar nutVis, const scalar nutLog, const scalar yPlus) const |
Return the blended nut according to the chosen blending treatment. 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 turbulent viscosity. More... | |
![]() | |
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 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) |
Estimate the y+ at the intersection of the two sublayers. More... | |
![]() | |
enum | blendingType { STEPWISE, MAX, BINOMIAL, EXPONENTIAL } |
Options for the blending treatment of viscous and inertial sublayers. More... | |
![]() | |
enum blendingType | blending_ |
Blending treatment (default = blendingType::STEPWISE) More... | |
const scalar | n_ |
word | UName_ |
Name of velocity field. More... | |
scalar | Cmu_ |
Empirical model coefficient. More... | |
scalar | kappa_ |
von Kármán constant More... | |
scalar | E_ |
Wall roughness parameter. More... | |
scalar | yPlusLam_ |
![]() | |
static const Enum< blendingType > | blendingTypeNames |
Names for blendingType. More... | |
This boundary condition provides a wall constraint on the turbulent viscosity, i.e. nut
, based on the turbulent kinetic energy, i.e. k
, for for low- and high-Reynolds number turbulence models.
<patchName> { // Mandatory entries (unmodifiable) type nutkWallFunction; // Optional (inherited) entries ... }
where the entries mean:
Property | Description | Type | Req'd | Dflt |
---|---|---|---|---|
type | Type name: nutkWallFunction | word | yes | - |
The inherited entries are elaborated in:
Definition at line 94 of file nutkWallFunctionFvPatchScalarField.H.
nutkWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 86 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 108 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 96 of file nutkWallFunctionFvPatchScalarField.C.
nutkWallFunctionFvPatchScalarField | ( | const nutkWallFunctionFvPatchScalarField & | wfpsf | ) |
Construct as copy.
Definition at line 119 of file nutkWallFunctionFvPatchScalarField.C.
nutkWallFunctionFvPatchScalarField | ( | const nutkWallFunctionFvPatchScalarField & | wfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 128 of file nutkWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Calculate the turbulent viscosity.
Implements nutWallFunctionFvPatchScalarField.
Reimplemented in atmNutWallFunctionFvPatchScalarField, atmNutkWallFunctionFvPatchScalarField, nutkRoughWallFunctionFvPatchScalarField, and nutkFilmWallFunctionFvPatchScalarField.
Definition at line 39 of file nutkWallFunctionFvPatchScalarField.C.
References nutWallFunctionFvPatchScalarField::blend(), nutWallFunctionFvPatchScalarField::Cmu_, Foam::constant::electromagnetic::e, nutWallFunctionFvPatchScalarField::E_, forAll, Foam::constant::atomic::group, IOobject::groupName(), k, turbulenceModel::k(), nutWallFunctionFvPatchScalarField::kappa_, Foam::log(), Foam::max(), turbulenceModel::nu(), nutWallFunctionFvPatchScalarField::nutw(), Foam::foamVersion::patch, Foam::pow025(), turbulenceModel::propertiesName, tmp< T >::ref(), Foam::sqrt(), y, turbulenceModel::y(), nutkWallFunctionFvPatchScalarField::yPlus(), and Foam::Zero.
TypeName | ( | "nutkWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Reimplemented in atmNutWallFunctionFvPatchScalarField, atmNutkWallFunctionFvPatchScalarField, nutkRoughWallFunctionFvPatchScalarField, and nutkFilmWallFunctionFvPatchScalarField.
Definition at line 147 of file nutkWallFunctionFvPatchScalarField.H.
References nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField().
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Reimplemented in atmNutWallFunctionFvPatchScalarField, atmNutkWallFunctionFvPatchScalarField, nutkRoughWallFunctionFvPatchScalarField, and nutkFilmWallFunctionFvPatchScalarField.
Definition at line 164 of file nutkWallFunctionFvPatchScalarField.H.
References nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField().
|
virtual |
Calculate and return the yPlus at the boundary.
Implements nutWallFunctionFvPatchScalarField.
Reimplemented in nutkFilmWallFunctionFvPatchScalarField.
Definition at line 140 of file nutkWallFunctionFvPatchScalarField.C.
References forAll, Foam::constant::atomic::group, IOobject::groupName(), k, Foam::mag(), tmp< T >::New(), Foam::foamVersion::patch, Foam::pow025(), turbulenceModel::propertiesName, fvPatchField< Type >::snGrad(), Foam::sqrt(), U, y, yPlus, and Foam::Zero.
Referenced by nutkWallFunctionFvPatchScalarField::calcNut(), nutkRoughWallFunctionFvPatchScalarField::calcNut(), and atmNutkWallFunctionFvPatchScalarField::calcNut().