Go to the documentation of this file.
166 #ifndef nutWallFunctionFvPatchScalarField_H
167 #define nutWallFunctionFvPatchScalarField_H
182 class nutWallFunctionFvPatchScalarField
184 public fixedValueFvPatchScalarField
241 virtual tmp<scalarField>
calcNut()
const = 0;
340 virtual tmp<scalarField>
yPlus()
const = 0;
352 virtual void write(Ostream&)
const;
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
blendingType
Options for the blending treatment of viscous and inertial sublayers.
scalar Cmu() const
Return Cmu.
A class for handling words, derived from Foam::string.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for managing temporary objects.
virtual tmp< scalarField > yPlus() const =0
Calculate and return the yPlus at the boundary.
virtual const volVectorField & U(const turbulenceModel &turb) const
"Exponential blending (smooth)"
scalar E() const
Return E.
scalar kappa_
von Kármán constant
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
virtual void write(Ostream &) const
Write.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual tmp< scalarField > calcNut() const =0
Calculate the turbulent viscosity.
enum blendingType blending_
Blending treatment (default = blendingType::STEPWISE)
TypeName("nutWallFunction")
Runtime type information.
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
"Binomial blending (smooth)"
scalar E_
Wall roughness parameter.
scalar Cmu_
Empirical model coefficient.
word UName_
Name of velocity field.
"Maximum value switch (discontinuous)"
GeometricField< vector, fvPatchField, volMesh > volVectorField
scalar blend(const scalar nutVis, const scalar nutLog, const scalar yPlus) const
Return the blended nut according to the chosen blending treatment.
static const Enum< blendingType > blendingTypeNames
Names for blendingType.
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
scalar kappa() const
Return kappa.
Foam::fvPatchFieldMapper.
"Stepwise switch (discontinuous)"
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
compressible::turbulenceModel & turb
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
scalar yPlusLam() const
Return the estimated y+ at the two-sublayer intersection.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void checkType()
Check the type of the patch.