This boundary condition provides a turbulent kinematic viscosity condition when using wall functions. 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. More...


Public Member Functions | |
| TypeName ("nutTabulatedWallFunction") | |
| Runtime type information. More... | |
| nutUTabulatedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
| Construct from patch and internal field. More... | |
| nutUTabulatedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
| Construct from patch, internal field and dictionary. More... | |
| nutUTabulatedWallFunctionFvPatchScalarField (const nutUTabulatedWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
| Construct by mapping given. More... | |
| nutUTabulatedWallFunctionFvPatchScalarField (const nutUTabulatedWallFunctionFvPatchScalarField &) | |
| Construct as copy. More... | |
| virtual tmp< fvPatchScalarField > | clone () const |
| Construct and return a clone. More... | |
| nutUTabulatedWallFunctionFvPatchScalarField (const nutUTabulatedWallFunctionFvPatchScalarField &, 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... | |
| 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. More... | |
| nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &) | |
| Construct as copy. More... | |
| nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
| Construct as copy setting internal field reference. More... | |
| scalar | yPlusLam () const |
| Return the Y+ at the edge of the laminar sublayer. More... | |
| virtual void | updateCoeffs () |
| Update the coefficients associated with the patch field. More... | |
Protected Member Functions | |
| virtual tmp< scalarField > | calcNut () const |
| Calculate the turbulence viscosity. More... | |
| virtual tmp< scalarField > | calcUPlus (const scalarField &Rey) const |
| Calculate wall u+ from table. More... | |
Protected Member Functions inherited from nutWallFunctionFvPatchScalarField | |
| virtual void | checkType () |
| Check the type of the patch. More... | |
| virtual void | writeLocalEntries (Ostream &) const |
| Write local wall function variables. More... | |
Protected Attributes | |
| word | uPlusTableName_ |
| Name of u+ table. More... | |
| uniformInterpolationTable< scalar > | uPlusTable_ |
| U+ table. More... | |
Protected Attributes inherited from nutWallFunctionFvPatchScalarField | |
| scalar | Cmu_ |
| Cmu coefficient. More... | |
| scalar | kappa_ |
| Von Karman constant. More... | |
| scalar | E_ |
| E coefficient. More... | |
| scalar | yPlusLam_ |
| Y+ at the edge of the laminar sublayer. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from nutWallFunctionFvPatchScalarField | |
| static scalar | yPlusLam (const scalar kappa, const scalar E) |
| Calculate the Y+ at the edge of the laminar sublayer. More... | |
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions. 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.
| Property | Description | Required | Default value |
|---|---|---|---|
uPlusTable | U+ as a function of Re table name | yes |
Example of the boundary condition specification:
<patchName>
{
type nutTabulatedWallFunction;
uPlusTable myUPlusTable;
}Definition at line 87 of file nutUTabulatedWallFunctionFvPatchScalarField.H.
| nutUTabulatedWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF | ||
| ) |
Construct from patch and internal field.
Definition at line 90 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
Referenced by nutUTabulatedWallFunctionFvPatchScalarField::calcUPlus(), nutUTabulatedWallFunctionFvPatchScalarField::clone(), and nutUTabulatedWallFunctionFvPatchScalarField::nutUTabulatedWallFunctionFvPatchScalarField().

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

| 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 115 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References nutUTabulatedWallFunctionFvPatchScalarField::nutUTabulatedWallFunctionFvPatchScalarField().

| nutUTabulatedWallFunctionFvPatchScalarField | ( | const nutUTabulatedWallFunctionFvPatchScalarField & | wfpsf | ) |
Construct as copy.
Definition at line 156 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References nutUTabulatedWallFunctionFvPatchScalarField::nutUTabulatedWallFunctionFvPatchScalarField().

| nutUTabulatedWallFunctionFvPatchScalarField | ( | const nutUTabulatedWallFunctionFvPatchScalarField & | wfpsf, |
| const DimensionedField< scalar, volMesh > & | iF | ||
| ) |
Construct as copy setting internal field reference.
Definition at line 168 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Calculate the turbulence viscosity.
Implements nutWallFunctionFvPatchScalarField.
Definition at line 39 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), nutUTabulatedWallFunctionFvPatchScalarField::calcUPlus(), Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), magUp, Foam::max(), turbulenceModel::nu(), patchi, fvPatchField< Type >::patchInternalField(), turbulenceModel::propertiesName, fvPatchField< Type >::snGrad(), Foam::sqr(), turbulenceModel::U(), y, and turbulenceModel::y().

|
protectedvirtual |
Calculate wall u+ from table.
Definition at line 70 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References forAll, uniformInterpolationTable< Type >::interpolateLog10(), nutUTabulatedWallFunctionFvPatchScalarField::nutUTabulatedWallFunctionFvPatchScalarField(), tmp< T >::ref(), uPlus, and nutUTabulatedWallFunctionFvPatchScalarField::uPlusTable_.
Referenced by nutUTabulatedWallFunctionFvPatchScalarField::calcNut(), and nutUTabulatedWallFunctionFvPatchScalarField::yPlus().


| TypeName | ( | "nutTabulatedWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Definition at line 152 of file nutUTabulatedWallFunctionFvPatchScalarField.H.
References nutUTabulatedWallFunctionFvPatchScalarField::nutUTabulatedWallFunctionFvPatchScalarField().

|
inlinevirtual |
Construct and return a clone setting internal field reference.
Definition at line 169 of file nutUTabulatedWallFunctionFvPatchScalarField.H.
References nutUTabulatedWallFunctionFvPatchScalarField::nutUTabulatedWallFunctionFvPatchScalarField(), nutUTabulatedWallFunctionFvPatchScalarField::write(), and nutUTabulatedWallFunctionFvPatchScalarField::yPlus().

|
virtual |
Calculate and return the yPlus at the boundary.
Implements nutWallFunctionFvPatchScalarField.
Definition at line 181 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), nutUTabulatedWallFunctionFvPatchScalarField::calcUPlus(), Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), magUp, turbulenceModel::nu(), patchi, fvPatchField< Type >::patchInternalField(), turbulenceModel::propertiesName, Rey, turbulenceModel::U(), y, and turbulenceModel::y().
Referenced by nutUTabulatedWallFunctionFvPatchScalarField::clone().


|
virtual |
Write.
Reimplemented from nutWallFunctionFvPatchScalarField.
Definition at line 204 of file nutUTabulatedWallFunctionFvPatchScalarField.C.
References token::END_STATEMENT, Foam::makePatchTypeField(), Foam::nl, nutUTabulatedWallFunctionFvPatchScalarField::uPlusTableName_, fvPatchField< Type >::write(), and Ostream::writeKeyword().
Referenced by nutUTabulatedWallFunctionFvPatchScalarField::clone().


|
protected |
Name of u+ table.
Definition at line 96 of file nutUTabulatedWallFunctionFvPatchScalarField.H.
Referenced by nutUTabulatedWallFunctionFvPatchScalarField::write().
|
protected |
U+ table.
Definition at line 99 of file nutUTabulatedWallFunctionFvPatchScalarField.H.
Referenced by nutUTabulatedWallFunctionFvPatchScalarField::calcUPlus().
1.8.13