Loading...
Searching...
No Matches
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 applications. More...

#include <nutUTabulatedWallFunctionFvPatchScalarField.H>

Inheritance diagram for nutUTabulatedWallFunctionFvPatchScalarField:
Collaboration diagram for nutUTabulatedWallFunctionFvPatchScalarField:

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< scalarFieldyPlus () 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 wallFunctionCoefficientswallCoeffs () const noexcept
 Return wallFunctionCoefficients.
virtual void updateCoeffs ()
 Update the coefficients associated with the patch field.

Protected Member Functions

virtual tmp< scalarFieldcalcNut () const
 Calculate the turbulent viscosity.
tmp< scalarFieldcalcUPlus (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 volVectorFieldU (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 nutWallFunctionFvPatchScalarFieldnutw (const turbulenceModel &turbModel, const label patchi)
 Return the nut patchField for the given wall patch.

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

Usage
Example of the boundary condition specification:
<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:

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

Definition at line 105 of file nutUTabulatedWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ nutUTabulatedWallFunctionFvPatchScalarField() [1/5]

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().

Here is the call graph for this function:
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 133 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

References dict, mesh, nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField(), p, uPlusTable_, and uPlusTableName_.

Here is the call graph for this function:

◆ nutUTabulatedWallFunctionFvPatchScalarField() [3/5]

nutUTabulatedWallFunctionFvPatchScalarField ( const nutUTabulatedWallFunctionFvPatchScalarField & ptf,
const fvPatch & p,
const DimensionedField< scalar, volMesh > & iF,
const fvPatchFieldMapper & mapper )

◆ nutUTabulatedWallFunctionFvPatchScalarField() [4/5]

nutUTabulatedWallFunctionFvPatchScalarField ( const nutUTabulatedWallFunctionFvPatchScalarField & wfpsf)

◆ nutUTabulatedWallFunctionFvPatchScalarField() [5/5]

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_.

Here is the call graph for this function:

Member Function Documentation

◆ calcNut()

Foam::tmp< Foam::scalarField > calcNut ( ) const
protectedvirtual

◆ calcUPlus()

Foam::tmp< Foam::scalarField > calcUPlus ( const scalarField & Rey) const
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().

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

◆ writeLocalEntries()

void writeLocalEntries ( Ostream & os) const
protected

Write local wall function variables.

Definition at line 82 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

References os(), and uPlusTableName_.

Referenced by write().

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

◆ TypeName()

TypeName ( "nutTabulatedWallFunction" )

Runtime type information.

References nutUTabulatedWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ clone() [1/2]

virtual tmp< fvPatchField< scalar > > clone ( ) const
inlinevirtual

Return a clone.

Definition at line 204 of file nutUTabulatedWallFunctionFvPatchScalarField.H.

References fvPatchField< Type >::Clone().

Here is the call graph for this function:

◆ clone() [2/2]

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

Clone with an internal field reference.

Definition at line 212 of file nutUTabulatedWallFunctionFvPatchScalarField.H.

References fvPatchField< Type >::Clone().

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 187 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

References calcUPlus(), IOobject::groupName(), Foam::mag(), magUp, fvPatchField< Type >::patchInternalField(), turbulenceModel::propertiesName, Rey, U, and y.

Here is the call graph for this function:

◆ write()

void write ( Ostream & os) const
virtual

Write.

Reimplemented from nutWallFunctionFvPatchScalarField.

Definition at line 214 of file nutUTabulatedWallFunctionFvPatchScalarField.C.

References os(), nutWallFunctionFvPatchScalarField::write(), writeLocalEntries(), and fvPatchField< Type >::writeValueEntry().

Here is the call graph for this function:

Member Data Documentation

◆ uPlusTableName_

◆ uPlusTable_


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