Go to the documentation of this file.
41 const label patchi =
patch().index();
48 internalField().
group()
63 sqr(
calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
84 nutw[facei] = this->operator[](facei);
99 return calcUTau(magGradU, maxIter_, err);
111 const label patchi =
patch().index();
118 internalField().
group()
134 err.setSize(
uTau.size());
139 scalar ut =
sqrt((nutw[facei] + nuw[facei])*magGradU[facei]);
149 scalar kUu =
min(kappa_*
magUp[facei]/ut, 50);
150 scalar fkUu =
exp(kUu) - 1 - kUu*(1 + 0.5*kUu);
153 - ut*
y[facei]/nuw[facei]
155 + 1/E_*(fkUu - 1.0/6.0*kUu*
sqr(kUu));
162 scalar uTauNew = ut +
f/df;
163 err[facei] =
mag((ut - uTauNew)/ut);
171 && err[facei] > tolerance_
320 const label patchi =
patch().index();
327 internalField().
group()
345 writeLocalEntries(
os);
346 writeEntry(
"value",
os);
const scalar tolerance_
Convergence tolerance.
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual ~nutUSpaldingWallFunctionFvPatchScalarField()
Destructor.
constexpr const char *const group
Group name for atomic constants.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual const volVectorField & U(const turbulenceModel &turb) const
const label maxIter_
Max iterations in calcNut.
dimensionedScalar exp(const dimensionedScalar &ds)
const nearWallDist & y() const
Return the near wall distances.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
#define forAll(list, i)
Loop across all elements in list.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
virtual tmp< scalarField > calcNut() const
Uncomment in case of intrumentation.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
virtual void write(Ostream &os) const
Write.
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,...
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const std::string patch
OpenFOAM patch number as a std::string.
nutUSpaldingWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fvPatchFieldMapper.
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,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
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.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...