nutLowReWallFunctionFvPatchScalarField Class Reference

This boundary condition provides a wall constraint on the turbulent viscosity, i.e. nut, for low Reynolds number models. It sets nut to zero, and provides an access function to calculate y+. More...

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

Public Member Functions

 TypeName ("nutLowReWallFunction")
 Runtime type information. More...
 
 nutLowReWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Construct from patch and internal field. More...
 
 nutLowReWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 Construct from patch, internal field and dictionary. More...
 
 nutLowReWallFunctionFvPatchScalarField (const nutLowReWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 
 nutLowReWallFunctionFvPatchScalarField (const nutLowReWallFunctionFvPatchScalarField &)
 Construct as copy. More...
 
virtual tmp< fvPatchScalarFieldclone () const
 Construct and return a clone. More...
 
 nutLowReWallFunctionFvPatchScalarField (const nutLowReWallFunctionFvPatchScalarField &, 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...
 
- 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...
 
virtual void write (Ostream &) const
 Write. More...
 

Protected Member Functions

virtual tmp< scalarFieldcalcNut () const
 Calculate the turbulent viscosity. 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...
 

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...
 
- 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_
 
- 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, for low Reynolds number models. It sets nut to zero, and provides an access function to calculate y+.

Usage
Example of the boundary condition specification:
<patchName>
{
    // Mandatory entries (unmodifiable)
    type            nutLowReWallFunction;

    // Optional (inherited) entries
    ...
}

where the entries mean:

Property Description Type Req'd Dflt
type Type name: nutLowReWallFunction word yes -

The inherited entries are elaborated in:

Source files

Definition at line 90 of file nutLowReWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ nutLowReWallFunctionFvPatchScalarField() [1/5]

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

Construct from patch and internal field.

Definition at line 49 of file nutLowReWallFunctionFvPatchScalarField.C.

Referenced by nutLowReWallFunctionFvPatchScalarField::clone().

Here is the caller graph for this function:

◆ nutLowReWallFunctionFvPatchScalarField() [2/5]

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

Construct from patch, internal field and dictionary.

Definition at line 73 of file nutLowReWallFunctionFvPatchScalarField.C.

◆ nutLowReWallFunctionFvPatchScalarField() [3/5]

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

Construct by mapping given nutLowReWallFunctionFvPatchScalarField onto a new patch

Definition at line 60 of file nutLowReWallFunctionFvPatchScalarField.C.

◆ nutLowReWallFunctionFvPatchScalarField() [4/5]

◆ nutLowReWallFunctionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 95 of file nutLowReWallFunctionFvPatchScalarField.C.

Member Function Documentation

◆ calcNut()

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

Calculate the turbulent viscosity.

Implements nutWallFunctionFvPatchScalarField.

Definition at line 39 of file nutLowReWallFunctionFvPatchScalarField.C.

References tmp< T >::New(), Foam::foamVersion::patch, and Foam::Zero.

Here is the call graph for this function:

◆ TypeName()

TypeName ( "nutLowReWallFunction"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchScalarField> clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 143 of file nutLowReWallFunctionFvPatchScalarField.H.

References nutLowReWallFunctionFvPatchScalarField::nutLowReWallFunctionFvPatchScalarField().

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 160 of file nutLowReWallFunctionFvPatchScalarField.H.

References nutLowReWallFunctionFvPatchScalarField::nutLowReWallFunctionFvPatchScalarField().

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 107 of file nutLowReWallFunctionFvPatchScalarField.C.

References Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), turbulenceModel::nu(), turbulenceModel::nuEff(), Foam::foamVersion::patch, turbulenceModel::propertiesName, fvPatchField< Type >::snGrad(), Foam::sqrt(), U, y, and turbulenceModel::y().

Here is the call graph for this function:

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