Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
nutWallFunctionFvPatchScalarField Class Referenceabstract

This boundary condition provides a turbulent kinematic viscosity condition when using wall functions, based on turbulence kinetic energy. More...

Inheritance diagram for nutWallFunctionFvPatchScalarField:
Inheritance graph
[legend]
Collaboration diagram for nutWallFunctionFvPatchScalarField:
Collaboration graph
[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 &)
 Construct by mapping given. More...
 
 nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &)
 Construct as copy. More...
 
 nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 Construct as copy setting internal field reference. 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 scalar yPlusLam (const scalar kappa, const scalar E)
 Calculate the Y+ at the edge of the laminar sublayer. More...
 

Protected Member Functions

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

scalar Cmu_
 Cmu coefficient. More...
 
scalar kappa_
 Von Karman constant. More...
 
scalar E_
 E coefficient. More...
 
scalar yPlusLam_
 Y+ at the edge of the laminar sublayer. More...
 

Detailed Description

This boundary condition provides a turbulent kinematic viscosity condition when using wall functions, based on turbulence kinetic energy.

- replicates OpenFOAM v1.5 (and earlier) behaviour

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>
    {
        type            nutWallFunction;
        value           uniform 0.0;
    }

Reference for the default model coefficients:

        H. Versteeg, W. Malalasekera
        An Introduction to Computational Fluid Dynamics: The Finite Volume
        Method, subsection "3.5.2 k-epsilon model"
See also
Foam::fixedValueFvPatchField
Source files

Definition at line 100 of file nutWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

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

Construct from patch and internal field.

Definition at line 69 of file nutWallFunctionFvPatchScalarField.C.

References nutWallFunctionFvPatchScalarField::checkType().

Referenced by nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField(), and nutWallFunctionFvPatchScalarField::writeLocalEntries().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Construct from patch, internal field and dictionary.

Definition at line 103 of file nutWallFunctionFvPatchScalarField.C.

References nutWallFunctionFvPatchScalarField::checkType(), and nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField().

Here is the call graph for this function:

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 85 of file nutWallFunctionFvPatchScalarField.C.

References nutWallFunctionFvPatchScalarField::checkType(), and nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField().

Here is the call graph for this function:

Construct as copy.

Definition at line 120 of file nutWallFunctionFvPatchScalarField.C.

References nutWallFunctionFvPatchScalarField::checkType(), and nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField().

Here is the call graph for this function:

Construct as copy setting internal field reference.

Definition at line 135 of file nutWallFunctionFvPatchScalarField.C.

References nutWallFunctionFvPatchScalarField::checkType(), and nutWallFunctionFvPatchScalarField::yPlusLam().

Here is the call graph for this function:

Member Function Documentation

void checkType ( )
protectedvirtual

Check the type of the patch.

Definition at line 44 of file nutWallFunctionFvPatchScalarField.C.

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

Referenced by nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField().

Here is the call graph for this function:

Here is the caller graph for this function:

virtual tmp<scalarField> calcNut ( ) const
protectedpure virtual
void writeLocalEntries ( Ostream os) const
protectedvirtual
TypeName ( "nutWallFunction"  )

Runtime type information.

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

Calculate the Y+ at the edge of the laminar sublayer.

Definition at line 153 of file nutWallFunctionFvPatchScalarField.C.

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

Referenced by nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField().

Here is the call graph for this function:

Here is the caller graph for this function:

virtual tmp<scalarField> yPlus ( ) const
pure virtual
void updateCoeffs ( )
virtual

Update the coefficients associated with the patch field.

Definition at line 169 of file nutWallFunctionFvPatchScalarField.C.

References nutWallFunctionFvPatchScalarField::calcNut(), and Foam::operator==().

Here is the call graph for this function:

void write ( Ostream os) const
virtual

Member Data Documentation

scalar Cmu_
protected
scalar kappa_
protected
scalar E_
protected
scalar yPlusLam_
protected

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