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

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

Public Member Functions

 TypeName ("nutUBlendedWallFunction")
 Runtime type information. More...
 
 nutUBlendedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Construct from patch and internal field. More...
 
 nutUBlendedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 Construct from patch, internal field and dictionary. More...
 
 nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 
 nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &)
 Construct as copy. More...
 
virtual tmp< fvPatchScalarFieldclone () const
 Construct and return a clone. More...
 
 nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &, 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...
 
virtual tmp< scalarFieldyPlus () const =0
 Calculate and return the yPlus at the boundary. More...
 
const wallFunctionCoefficientswallCoeffs () const noexcept
 Return wallFunctionCoefficients. More...
 
virtual void updateCoeffs ()
 Update the coefficients associated with the patch field. More...
 
virtual void write (Ostream &) const
 Write. More...
 

Protected Member Functions

virtual tmp< scalarFieldcalcNut () const
 Calculate the turbulent viscosity. More...
 
tmp< scalarFieldcalcUTau (const scalarField &magGradU) const
 Calculate the friction velocity. More...
 
void writeLocalEntries (Ostream &) const
 Write local wall function variables. 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 tmp< scalarFieldcalcNut () const =0
 Calculate the turbulent viscosity. More...
 
void writeLocalEntries (Ostream &) const
 Write local wall function variables. More...
 

Protected Attributes

scalar n_
 Model coefficient; default = 4. More...
 
- Protected Attributes inherited from nutWallFunctionFvPatchScalarField
word UName_
 Name of velocity field. More...
 
wallFunctionCoefficients wallCoeffs_
 Wall-function coefficients. More...
 

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

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 143 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 138 of file nutUBlendedWallFunctionFvPatchScalarField.C.

◆ 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 164 of file nutUBlendedWallFunctionFvPatchScalarField.C.

◆ nutUBlendedWallFunctionFvPatchScalarField() [3/5]

Construct by mapping given nutUBlendedWallFunctionFvPatchScalarField onto a new patch

Definition at line 150 of file nutUBlendedWallFunctionFvPatchScalarField.C.

◆ nutUBlendedWallFunctionFvPatchScalarField() [4/5]

◆ nutUBlendedWallFunctionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 188 of file nutUBlendedWallFunctionFvPatchScalarField.C.

Member Function Documentation

◆ calcNut()

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

Calculate the turbulent viscosity.

Implements nutWallFunctionFvPatchScalarField.

Definition at line 37 of file nutUBlendedWallFunctionFvPatchScalarField.C.

References nutUBlendedWallFunctionFvPatchScalarField::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

Calculate the friction velocity.

Definition at line 65 of file nutUBlendedWallFunctionFvPatchScalarField.C.

References e, forAll, IOobject::groupName(), Foam::log(), Foam::mag(), magUp, Foam::max(), n, Time::New(), fvPatchField< Type >::patchInternalField(), Foam::pow(), turbulenceModel::propertiesName, Foam::sqrt(), U, y, yPlus, and Foam::Zero.

Referenced by nutUBlendedWallFunctionFvPatchScalarField::calcNut().

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 127 of file nutUBlendedWallFunctionFvPatchScalarField.C.

References os(), and Ostream::writeEntryIfDifferent().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "nutUBlendedWallFunction"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp< fvPatchScalarField > clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 208 of file nutUBlendedWallFunctionFvPatchScalarField.H.

◆ 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 224 of file nutUBlendedWallFunctionFvPatchScalarField.H.

◆ yPlus()

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

Calculate and return the yPlus at the boundary.

Implements nutWallFunctionFvPatchScalarField.

Definition at line 203 of file nutUBlendedWallFunctionFvPatchScalarField.C.

References IOobject::groupName(), Foam::mag(), turbulenceModel::propertiesName, fvPatchField< Type >::snGrad(), 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 228 of file nutUBlendedWallFunctionFvPatchScalarField.C.

References os(), and ObukhovLength::write().

Here is the call graph for this function:

Member Data Documentation

◆ n_

scalar n_
protected

Model coefficient; default = 4.

Definition at line 152 of file nutUBlendedWallFunctionFvPatchScalarField.H.


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