nutUTabulatedWallFunctionFvPatchScalarField Class Reference

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 turbulence models. More...

Inheritance diagram for nutUTabulatedWallFunctionFvPatchScalarField:
[legend]
Collaboration diagram for nutUTabulatedWallFunctionFvPatchScalarField:
[legend]

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 &)
 
 nutUTabulatedWallFunctionFvPatchScalarField (const nutUTabulatedWallFunctionFvPatchScalarField &)
 Construct as copy. More...
 
virtual tmp< fvPatchScalarFieldclone () const
 Construct and return a clone. More...
 
 nutUTabulatedWallFunctionFvPatchScalarField (const nutUTabulatedWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 Construct as copy setting internal field reference. More...
 
virtual tmp< fvPatchScalarFieldclone (const DimensionedField< scalar, volMesh > &iF) const
 Construct and return a clone setting internal field reference. More...
 
virtual tmp< scalarFieldyPlus () 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 &)
 
 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...
 

Protected Member Functions

virtual tmp< scalarFieldcalcNut () const
 Calculate the turbulent viscosity. More...
 
virtual tmp< scalarFieldcalcUPlus (const scalarField &Rey) const
 Calculate wall u+ from table. More...
 
- Protected Member Functions inherited from nutWallFunctionFvPatchScalarField
virtual const volVectorFieldU (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...
 

Protected Attributes

word uPlusTableName_
 Name of u+ table. More...
 
uniformInterpolationTable< scalar > uPlusTable_
 u+ table More...
 
- Protected Attributes inherited from nutWallFunctionFvPatchScalarField
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_
 

Additional Inherited Members

- Static Public Member Functions inherited from nutWallFunctionFvPatchScalarField
static const nutWallFunctionFvPatchScalarFieldnutw (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...
 
- Protected Types inherited from nutWallFunctionFvPatchScalarField
enum  blendingType { STEPWISE, MAX, BINOMIAL, EXPONENTIAL }
 Options for the blending treatment of viscous and inertial sublayers. More...
 
- Static Protected Attributes inherited from nutWallFunctionFvPatchScalarField
static const Enum< blendingTypeblendingTypeNames
 Names for blendingType. More...
 

Detailed Description

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 turbulence models.

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.

Usage
Example of the boundary condition specification:
<patchName>
{
    // Mandatory entries (unmodifiable)
    type            nutTabulatedWallFunction;
    uPlusTable      myUPlusTable;

    // Optional (inherited) entries
    ...
}

where the entries mean:

Property Description Type Req'd Dflt
type Type name: nutUTabulatedWallFunction word yes -
uPlusTable u+ as a function of Re table name word yes -

The inherited entries are elaborated in:

Note
  • The tables are not registered since the same table object may be used for more than one patch.
Source files

Definition at line 108 of file nutUTabulatedWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ nutUTabulatedWallFunctionFvPatchScalarField() [1/5]

Construct from patch and internal field.

Definition at line 91 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

Referenced by nutUTabulatedWallFunctionFvPatchScalarField::clone().

Here is the caller graph for this function:

◆ nutUTabulatedWallFunctionFvPatchScalarField() [2/5]

nutUTabulatedWallFunctionFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF,
const dictionary dict 
)

Construct from patch, internal field and dictionary.

Definition at line 131 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

◆ nutUTabulatedWallFunctionFvPatchScalarField() [3/5]

Construct by mapping given nutUTabulatedWallFunctionFvPatchScalarField onto a new patch

Definition at line 116 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

◆ nutUTabulatedWallFunctionFvPatchScalarField() [4/5]

◆ nutUTabulatedWallFunctionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 169 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

Member Function Documentation

◆ calcNut()

◆ calcUPlus()

Foam::tmp< Foam::scalarField > calcUPlus ( const scalarField Rey) const
protectedvirtual

Calculate wall u+ from table.

Definition at line 71 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

References forAll, Foam::foamVersion::patch, tmp< T >::ref(), Rey, uPlus, and Foam::Zero.

Referenced by nutUTabulatedWallFunctionFvPatchScalarField::calcNut().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ TypeName()

TypeName ( "nutTabulatedWallFunction"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchScalarField> clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 173 of file nutUTabulatedWallFunctionFvPatchScalarField.H.

References nutUTabulatedWallFunctionFvPatchScalarField::nutUTabulatedWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ clone() [2/2]

virtual tmp<fvPatchScalarField> clone ( const DimensionedField< scalar, volMesh > &  iF) const
inlinevirtual

Construct and return a clone setting internal field reference.

Definition at line 190 of file nutUTabulatedWallFunctionFvPatchScalarField.H.

References nutUTabulatedWallFunctionFvPatchScalarField::nutUTabulatedWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ yPlus()

Foam::tmp< Foam::scalarField > yPlus ( ) const
virtual

Calculate and return the yPlus at the boundary.

Implements nutWallFunctionFvPatchScalarField.

Definition at line 183 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

References Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), magUp, turbulenceModel::nu(), Foam::foamVersion::patch, fvPatchField< Type >::patchInternalField(), turbulenceModel::propertiesName, Rey, U, y, and turbulenceModel::y().

Here is the call graph for this function:

◆ write()

void write ( Ostream os) const
virtual

Write.

Reimplemented from nutWallFunctionFvPatchScalarField.

Definition at line 207 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

References os(), fvPatchField< Type >::write(), and Ostream::writeEntry().

Here is the call graph for this function:

Member Data Documentation

◆ uPlusTableName_

word uPlusTableName_
protected

Name of u+ table.

Definition at line 117 of file nutUTabulatedWallFunctionFvPatchScalarField.H.

◆ uPlusTable_

uniformInterpolationTable<scalar> uPlusTable_
protected

u+ table

Definition at line 120 of file nutUTabulatedWallFunctionFvPatchScalarField.H.


The documentation for this class was generated from the following files: