nutWallFunctionFvPatchScalarField Class Referenceabstract

The class nutWallFunction is a base class that parents the derived boundary conditions which provide a wall constraint on various fields, such as turbulent viscosity, i.e. nut, or turbulent kinetic energy dissipation rate, i.e. epsilon, for low- and high-Reynolds number turbulence models. The class is not an executable itself, yet a provider for common entries to its derived boundary conditions. More...

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

Public Member Functions

 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 tmp< scalarFieldyPlus () const =0
 Calculate and return the yPlus at the boundary. More...
 
virtual void updateCoeffs ()
 Update the coefficients associated with the patch field. More...
 
virtual void write (Ostream &) const
 Write. More...
 

Static Public Member Functions

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

enum  blendingType { STEPWISE, MAX, BINOMIAL, EXPONENTIAL }
 Options for the blending treatment of viscous and inertial sublayers. More...
 

Protected Member Functions

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...
 
virtual void writeLocalEntries (Ostream &) const
 Write local wall function variables. More...
 

Protected Attributes

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_
 

Static Protected Attributes

static const Enum< blendingTypeblendingTypeNames
 Names for blendingType. More...
 

Detailed Description

The class nutWallFunction is a base class that parents the derived boundary conditions which provide a wall constraint on various fields, such as turbulent viscosity, i.e. nut, or turbulent kinetic energy dissipation rate, i.e. epsilon, for low- and high-Reynolds number turbulence models. The class is not an executable itself, yet a provider for common entries to its derived boundary conditions.

Reference:

    Default model coefficients (tag:VM):
        Versteeg, H. K., & Malalasekera, W. (2011).
        An introduction to computational fluid dynamics: The finite
        volume method. Harlow: Pearson Education.
        Subsection "3.5.2 k-epsilon model".

    Binomial blending of the viscous and inertial sublayers (tag:ME):
        Menter, F., & Esch, T. (2001).
        Elements of industrial heat transfer prediction.
        In Proceedings of the 16th Brazilian Congress of Mechanical
        Engineering (COBEM), November 2001. vol. 20, p. 117-127.

    Exponential/Max blending of the viscous and inertial sublayers (tag:PH):
        Popovac, M., & Hanjalić, K. (2007).
        Compound wall treatment for RANS computation of complex
        turbulent flows and heat transfer.
        Flow, turbulence and combustion, 78(2), 177-202.
        DOI:10.1007/s10494-006-9067-x
Usage
Example of the boundary condition specification:
<patchName>
{
    // Mandatory and other optional entries
    ...

    // Optional (inherited) entries
    Cmu             0.09;
    kappa           0.41;
    E               9.8;
    blending        stepwise;
    n               4.0;
    U               U;

    // Optional (inherited) entries
    ...
}

where the entries mean:

Property Description Type Req'd Dflt
Cmu Empirical model coefficient scalar no 0.09
kappa von Kármán constant scalar no 0.41
E Wall roughness parameter scalar no 9.8
blending Viscous/inertial sublayer blending word no stepwise
n Binomial blending exponent scalar no 2.0
U Name of the velocity field word no U

The inherited entries are elaborated in:

Options for the blending entry:

      stepwise    | Stepwise switch (discontinuous)
      max         | Maximum value switch (discontinuous)
      binomial    | Binomial blending (smooth)
      exponential | Exponential blending (smooth)

wherein nut predictions for the viscous and inertial sublayers are blended according to the following expressions:

  • stepwise (default):

\[ \nu_t = {\nu_t}_{log} \qquad if \quad y^+ > y^+_{lam} \]

\[ \nu_t = {\nu_t}_{vis} \qquad if \quad y^+ <= y^+_{lam} \]

where

\( {\nu_t}_{vis} \) = \(\nu_t\) prediction in the viscous sublayer
\( {\nu_t}_{log} \) = \(\nu_t\) prediction in the inertial sublayer
\( y^+ \) = estimated wall-normal height of the cell centre in wall units
\( y^+_{lam} \) = estimated intersection of the viscous and inertial sublayers
  • max (PH:Eq. 27):

\[ \nu_t = max({\nu_t}_{vis}, {\nu_t}_{log}) \]

  • binomial (ME:Eqs. 15-16):

\[ \nu_t = (({\nu_t}_{vis})^n + ({\nu_t}_{log})^n)^{1/n} \]

where

\( n \) = Binomial blending exponent
  • exponential (PH:Eq. 32):

\[ \nu_t = {\nu_t}_{vis} \exp[-\Gamma] + {\nu_t}_{log} \exp[-1/\Gamma] \]

where (PH:Eq. 31)

\( \Gamma \) = Blending expression
\( \Gamma \) = \(0.01 (y^+)^4 / (1.0 + 5.0 y^+)\)
See also
Source files

Definition at line 251 of file nutWallFunctionFvPatchScalarField.H.

Member Enumeration Documentation

◆ blendingType

enum blendingType
protected

Options for the blending treatment of viscous and inertial sublayers.

Enumerator
STEPWISE 

"Stepwise switch (discontinuous)"

MAX 

"Maximum value switch (discontinuous)"

BINOMIAL 

"Binomial blending (smooth)"

EXPONENTIAL 

"Exponential blending (smooth)"

Definition at line 260 of file nutWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ nutWallFunctionFvPatchScalarField() [1/5]

nutWallFunctionFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF 
)

Construct from patch and internal field.

Definition at line 116 of file nutWallFunctionFvPatchScalarField.C.

◆ nutWallFunctionFvPatchScalarField() [2/5]

nutWallFunctionFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF,
const dictionary dict 
)

Construct from patch, internal field and dictionary.

Definition at line 156 of file nutWallFunctionFvPatchScalarField.C.

◆ nutWallFunctionFvPatchScalarField() [3/5]

nutWallFunctionFvPatchScalarField ( const nutWallFunctionFvPatchScalarField ptf,
const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF,
const fvPatchFieldMapper mapper 
)

Construct by mapping given nutWallFunctionFvPatchScalarField onto a new patch

Definition at line 135 of file nutWallFunctionFvPatchScalarField.C.

◆ nutWallFunctionFvPatchScalarField() [4/5]

Construct as copy.

Definition at line 195 of file nutWallFunctionFvPatchScalarField.C.

◆ nutWallFunctionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 213 of file nutWallFunctionFvPatchScalarField.C.

Member Function Documentation

◆ U()

const Foam::volVectorField & U ( const turbulenceModel turb) const
protectedvirtual

Helper to return the velocity field either from the turbulence model (default) or the mesh database

Definition at line 75 of file nutWallFunctionFvPatchScalarField.C.

References word::null, and turb.

Referenced by nutUWallFunctionFvPatchScalarField::calcNut(), nutUTabulatedWallFunctionFvPatchScalarField::calcNut(), nutUBlendedWallFunctionFvPatchScalarField::calcNut(), nutUSpaldingWallFunctionFvPatchScalarField::calcNut(), and atmNutUWallFunctionFvPatchScalarField::calcNut().

Here is the caller graph for this function:

◆ checkType()

void checkType ( )
protectedvirtual

Check the type of the patch.

Definition at line 60 of file nutWallFunctionFvPatchScalarField.C.

References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, Foam::nl, and Foam::foamVersion::patch.

Here is the call graph for this function:

◆ calcNut()

◆ writeLocalEntries()

void writeLocalEntries ( Ostream os) const
protectedvirtual

Write local wall function variables.

Reimplemented in nutUSpaldingWallFunctionFvPatchScalarField.

Definition at line 91 of file nutWallFunctionFvPatchScalarField.C.

References word::null, os(), and Ostream::writeEntry().

Referenced by nutkFilmWallFunctionFvPatchScalarField::write(), atmNutUWallFunctionFvPatchScalarField::write(), atmNutkWallFunctionFvPatchScalarField::write(), atmNutWallFunctionFvPatchScalarField::write(), and nutUSpaldingWallFunctionFvPatchScalarField::writeLocalEntries().

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

◆ TypeName()

TypeName ( "nutWallFunction"  )

Runtime type information.

◆ Cmu()

scalar Cmu ( ) const
inline

Return Cmu.

Definition at line 369 of file nutWallFunctionFvPatchScalarField.H.

References nutWallFunctionFvPatchScalarField::Cmu_.

Referenced by atmOmegaWallFunctionFvPatchScalarField::calculate(), atmEpsilonWallFunctionFvPatchScalarField::calculate(), epsilonWallFunctionFvPatchScalarField::calculate(), and omegaWallFunctionFvPatchScalarField::calculate().

Here is the caller graph for this function:

◆ kappa()

scalar kappa ( ) const
inline

Return kappa.

Definition at line 375 of file nutWallFunctionFvPatchScalarField.H.

References nutWallFunctionFvPatchScalarField::kappa_.

Referenced by atmOmegaWallFunctionFvPatchScalarField::calculate(), atmEpsilonWallFunctionFvPatchScalarField::calculate(), epsilonWallFunctionFvPatchScalarField::calculate(), and omegaWallFunctionFvPatchScalarField::calculate().

Here is the caller graph for this function:

◆ E()

scalar E ( ) const
inline

Return E.

Definition at line 381 of file nutWallFunctionFvPatchScalarField.H.

References nutWallFunctionFvPatchScalarField::E_.

◆ nutw()

const Foam::nutWallFunctionFvPatchScalarField & nutw ( const turbulenceModel turbModel,
const label  patchi 
)
static

Return the nut patchField for the given wall patch.

Definition at line 235 of file nutWallFunctionFvPatchScalarField.C.

References turbulenceModel::nut().

Referenced by nutkWallFunctionFvPatchScalarField::calcNut(), nutUWallFunctionFvPatchScalarField::calcNut(), nutkRoughWallFunctionFvPatchScalarField::calcNut(), nutUSpaldingWallFunctionFvPatchScalarField::calcNut(), atmNutUWallFunctionFvPatchScalarField::calcNut(), atmNutkWallFunctionFvPatchScalarField::calcNut(), atmNutWallFunctionFvPatchScalarField::calcNut(), epsilonWallFunctionFvPatchScalarField::calculate(), omegaWallFunctionFvPatchScalarField::calculate(), and kLowReWallFunctionFvPatchScalarField::updateCoeffs().

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

◆ yPlusLam() [1/2]

Foam::scalar yPlusLam ( const scalar  kappa,
const scalar  E 
)
static

Estimate the y+ at the intersection of the two sublayers.

Definition at line 250 of file nutWallFunctionFvPatchScalarField.C.

References Foam::constant::electromagnetic::kappa, Foam::log(), and Foam::max().

Referenced by atmEpsilonWallFunctionFvPatchScalarField::calculate(), epsilonWallFunctionFvPatchScalarField::calculate(), and omegaWallFunctionFvPatchScalarField::calculate().

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

◆ yPlusLam() [2/2]

Foam::scalar yPlusLam ( ) const

Return the estimated y+ at the two-sublayer intersection.

Definition at line 324 of file nutWallFunctionFvPatchScalarField.C.

◆ blend()

Foam::scalar blend ( const scalar  nutVis,
const scalar  nutLog,
const scalar  yPlus 
) const

Return the blended nut according to the chosen blending treatment.

Definition at line 267 of file nutWallFunctionFvPatchScalarField.C.

References Foam::exp(), Foam::max(), Foam::pow(), Foam::pow4(), and yPlus.

Referenced by nutkWallFunctionFvPatchScalarField::calcNut(), and nutUWallFunctionFvPatchScalarField::calcNut().

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

◆ yPlus()

virtual tmp<scalarField> yPlus ( ) const
pure virtual

Calculate and return the yPlus at the boundary.

yPlus is the first-cell-centre height from boundary in wall units

Implemented in nutURoughWallFunctionFvPatchScalarField, nutUSpaldingWallFunctionFvPatchScalarField, nutUBlendedWallFunctionFvPatchScalarField, nutUTabulatedWallFunctionFvPatchScalarField, nutUWallFunctionFvPatchScalarField, nutkWallFunctionFvPatchScalarField, nutLowReWallFunctionFvPatchScalarField, and nutkFilmWallFunctionFvPatchScalarField.

Referenced by yPlus::execute(), and alphatJayatillekeWallFunctionFvPatchScalarField::yPlus().

Here is the caller graph for this function:

◆ updateCoeffs()

void updateCoeffs ( )
virtual

Update the coefficients associated with the patch field.

Definition at line 330 of file nutWallFunctionFvPatchScalarField.C.

References Foam::operator==().

Here is the call graph for this function:

◆ write()

Member Data Documentation

◆ blendingTypeNames

const Foam::Enum< Foam::nutWallFunctionFvPatchScalarField::blendingType > blendingTypeNames
staticprotected

Names for blendingType.

Definition at line 269 of file nutWallFunctionFvPatchScalarField.H.

◆ blending_

enum blendingType blending_
protected

Blending treatment (default = blendingType::STEPWISE)

Definition at line 275 of file nutWallFunctionFvPatchScalarField.H.

◆ n_

const scalar n_
protected

Binomial blending exponent being used when blendingType is blendingType::BINOMIAL (default = 4)

Definition at line 279 of file nutWallFunctionFvPatchScalarField.H.

◆ UName_

word UName_
protected

Name of velocity field.

Default is null (not specified) in which case the velocity is retrieved from the turbulence model

Definition at line 284 of file nutWallFunctionFvPatchScalarField.H.

◆ Cmu_

◆ kappa_

◆ E_

◆ yPlusLam_

scalar yPlusLam_
protected

Estimated y+ value at the intersection of the viscous and inertial sublayers

Definition at line 297 of file nutWallFunctionFvPatchScalarField.H.


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