Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
nutUBlendedWallFunctionFvPatchScalarField Class Reference

This boundary condition provides a wall constraint on 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 turbulence models. More...

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

Public Member Functions

 TypeName ("nutUBlendedWallFunction")
 
 nutUBlendedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 
 nutUBlendedWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 
 nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 
 nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &)
 
virtual tmp< fvPatchScalarFieldclone () const
 
 nutUBlendedWallFunctionFvPatchScalarField (const nutUBlendedWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 
virtual tmp< fvPatchScalarFieldclone (const DimensionedField< scalar, volMesh > &iF) const
 
virtual tmp< scalarFieldyPlus () const
 
virtual void write (Ostream &os) const
 
- Public Member Functions inherited from nutWallFunctionFvPatchScalarField
 TypeName ("nutWallFunction")
 
 nutWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 
 nutWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 
 nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 
 nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &)
 
 nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 
scalar Cmu () const
 
scalar kappa () const
 
scalar E () const
 
scalar yPlusLam () const
 
scalar blend (const scalar nutVis, const scalar nutLog, const scalar yPlus) const
 
virtual void updateCoeffs ()
 

Protected Member Functions

virtual tmp< scalarFieldcalcNut () const
 
virtual tmp< scalarFieldcalcUTau (const scalarField &magGradU) const
 
- Protected Member Functions inherited from nutWallFunctionFvPatchScalarField
virtual const volVectorFieldU (const turbulenceModel &turb) const
 
virtual void checkType ()
 
virtual void writeLocalEntries (Ostream &) const
 

Protected Attributes

scalar n_
 
- Protected Attributes inherited from nutWallFunctionFvPatchScalarField
enum blendingType blending_
 
const scalar n_
 
word UName_
 
scalar Cmu_
 
scalar kappa_
 
scalar E_
 
scalar yPlusLam_
 

Additional Inherited Members

- Static Public Member Functions inherited from nutWallFunctionFvPatchScalarField
static const nutWallFunctionFvPatchScalarFieldnutw (const turbulenceModel &turbModel, const label patchi)
 
static scalar yPlusLam (const scalar kappa, const scalar E)
 
- Protected Types inherited from nutWallFunctionFvPatchScalarField
enum  blendingType { STEPWISE, MAX, BINOMIAL, EXPONENTIAL }
 
- Static Protected Attributes inherited from nutWallFunctionFvPatchScalarField
static const Enum< blendingTypeblendingTypeNames
 

Detailed Description

This boundary condition provides a wall constraint on 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 turbulence models.

\[ 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 sublayer
$ u_{\tau,l} $ = 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 (unmodifiable)
    type            nutUBlendedWallFunction;

    // Optional entries (unmodifiable)
    n               4.0;

    // Optional (inherited) entries
    ...
 }

where the entries mean:

Property Description Type Req'd Dflt
type Type name: nutUBlendedWallFunction word yes -
n Blending factor scalar no 4.0

The inherited entries are elaborated in:

Note
See also
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 
)

Definition at line 118 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 
)

Definition at line 144 of file nutUBlendedWallFunctionFvPatchScalarField.C.

◆ nutUBlendedWallFunctionFvPatchScalarField() [3/5]

Definition at line 130 of file nutUBlendedWallFunctionFvPatchScalarField.C.

◆ nutUBlendedWallFunctionFvPatchScalarField() [4/5]

Definition at line 157 of file nutUBlendedWallFunctionFvPatchScalarField.C.

◆ nutUBlendedWallFunctionFvPatchScalarField() [5/5]

Definition at line 168 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

◆ TypeName()

TypeName ( "nutUBlendedWallFunction"  )

◆ clone() [1/2]

virtual tmp<fvPatchScalarField> clone ( ) const
inlinevirtual

◆ clone() [2/2]

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

◆ yPlus()

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

◆ write()

void write ( Ostream os) const
virtual

Reimplemented from nutWallFunctionFvPatchScalarField.

Definition at line 204 of file nutUBlendedWallFunctionFvPatchScalarField.C.

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

Here is the call graph for this function:

Member Data Documentation

◆ n_

scalar n_
protected

Definition at line 152 of file nutUBlendedWallFunctionFvPatchScalarField.H.


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