This boundary condition provides a wall function for the turbulent viscosity (i.e. nut) based on velocity (i.e. U) using a binomial-function wall-function blending method between the viscous and inertial sublayer predictions of nut for low- and high-Reynolds number applications.
More...


Public Member Functions | |
| TypeName ("nutUBlendedWallFunction") | |
| Runtime type information. More... | |
| nutUBlendedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
| Construct from patch and internal field. More... | |
| nutUBlendedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
| Construct from patch, internal field and dictionary. More... | |
| nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
| Construct by mapping given nutUBlendedWallFunctionFvPatchScalarField onto a new patch. More... | |
| nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &) | |
| Construct as copy. More... | |
| nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
| Construct as copy setting internal field reference. More... | |
| virtual tmp< fvPatchField< scalar > > | clone () const |
| Return a clone. More... | |
| virtual tmp< fvPatchField< scalar > > | clone (const DimensionedField< scalar, volMesh > &iF) const |
| Clone with an internal field reference. More... | |
| virtual tmp< scalarField > | yPlus () const |
| Calculate and return the yPlus at the boundary. More... | |
| virtual void | write (Ostream &os) const |
| Write. 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 &) | |
| Construct by mapping given nutWallFunctionFvPatchScalarField onto a new patch. More... | |
| nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &) | |
| Construct as copy. More... | |
| nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
| Construct as copy setting internal field reference. More... | |
| const wallFunctionCoefficients & | wallCoeffs () const noexcept |
| Return wallFunctionCoefficients. More... | |
| virtual void | updateCoeffs () |
| Update the coefficients associated with the patch field. More... | |
Protected Member Functions | |
| virtual tmp< scalarField > | calcNut () const |
| Calculate the turbulent viscosity. More... | |
| tmp< scalarField > | calcUTau (const scalarField &magGradU) const |
| Calculate the friction velocity. More... | |
| void | writeLocalEntries (Ostream &) const |
| Write local wall function variables. More... | |
Protected Member Functions inherited from nutWallFunctionFvPatchScalarField | |
| virtual const volVectorField & | U (const turbulenceModel &turb) const |
| Helper to return the velocity field either from the turbulence model (default) or the mesh database. More... | |
| virtual void | checkType () |
| Check the type of the patch. More... | |
| void | writeLocalEntries (Ostream &) const |
| Write local wall function variables. More... | |
Protected Attributes | |
| scalar | n_ |
| Model coefficient; default = 4. More... | |
Protected Attributes inherited from nutWallFunctionFvPatchScalarField | |
| word | UName_ |
| Name of velocity field. More... | |
| wallFunctionCoefficients | wallCoeffs_ |
| Wall-function coefficients. 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... | |
This boundary condition provides a wall function for the turbulent viscosity (i.e. nut) based on velocity (i.e. U) using a binomial-function wall-function blending method between the viscous and inertial sublayer predictions of nut for low- and high-Reynolds number applications.
where
| = | Friction velocity |
| = | Friction velocity in the viscous sublayer |
| = | Friction velocity in the inertial sublayer |
Reference:
See the section that describes 'automatic wall treatment':
Menter, F., Ferreira, J. C., Esch, T., Konno, B. (2003).
The SST turbulence model with improved wall treatment
for heat transfer predictions in gas turbines.
In Proceedings of the International Gas Turbine Congress.
November, 2003. Tokyo, Japan. pp. 2-7.<patchName>
{
// Mandatory entries
type nutUBlendedWallFunction;
// Optional entries
n 4.0;
// Inherited entries
...
}
where the entries mean:
| Property | Description | Type | Reqd | Deflt |
|---|---|---|---|---|
type | Type name: nutUBlendedWallFunction | word | yes | - |
n | Blending factor | scalar | no | 4.0 |
The inherited entries are elaborated in:
blending option binomial or with the deprecated blended flag set to on.correctNut() (called through turbulence->validate) returns a slightly different value every time it is called. See nutUSpaldingWallFunctionFvPatchScalarField.C.Definition at line 139 of file nutUBlendedWallFunctionFvPatchScalarField.H.
| nutUBlendedWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF | ||
| ) |
Construct from patch and internal field.
Definition at line 133 of file nutUBlendedWallFunctionFvPatchScalarField.C.
| nutUBlendedWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF, | ||
| const dictionary & | dict | ||
| ) |
Construct from patch, internal field and dictionary.
Definition at line 159 of file nutUBlendedWallFunctionFvPatchScalarField.C.
| nutUBlendedWallFunctionFvPatchScalarField | ( | const nutUBlendedWallFunctionFvPatchScalarField & | ptf, |
| const fvPatch & | p, | ||
| const DimensionedField< scalar, volMesh > & | iF, | ||
| const fvPatchFieldMapper & | mapper | ||
| ) |
Construct by mapping given nutUBlendedWallFunctionFvPatchScalarField onto a new patch.
Definition at line 145 of file nutUBlendedWallFunctionFvPatchScalarField.C.
| nutUBlendedWallFunctionFvPatchScalarField | ( | const nutUBlendedWallFunctionFvPatchScalarField & | wfpsf | ) |
Construct as copy.
Definition at line 172 of file nutUBlendedWallFunctionFvPatchScalarField.C.
| nutUBlendedWallFunctionFvPatchScalarField | ( | const nutUBlendedWallFunctionFvPatchScalarField & | wfpsf, |
| const DimensionedField< scalar, volMesh > & | iF | ||
| ) |
Construct as copy setting internal field reference.
Definition at line 183 of file nutUBlendedWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Calculate the turbulent viscosity.
Implements nutWallFunctionFvPatchScalarField.
Definition at line 30 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), nutUBlendedWallFunctionFvPatchScalarField::calcUTau(), Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), Foam::max(), Foam::foamVersion::patch, turbulenceModel::propertiesName, Foam::sqr(), and nutWallFunctionFvPatchScalarField::U().

|
protected |
Calculate the friction velocity.
Definition at line 59 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References Foam::constant::electromagnetic::e, forAll, Foam::constant::atomic::group, IOobject::groupName(), Foam::constant::electromagnetic::kappa, Foam::log(), Foam::mag(), magUp, Foam::max(), n, tmp< T >::New(), Foam::foamVersion::patch, fvPatchField< Type >::patchInternalField(), Foam::pow(), turbulenceModel::propertiesName, Foam::sqrt(), U, y, yPlus, and Foam::Zero.
Referenced by nutUBlendedWallFunctionFvPatchScalarField::calcNut().


|
protected |
Write local wall function variables.
Definition at line 121 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References os(), and Ostream::writeEntryIfDifferent().

| TypeName | ( | "nutUBlendedWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Return a clone.
Definition at line 233 of file nutUBlendedWallFunctionFvPatchScalarField.H.
References fvPatchField< Type >::Clone().

|
inlinevirtual |
Clone with an internal field reference.
Definition at line 242 of file nutUBlendedWallFunctionFvPatchScalarField.H.
References fvPatchField< Type >::Clone().

|
virtual |
Calculate and return the yPlus at the boundary.
Implements nutWallFunctionFvPatchScalarField.
Definition at line 196 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), Foam::foamVersion::patch, turbulenceModel::propertiesName, fvPatchField< Type >::snGrad(), U, and y.

|
virtual |
Write.
Reimplemented from nutWallFunctionFvPatchScalarField.
Definition at line 222 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References os(), nutWallFunctionFvPatchScalarField::write(), and fvPatchField< Type >::writeValueEntry().

|
protected |
Model coefficient; default = 4.
Definition at line 150 of file nutUBlendedWallFunctionFvPatchScalarField.H.