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...
#include <nutUBlendedWallFunctionFvPatchScalarField.H>


Public Member Functions | |
| TypeName ("nutUBlendedWallFunction") | |
| Runtime type information. | |
| nutUBlendedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
| Construct from patch and internal field. | |
| nutUBlendedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
| Construct from patch, internal field and dictionary. | |
| nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
| Construct by mapping given nutUBlendedWallFunctionFvPatchScalarField onto a new patch. | |
| nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &) | |
| Construct as copy. | |
| nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
| Construct as copy setting internal field reference. | |
| virtual tmp< fvPatchField< scalar > > | clone () const |
| Return a clone. | |
| virtual tmp< fvPatchField< scalar > > | clone (const DimensionedField< scalar, volMesh > &iF) const |
| Clone with an internal field reference. | |
| virtual tmp< scalarField > | yPlus () const |
| Calculate and return the yPlus at the boundary. | |
| virtual void | write (Ostream &os) const |
| Write. | |
| Public Member Functions inherited from nutWallFunctionFvPatchScalarField | |
| TypeName ("nutWallFunction") | |
| Runtime type information. | |
| nutWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
| Construct from patch and internal field. | |
| nutWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
| Construct from patch, internal field and dictionary. | |
| nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
| Construct by mapping given nutWallFunctionFvPatchScalarField onto a new patch. | |
| nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &) | |
| Construct as copy. | |
| nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
| Construct as copy setting internal field reference. | |
| const wallFunctionCoefficients & | wallCoeffs () const noexcept |
| Return wallFunctionCoefficients. | |
| virtual void | updateCoeffs () |
| Update the coefficients associated with the patch field. | |
Protected Member Functions | |
| virtual tmp< scalarField > | calcNut () const |
| Calculate the turbulent viscosity. | |
| tmp< scalarField > | calcUTau (const scalarField &magGradU) const |
| Calculate the friction velocity. | |
| void | writeLocalEntries (Ostream &) const |
| Write local wall function variables. | |
| 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. | |
| virtual void | checkType () |
| Check the type of the patch. | |
| void | writeLocalEntries (Ostream &) const |
| Write local wall function variables. | |
Protected Attributes | |
| scalar | n_ |
| Model coefficient; default = 4. | |
| Protected Attributes inherited from nutWallFunctionFvPatchScalarField | |
| word | UName_ |
| Name of velocity field. | |
| wallFunctionCoefficients | wallCoeffs_ |
| Wall-function coefficients. | |
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. | |
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.
![\[ u_\tau = (u_{\tau,vis}^n + u_{\tau,log}^n)^{1/n}
\]](form_1038.png)
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 131 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References n_, nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField(), and p.
Referenced by nutUBlendedWallFunctionFvPatchScalarField(), nutUBlendedWallFunctionFvPatchScalarField(), nutUBlendedWallFunctionFvPatchScalarField(), and TypeName().


| nutUBlendedWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF, | ||
| const dictionary & | dict ) |
Construct from patch, internal field and dictionary.
Definition at line 157 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References dict, n_, nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField(), and p.

| 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 143 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References n_, nutUBlendedWallFunctionFvPatchScalarField(), nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField(), and p.

| nutUBlendedWallFunctionFvPatchScalarField | ( | const nutUBlendedWallFunctionFvPatchScalarField & | wfpsf | ) |
Construct as copy.
Definition at line 170 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References n_, nutUBlendedWallFunctionFvPatchScalarField(), and nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField().

| nutUBlendedWallFunctionFvPatchScalarField | ( | const nutUBlendedWallFunctionFvPatchScalarField & | wfpsf, |
| const DimensionedField< scalar, volMesh > & | iF ) |
Construct as copy setting internal field reference.
Definition at line 181 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References n_, nutUBlendedWallFunctionFvPatchScalarField(), and nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField().

|
protectedvirtual |
Calculate the turbulent viscosity.
Implements nutWallFunctionFvPatchScalarField.
Definition at line 30 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References calcUTau(), IOobject::groupName(), Foam::mag(), Foam::max(), turbulenceModel::propertiesName, Foam::sqr(), and U.

|
protected |
Calculate the friction velocity.
Definition at line 58 of file nutUBlendedWallFunctionFvPatchScalarField.C.
References e, forAll, IOobject::groupName(), Foam::log(), Foam::mag(), magUp, Foam::max(), n, n_, tmp< T >::New(), nutWallFunctionFvPatchScalarField::nutw(), fvPatchField< Type >::patchInternalField(), Foam::pow(), turbulenceModel::propertiesName, Foam::sqrt(), U, nutWallFunctionFvPatchScalarField::wallCoeffs_, y, yPlus(), and Foam::Zero.
Referenced by calcNut(), and yPlus().


|
protected |
Write local wall function variables.
Definition at line 120 of file nutUBlendedWallFunctionFvPatchScalarField.C.
Referenced by write().


| TypeName | ( | "nutUBlendedWallFunction" | ) |
Runtime type information.
References nutUBlendedWallFunctionFvPatchScalarField().

|
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 241 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 calcUTau(), IOobject::groupName(), Foam::mag(), turbulenceModel::propertiesName, fvPatchField< Type >::snGrad(), U, and y.
Referenced by calcUTau().


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

|
protected |
Model coefficient; default = 4.
Definition at line 150 of file nutUBlendedWallFunctionFvPatchScalarField.H.
Referenced by calcUTau(), nutUBlendedWallFunctionFvPatchScalarField(), nutUBlendedWallFunctionFvPatchScalarField(), nutUBlendedWallFunctionFvPatchScalarField(), nutUBlendedWallFunctionFvPatchScalarField(), nutUBlendedWallFunctionFvPatchScalarField(), and writeLocalEntries().