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:
[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 estimated y+ at the two-sublayer intersection. More...
 
scalar blend (const scalar nutVis, const scalar nutLog, const scalar yPlus) const
 Return the blended nut according to the chosen blending treatment. More...
 
virtual void updateCoeffs ()
 Update the coefficients associated with the patch field. More...
 

Protected Member Functions

virtual tmp< scalarFieldcalcNut () const
 Calculate the turbulent 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
enum blendingType blending_
 Blending treatment (default = blendingType::STEPWISE) More...
 
const scalar n_
 
word UName_
 Name of velocity field. More...
 
scalar Cmu_
 Empirical model coefficient. More...
 
scalar kappa_
 von Kármán constant More...
 
scalar E_
 Wall roughness parameter. More...
 
scalar yPlusLam_
 

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)
 Estimate the y+ at the intersection of the two sublayers. More...
 
- Protected Types inherited from nutWallFunctionFvPatchScalarField
enum  blendingType { STEPWISE, MAX, BINOMIAL, EXPONENTIAL }
 Options for the blending treatment of viscous and inertial sublayers. More...
 
- Static Protected Attributes inherited from nutWallFunctionFvPatchScalarField
static const Enum< blendingTypeblendingTypeNames
 Names for blendingType. More...
 

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 147 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 209 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 226 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 os(), 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 156 of file nutUBlendedWallFunctionFvPatchScalarField.H.


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