nutUBlendedWallFunctionFvPatchScalarField Class Reference

This boundary condition provides a wall constraint on the turbulent kinematic viscosity, i.e. nut, when using wall functions based on a blending of laminar sub-layer and log region contributions. 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...
 
scalar Cmu () const
 Return Cmu. More...
 
scalar kappa () const
 Return kappa. More...
 
scalar E () const
 Return E. More...
 
scalar yPlusLam () const
 Return the y+ at the edge of the viscous sublayer. More...
 
virtual void updateCoeffs ()
 Update the coefficients associated with the patch field. More...
 

Protected Member Functions

virtual tmp< scalarFieldcalcNut () const
 Calculate the turbulence viscosity. More...
 
virtual tmp< scalarFieldcalcUTau (const scalarField &magGradU) const
 Calculate the friction velocity. 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

scalar n_
 Model coefficient; default = 4. More...
 
- Protected Attributes inherited from nutWallFunctionFvPatchScalarField
word UName_
 Name of velocity field. More...
 
scalar Cmu_
 Cmu coefficient. More...
 
scalar kappa_
 Von Karman constant. More...
 
scalar E_
 E coefficient. More...
 
scalar yPlusLam_
 Estimated y+ value at the edge of the viscous sublayer. 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...
 
static scalar yPlusLam (const scalar kappa, const scalar E)
 Calculate the y+ at the edge of the viscous sublayer. More...
 

Detailed Description

This boundary condition provides a wall constraint on the turbulent kinematic viscosity, i.e. nut, when using wall functions based on a blending of laminar sub-layer and log region contributions.

\[ u_\tau = (u_{\tau,v}^n + u_{\tau,l}^n)^{1/n} \]

where

\( u_\tau \) = friction velocity
\( u_{\tau,v} \) = friction velocity in the viscous region
\( u_{\tau,l} \) = friction velocity in the log region

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.
        Proceedings of the International Gas Turbine Congress 2003 Tokyo
Usage
Property Description Required Default value
Cmu Model coefficient no 0.09
kappa Von Karman constant no 0.41
E Model coefficient no 9.8

Example of the boundary condition specification:

    <patchName>
    {
        // Mandatory entries
        type            nutUBlendedWallFunction;

        // Optional entries
     }
Note
The full 'automatic wall treatment' description also requires use of the Foam::omegaWallFunction with the blended flag set to 'on'

Suffers from non-exact restart since correctNut() (called through turbulence->validate) returns a slightly different value every time it is called. See nutUSpaldingWallFunctionFvPatchScalarField.C

See also
Foam::nutWallFunctionFvPatchScalarField Foam::omegaWallFunctionFvPatchScalarField
Source files

Definition at line 138 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 125 of file nutUBlendedWallFunctionFvPatchScalarField.C.

Referenced by nutUBlendedWallFunctionFvPatchScalarField::clone().

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

◆ nutUBlendedWallFunctionFvPatchScalarField() [3/5]

Construct by mapping given nutUBlendedWallFunctionFvPatchScalarField onto a new patch

Definition at line 137 of file nutUBlendedWallFunctionFvPatchScalarField.C.

◆ nutUBlendedWallFunctionFvPatchScalarField() [4/5]

◆ nutUBlendedWallFunctionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 175 of file nutUBlendedWallFunctionFvPatchScalarField.C.

Member Function Documentation

◆ calcNut()

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

◆ calcUTau()

Foam::tmp< Foam::scalarField > calcUTau ( const scalarField magGradU) const
protectedvirtual

Calculate the friction velocity.

Definition at line 64 of file nutUBlendedWallFunctionFvPatchScalarField.C.

References forAll, Foam::constant::atomic::group, IOobject::groupName(), Foam::log(), Foam::mag(), magUp, n, turbulenceModel::nu(), Foam::foamVersion::patch, fvPatchField< Type >::patchInternalField(), Foam::pow(), turbulenceModel::propertiesName, tmp< T >::ref(), Foam::sqrt(), U, y, turbulenceModel::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:

◆ TypeName()

TypeName ( "nutUBlendedWallFunction"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchScalarField> clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 200 of file nutUBlendedWallFunctionFvPatchScalarField.H.

References nutUBlendedWallFunctionFvPatchScalarField::nutUBlendedWallFunctionFvPatchScalarField().

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

References nutUBlendedWallFunctionFvPatchScalarField::nutUBlendedWallFunctionFvPatchScalarField().

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

References Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), turbulenceModel::nu(), Foam::foamVersion::patch, turbulenceModel::propertiesName, fvPatchField< Type >::snGrad(), 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 211 of file nutUBlendedWallFunctionFvPatchScalarField.C.

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

Here is the call graph for this function:

Member Data Documentation

◆ n_

scalar n_
protected

Model coefficient; default = 4.

Definition at line 147 of file nutUBlendedWallFunctionFvPatchScalarField.H.


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