nutUWallFunctionFvPatchScalarField Class Reference

This boundary condition provides a wall constraint on the turbulent viscosity, i.e. nut, based on velocity, i.e. U, for low- and high-Reynolds number turbulence models. More...

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

Public Member Functions

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

Protected Member Functions

virtual tmp< scalarFieldcalcYPlus (const scalarField &magUp) const
 Calculate yPlus. More...
 
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, based on velocity, i.e. U, for low- and high-Reynolds number turbulence models.

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

    // Optional (inherited) entries
    ...
}

where the entries mean:

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

The inherited entries are elaborated in:

Note
Source files

Definition at line 98 of file nutUWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ nutUWallFunctionFvPatchScalarField() [1/5]

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

Construct from patch and internal field.

Definition at line 129 of file nutUWallFunctionFvPatchScalarField.C.

Referenced by nutUWallFunctionFvPatchScalarField::clone().

Here is the caller graph for this function:

◆ nutUWallFunctionFvPatchScalarField() [2/5]

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

Construct from patch, internal field and dictionary.

Definition at line 151 of file nutUWallFunctionFvPatchScalarField.C.

◆ nutUWallFunctionFvPatchScalarField() [3/5]

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

Construct by mapping given nutUWallFunctionFvPatchScalarField onto a new patch

Definition at line 139 of file nutUWallFunctionFvPatchScalarField.C.

◆ nutUWallFunctionFvPatchScalarField() [4/5]

Construct as copy.

Definition at line 162 of file nutUWallFunctionFvPatchScalarField.C.

◆ nutUWallFunctionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 171 of file nutUWallFunctionFvPatchScalarField.C.

Member Function Documentation

◆ calcYPlus()

Foam::tmp< Foam::scalarField > calcYPlus ( const scalarField magUp) const
protectedvirtual

Calculate yPlus.

Definition at line 81 of file nutUWallFunctionFvPatchScalarField.C.

References forAll, Foam::constant::atomic::group, IOobject::groupName(), Foam::log(), Foam::mag(), magUp, Foam::max(), turbulenceModel::nu(), Foam::foamVersion::patch, turbulenceModel::propertiesName, tmp< T >::ref(), y, turbulenceModel::y(), yPlus, and Foam::Zero.

Referenced by nutUWallFunctionFvPatchScalarField::calcNut().

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

◆ calcNut()

◆ TypeName()

TypeName ( "nutUWallFunction"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchScalarField> clone ( ) const
inlinevirtual

Construct and return a clone.

Reimplemented in atmNutUWallFunctionFvPatchScalarField.

Definition at line 154 of file nutUWallFunctionFvPatchScalarField.H.

References nutUWallFunctionFvPatchScalarField::nutUWallFunctionFvPatchScalarField().

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.

Reimplemented in atmNutUWallFunctionFvPatchScalarField.

Definition at line 171 of file nutUWallFunctionFvPatchScalarField.H.

References nutUWallFunctionFvPatchScalarField::nutUWallFunctionFvPatchScalarField().

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 183 of file nutUWallFunctionFvPatchScalarField.C.

References forAll, Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), magUp, Foam::foamVersion::patch, fvPatchField< Type >::patchInternalField(), turbulenceModel::propertiesName, tmp< T >::ref(), fvPatchField< Type >::snGrad(), Foam::sqrt(), U, y, and yPlus.

Referenced by nutUWallFunctionFvPatchScalarField::calcNut().

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

◆ write()

void write ( Ostream os) const
virtual

Write.

Reimplemented from nutWallFunctionFvPatchScalarField.

Reimplemented in atmNutUWallFunctionFvPatchScalarField.

Definition at line 226 of file nutUWallFunctionFvPatchScalarField.C.

References os(), and fvPatchField< Type >::write().

Here is the call graph for this function:

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