Go to the documentation of this file.
39 const label patchi =
patch().index();
46 internalField().
group()
57 sqr(
calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
68 const label patchi =
patch().index();
75 internalField().
group()
97 scalar ut =
sqrt((nutw[facei] + nuw[facei])*magGradU[facei]);
98 if (
mag(ut) > ROOTVSMALL)
100 scalar
error = GREAT;
102 while (iter++ < 10 &&
error > 0.001)
104 const scalar
yPlus =
y[facei]*ut/nuw[facei];
109 pow(
pow(uTauVis, n_) +
pow(uTauLog, n_), 1.0/n_);
110 error =
mag(ut - utNew)/(ut + ROOTVSMALL);
111 ut = 0.5*(ut + utNew);
190 const label patchi =
patch().index();
197 internalField().
group()
206 return y*calcUTau(magGradU)/nuw;
216 writeLocalEntries(os);
218 writeEntry(
"value", os);
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
constexpr const char *const group
Group name for atomic constants.
A class for managing temporary objects.
virtual void write(Ostream &os) const
Write.
static constexpr const zero Zero
Global zero (0)
nutUBlendedWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual const volVectorField & U(const turbulenceModel &turb) const
const nearWallDist & y() const
Return the near wall distances.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
#define forAll(list, i)
Loop across all elements in list.
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
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,...
dimensionedScalar log(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
scalar n_
Model coefficient; default = 4.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const std::string patch
OpenFOAM patch number as a std::string.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
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 tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
Class to handle errors and exceptions in a simple, consistent stream-based manner.
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...