This boundary condition provides a wall constraint on the turbulent viscosity (i.e. nut) based on velocity (i.e. U), for low- and high-Reynolds number applications.
More...
#include <nutUTabulatedWallFunctionFvPatchScalarField.H>


Public Member Functions | |
| TypeName ("nutTabulatedWallFunction") | |
| Runtime type information. | |
| nutUTabulatedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
| Construct from patch and internal field. | |
| nutUTabulatedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
| Construct from patch, internal field and dictionary. | |
| nutUTabulatedWallFunctionFvPatchScalarField (const nutUTabulatedWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
| Construct by mapping given nutUTabulatedWallFunctionFvPatchScalarField onto a new patch. | |
| nutUTabulatedWallFunctionFvPatchScalarField (const nutUTabulatedWallFunctionFvPatchScalarField &) | |
| Construct as copy. | |
| nutUTabulatedWallFunctionFvPatchScalarField (const nutUTabulatedWallFunctionFvPatchScalarField &, 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 > | calcUPlus (const scalarField &Rey) const |
| Calculate wall u+ from table. | |
| 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 | |
| word | uPlusTableName_ |
| Name of u+ table. | |
| uniformInterpolationTable< scalar > | uPlusTable_ |
| u+ table | |
| 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 constraint on the turbulent viscosity (i.e. nut) based on velocity (i.e. U), for low- and high-Reynolds number applications.
As input, the user specifies a look-up table of u+ as a function of near-wall Reynolds number.
The table should be located in the $FOAM_CASE/constant directory.
<patchName>
{
// Mandatory entries
type nutTabulatedWallFunction;
uPlusTable myUPlusTable;
// Inherited entries
...
}
where the entries mean:
| Property | Description | Type | Reqd | Deflt |
|---|---|---|---|---|
type | Type name: nutUTabulatedWallFunction | word | yes | - |
uPlusTable | u+ as a function of Re table name | word | yes | - |
The inherited entries are elaborated in:
Definition at line 105 of file nutUTabulatedWallFunctionFvPatchScalarField.H.
| nutUTabulatedWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF ) |
Construct from patch and internal field.
Definition at line 93 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References mesh, nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField(), p, uPlusTable_, and uPlusTableName_.
Referenced by nutUTabulatedWallFunctionFvPatchScalarField(), nutUTabulatedWallFunctionFvPatchScalarField(), nutUTabulatedWallFunctionFvPatchScalarField(), and TypeName().


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

| nutUTabulatedWallFunctionFvPatchScalarField | ( | const nutUTabulatedWallFunctionFvPatchScalarField & | ptf, |
| const fvPatch & | p, | ||
| const DimensionedField< scalar, volMesh > & | iF, | ||
| const fvPatchFieldMapper & | mapper ) |
Construct by mapping given nutUTabulatedWallFunctionFvPatchScalarField onto a new patch.
Definition at line 118 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References nutUTabulatedWallFunctionFvPatchScalarField(), nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField(), p, uPlusTable_, and uPlusTableName_.

| nutUTabulatedWallFunctionFvPatchScalarField | ( | const nutUTabulatedWallFunctionFvPatchScalarField & | wfpsf | ) |
Construct as copy.
Definition at line 159 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References nutUTabulatedWallFunctionFvPatchScalarField(), nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField(), uPlusTable_, and uPlusTableName_.

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

|
protectedvirtual |
Calculate the turbulent viscosity.
Implements nutWallFunctionFvPatchScalarField.
Definition at line 31 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References calcUPlus(), IOobject::groupName(), Foam::mag(), magUp, Foam::max(), fvPatchField< Type >::patchInternalField(), turbulenceModel::propertiesName, fvPatchField< Type >::snGrad(), Foam::sqr(), U, and y.

|
protected |
Calculate wall u+ from table.
Definition at line 65 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References forAll, tmp< T >::New(), Rey, uPlus, uPlusTable_, and Foam::Zero.
Referenced by calcNut(), and yPlus().


|
protected |
Write local wall function variables.
Definition at line 82 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References os(), and uPlusTableName_.
Referenced by write().


| TypeName | ( | "nutTabulatedWallFunction" | ) |
Runtime type information.
References nutUTabulatedWallFunctionFvPatchScalarField().

|
inlinevirtual |
Return a clone.
Definition at line 204 of file nutUTabulatedWallFunctionFvPatchScalarField.H.
References fvPatchField< Type >::Clone().

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

|
virtual |
Calculate and return the yPlus at the boundary.
Implements nutWallFunctionFvPatchScalarField.
Definition at line 187 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References calcUPlus(), IOobject::groupName(), Foam::mag(), magUp, fvPatchField< Type >::patchInternalField(), turbulenceModel::propertiesName, Rey, U, and y.

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

|
protected |
Name of u+ table.
Definition at line 116 of file nutUTabulatedWallFunctionFvPatchScalarField.H.
Referenced by nutUTabulatedWallFunctionFvPatchScalarField(), nutUTabulatedWallFunctionFvPatchScalarField(), nutUTabulatedWallFunctionFvPatchScalarField(), nutUTabulatedWallFunctionFvPatchScalarField(), nutUTabulatedWallFunctionFvPatchScalarField(), and writeLocalEntries().
|
protected |
u+ table
Definition at line 121 of file nutUTabulatedWallFunctionFvPatchScalarField.H.
Referenced by calcUPlus(), nutUTabulatedWallFunctionFvPatchScalarField(), nutUTabulatedWallFunctionFvPatchScalarField(), nutUTabulatedWallFunctionFvPatchScalarField(), nutUTabulatedWallFunctionFvPatchScalarField(), and nutUTabulatedWallFunctionFvPatchScalarField().