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

This boundary condition provides a wall constraint on the turbulent viscosity, i.e. nut, based on velocity, i.e. U. Using Spalding's law gives a continuous nut profile to the wall. More...

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

Public Member Functions

 TypeName ("nutUSpaldingWallFunction")
 
 nutUSpaldingWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 
 nutUSpaldingWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 
 nutUSpaldingWallFunctionFvPatchScalarField (const nutUSpaldingWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 
 nutUSpaldingWallFunctionFvPatchScalarField (const nutUSpaldingWallFunctionFvPatchScalarField &)
 
virtual tmp< fvPatchScalarFieldclone () const
 
 nutUSpaldingWallFunctionFvPatchScalarField (const nutUSpaldingWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 
virtual tmp< fvPatchScalarFieldclone (const DimensionedField< scalar, volMesh > &iF) const
 
virtual ~nutUSpaldingWallFunctionFvPatchScalarField ()
 
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
 
virtual tmp< scalarFieldcalcUTau (const scalarField &magGradU, const label maxIter, scalarField &err) const
 
virtual void writeLocalEntries (Ostream &) const
 
- Protected Member Functions inherited from nutWallFunctionFvPatchScalarField
virtual const volVectorFieldU (const turbulenceModel &turb) const
 
virtual void checkType ()
 

Protected Attributes

const label maxIter_
 
const scalar tolerance_
 
- 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 Spalding's law gives a continuous nut profile to the wall.

\[ y^+ = u^+ + \frac{1}{E} \left[exp(\kappa u^+) - 1 - \kappa u^+\, - 0.5 (\kappa u^+)^2 - \frac{1}{6} (\kappa u^+)^3\right] \]

where

$ y^+ $ = Non-dimensional position
$ u^+ $ = Non-dimensional velocity
$ \kappa $ = von Kármán constant
Usage
Example of the boundary condition specification:
<patchName>
{
    // Mandatory entries (unmodifiable)
    type            nutUSpaldingWallFunction;

    // Optional entries (unmodifiable)
    maxIter         10;
    tolerance       0.0001;

    // Optional (inherited) entries
    ...
}

where the entries mean:

Property Description Type Req'd Dflt
type Type name: nutUBlendedWallFunction word yes -
maxIter Number of Newton-Raphson iterations label no 10
tolerance Convergence tolerance scalar no 0.0001

The inherited entries are elaborated in:

Note
  • Suffers from non-exact restart since correctNut() (called through turbulence->validate) returns a slightly different value every time it is called. This is since the seed for the Newton-Raphson iteration uses the current value of *this (=nut ).
  • This can be avoided by overriding the tolerance. This also switches on a pre-detection whether the current nut already satisfies the turbulence conditions and if so does not change it at all. This means that the nut only changes if it 'has really changed'. This probably should be used with a tight tolerance, to make sure to kick every iteration, e.g. maxIter 100; tolerance 1e-7;
Source files

Definition at line 141 of file nutUSpaldingWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ nutUSpaldingWallFunctionFvPatchScalarField() [1/5]

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

Definition at line 202 of file nutUSpaldingWallFunctionFvPatchScalarField.C.

Referenced by nutUSpaldingWallFunctionFvPatchScalarField::clone().

Here is the caller graph for this function:

◆ nutUSpaldingWallFunctionFvPatchScalarField() [2/5]

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

◆ nutUSpaldingWallFunctionFvPatchScalarField() [3/5]

◆ nutUSpaldingWallFunctionFvPatchScalarField() [4/5]

◆ nutUSpaldingWallFunctionFvPatchScalarField() [5/5]

◆ ~nutUSpaldingWallFunctionFvPatchScalarField()

Member Function Documentation

◆ calcNut()

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

◆ calcUTau() [1/2]

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

Definition at line 87 of file nutUSpaldingWallFunctionFvPatchScalarField.C.

Referenced by nutUSpaldingWallFunctionFvPatchScalarField::calcNut().

Here is the caller graph for this function:

◆ calcUTau() [2/2]

Foam::tmp< Foam::scalarField > calcUTau ( const scalarField magGradU,
const label  maxIter,
scalarField err 
) const
protectedvirtual

◆ writeLocalEntries()

void writeLocalEntries ( Ostream os) const
protectedvirtual

Reimplemented from nutWallFunctionFvPatchScalarField.

Definition at line 187 of file nutUSpaldingWallFunctionFvPatchScalarField.C.

References os(), Ostream::writeEntryIfDifferent(), and nutWallFunctionFvPatchScalarField::writeLocalEntries().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "nutUSpaldingWallFunction"  )

◆ 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 333 of file nutUSpaldingWallFunctionFvPatchScalarField.C.

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

Here is the call graph for this function:

Member Data Documentation

◆ maxIter_

const label maxIter_
protected

◆ tolerance_

const scalar tolerance_
protected

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