Loading...
Searching...
No Matches
nutUBlendedWallFunctionFvPatchScalarField Class Reference

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>

Inheritance diagram for nutUBlendedWallFunctionFvPatchScalarField:
Collaboration diagram for nutUBlendedWallFunctionFvPatchScalarField:

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< 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< scalarFieldcalcUTau (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 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

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 nutWallFunctionFvPatchScalarFieldnutw (const turbulenceModel &turbModel, const label patchi)
 Return the nut patchField for the given wall patch.

Detailed Description

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}
\‍]

where

$      u_\tau     $=Friction velocity
$      u_{\tau,vis} $=Friction velocity in the viscous sublayer
$      u_{\tau,log} $=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.
Usage
Example of the boundary condition specification:
<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:

Note
Source files

Definition at line 139 of file nutUBlendedWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ nutUBlendedWallFunctionFvPatchScalarField() [1/5]

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

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

◆ nutUBlendedWallFunctionFvPatchScalarField() [2/5]

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.

Here is the call graph for this function:

◆ nutUBlendedWallFunctionFvPatchScalarField() [3/5]

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.

Here is the call graph for this function:

◆ nutUBlendedWallFunctionFvPatchScalarField() [4/5]

nutUBlendedWallFunctionFvPatchScalarField ( const nutUBlendedWallFunctionFvPatchScalarField & wfpsf)

Construct as copy.

Definition at line 170 of file nutUBlendedWallFunctionFvPatchScalarField.C.

References n_, nutUBlendedWallFunctionFvPatchScalarField(), and nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ nutUBlendedWallFunctionFvPatchScalarField() [5/5]

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

Here is the call graph for this function:

Member Function Documentation

◆ calcNut()

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

Here is the call graph for this function:

◆ calcUTau()

Foam::tmp< Foam::scalarField > calcUTau ( const scalarField & magGradU) const
protected

◆ writeLocalEntries()

void writeLocalEntries ( Ostream & os) const
protected

Write local wall function variables.

Definition at line 120 of file nutUBlendedWallFunctionFvPatchScalarField.C.

References n_, and os().

Referenced by write().

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

◆ TypeName()

TypeName ( "nutUBlendedWallFunction" )

Runtime type information.

References nutUBlendedWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ clone() [1/2]

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

Return a clone.

Definition at line 233 of file nutUBlendedWallFunctionFvPatchScalarField.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 241 of file nutUBlendedWallFunctionFvPatchScalarField.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 196 of file nutUBlendedWallFunctionFvPatchScalarField.C.

References calcUTau(), IOobject::groupName(), Foam::mag(), turbulenceModel::propertiesName, fvPatchField< Type >::snGrad(), U, and y.

Referenced by calcUTau().

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

◆ write()

void write ( Ostream & os) const
virtual

Write.

Reimplemented from nutWallFunctionFvPatchScalarField.

Definition at line 221 of file nutUBlendedWallFunctionFvPatchScalarField.C.

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

Here is the call graph for this function:

Member Data Documentation

◆ n_


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