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 turbulence kinematic viscosity, i.e. nut, for low- and high-Reynolds number turbulence models. 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 y+ at the edge of the viscous sublayer. 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)
 Calculate the y+ at the edge of the viscous sublayer. 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 turbulence viscosity. More...
 
virtual void writeLocalEntries (Ostream &) const
 Write local wall function variables. More...
 

Protected Attributes

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

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 turbulence kinematic viscosity, i.e. nut, for low- and high-Reynolds number turbulence models.

Reference:

    Default model coefficients:
        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".
Usage
Property Description Required Default value
Cmu Cmu coefficient no 0.09
kappa Von Karman constant no 0.41
E E coefficient no 9.8

Examples of the boundary condition specification:

    <patchName>
    {
        // Mandatory entries
        type            nutWallFunction;

        // Optional entries
    }
See also
Foam::fixedValueFvPatchField
Source files

Definition at line 110 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 94 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 130 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 111 of file nutWallFunctionFvPatchScalarField.C.

◆ nutWallFunctionFvPatchScalarField() [4/5]

Construct as copy.

Definition at line 148 of file nutWallFunctionFvPatchScalarField.C.

◆ nutWallFunctionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 164 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 60 of file nutWallFunctionFvPatchScalarField.C.

References word::null, and turb.

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

Here is the caller graph for this function:

◆ checkType()

void checkType ( )
protectedvirtual

Check the type of the patch.

Definition at line 45 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 76 of file nutWallFunctionFvPatchScalarField.C.

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

Referenced by nutkFilmWallFunctionFvPatchScalarField::write(), nutkAtmRoughWallFunctionFvPatchScalarField::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 203 of file nutWallFunctionFvPatchScalarField.H.

References nutWallFunctionFvPatchScalarField::Cmu_.

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

Here is the caller graph for this function:

◆ kappa()

scalar kappa ( ) const
inline

Return kappa.

Definition at line 209 of file nutWallFunctionFvPatchScalarField.H.

References nutWallFunctionFvPatchScalarField::kappa_.

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

Here is the caller graph for this function:

◆ E()

scalar E ( ) const
inline

Return E.

Definition at line 215 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 184 of file nutWallFunctionFvPatchScalarField.C.

References turbulenceModel::nut().

Referenced by nutkWallFunctionFvPatchScalarField::calcNut(), nutUWallFunctionFvPatchScalarField::calcNut(), nutkAtmRoughWallFunctionFvPatchScalarField::calcNut(), nutkRoughWallFunctionFvPatchScalarField::calcNut(), nutUSpaldingWallFunctionFvPatchScalarField::calcNut(), epsilonWallFunctionFvPatchScalarField::calculate(), omegaWallFunctionFvPatchScalarField::calculate(), fWallFunctionFvPatchScalarField::updateCoeffs(), v2WallFunctionFvPatchScalarField::updateCoeffs(), 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

Calculate the y+ at the edge of the viscous sublayer.

Definition at line 199 of file nutWallFunctionFvPatchScalarField.C.

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

Referenced by 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 y+ at the edge of the viscous sublayer.

Definition at line 215 of file nutWallFunctionFvPatchScalarField.C.

◆ yPlus()

virtual tmp<scalarField> yPlus ( ) const
pure virtual

◆ updateCoeffs()

void updateCoeffs ( )
virtual

Update the coefficients associated with the patch field.

Definition at line 221 of file nutWallFunctionFvPatchScalarField.C.

References Foam::operator==().

Here is the call graph for this function:

◆ write()

void write ( Ostream os) const
virtual

Member Data Documentation

◆ UName_

word UName_
protected

Name of velocity field.

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

Definition at line 121 of file nutWallFunctionFvPatchScalarField.H.

◆ Cmu_

◆ kappa_

◆ E_

◆ yPlusLam_

scalar yPlusLam_
protected

Estimated y+ value at the edge of the viscous sublayer.

Definition at line 133 of file nutWallFunctionFvPatchScalarField.H.

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


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