61 const scalar
Re =
magUp[facei]*
y[facei]/nuw[facei];
94 static const scalar
c2 = 2.25/(90 - 2.25);
95 static const scalar c3 = 2.0*
atan(1.0)/
log(90/2.25);
96 static const scalar c4 = c3*
log(2.25);
106 const scalar
c1 = 1/(90 - 2.25) + Cs_[facei];
108 const scalar magUpara =
magUp[facei];
109 const scalar
Re = magUpara*
y[facei]/nuw[facei];
113 const scalar ryPlusLam = 1/yp;
116 scalar yPlusLast = 0.0;
118 const scalar dKsPlusdYPlus = Ks_[facei]/
y[facei];
125 scalar KsPlus = yp*dKsPlusdYPlus;
128 scalar yPlusGPrime = 0;
133 const scalar t1 = 1 + Cs_[facei]*KsPlus;
135 yPlusGPrime = Cs_[facei]*KsPlus/t1;
137 else if (KsPlus > 2.25)
139 const scalar t1 =
c1*KsPlus -
c2;
140 const scalar t2 = c3*
log(KsPlus) - c4;
141 const scalar sint2 =
sin(t2);
142 const scalar logt1 =
log(t1);
145 (
c1*sint2*KsPlus/t1) + (c3*logt1*
cos(t2));
154 (kappaRe + yp*(1 - yPlusGPrime))
155 /(1 +
log(
E*yp) - yPlusGPrime),
165 (kappaRe + yp*(1 - yPlusGPrime))
166 /(1 +
log(
E*yp) - yPlusGPrime),
175 }
while(
mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 20);
182 const scalar
Re =
magUp[facei]*
y[facei]/nuw[facei];
187 scalar yPlusLast = yp;
200 }
while(
mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 20);
234 Ks_(mapper(ptf.Ks_)),
262 refCast<const nutURoughWallFunctionFvPatchScalarField>(ptf);
264 mapper(Ks_, nrwfpsf.Ks_);
265 mapper(Cs_, nrwfpsf.Cs_);
274 nutUWallFunctionFvPatchScalarField::reset(ptf);
277 refCast<const nutURoughWallFunctionFvPatchScalarField>(ptf);
279 Ks_.
reset(nrwfpsf.Ks_);
280 Cs_.
reset(nrwfpsf.Cs_);
#define forAll(list, i)
Loop across all elements in list.
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...
void reset(const Field< Type > &)
Reset the field values to the given field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Abstract base class for field mapping.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Abstract base class for turbulence models (RAS, LES and laminar).
const volScalarField::Boundary & yb() const
Return the near wall distance.
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions ...
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
nutURoughWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void write(Ostream &os) const
Write.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
scalar kappa_
Von Karman constant.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
scalar E() const
Return E.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
const dimensionedScalar c2
Second radiation constant: default SI units: [m K].
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimless
dimensionedScalar sin(const dimensionedScalar &ds)
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensionedScalar atan(const dimensionedScalar &ds)
dimensionedScalar cos(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cf)