40 const label patchi = patch().index();
47 internalField().group()
65 /(magGradU + ROOTVSMALL)
78 auto&
uPlus = tuPlus.ref();
82 uPlus[facei] = uPlusTable_.interpolateLog10(
Rey[facei]);
100Foam::nutUTabulatedWallFunctionFvPatchScalarField::
101nutUTabulatedWallFunctionFvPatchScalarField
108 uPlusTableName_(
"undefined-uPlusTableName"),
125Foam::nutUTabulatedWallFunctionFvPatchScalarField::
126nutUTabulatedWallFunctionFvPatchScalarField
135 uPlusTableName_(ptf.uPlusTableName_),
136 uPlusTable_(ptf.uPlusTable_)
140Foam::nutUTabulatedWallFunctionFvPatchScalarField::
141nutUTabulatedWallFunctionFvPatchScalarField
149 uPlusTableName_(
dict.get<
word>(
"uPlusTable")),
166Foam::nutUTabulatedWallFunctionFvPatchScalarField::
167nutUTabulatedWallFunctionFvPatchScalarField
173 uPlusTableName_(wfpsf.uPlusTableName_),
174 uPlusTable_(wfpsf.uPlusTable_)
178Foam::nutUTabulatedWallFunctionFvPatchScalarField::
179nutUTabulatedWallFunctionFvPatchScalarField
186 uPlusTableName_(wfpsf.uPlusTableName_),
187 uPlusTable_(wfpsf.uPlusTable_)
196 const label patchi = patch().index();
203 internalField().group()
217 return Rey/(calcUPlus(
Rey) + ROOTVSMALL);
227 writeLocalEntries(
os);
228 writeEntry(
"value",
os);
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
This boundary condition provides a wall constraint on the turbulent viscosity (i.e....
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
tmp< scalarField > calcUPlus(const scalarField &Rey) const
Calculate wall u+ from table.
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
constant condensation/saturation model.
A class for managing temporary objects.
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
#define forAll(list, i)
Loop across all elements in list.