Go to the documentation of this file.
69 const scalar
Re =
magUp[facei]*
y[facei]/nuw[facei] + ROOTVSMALL;
70 nutw[facei] = nuw[facei]*(
sqr(
yPlus[facei])/
Re - 1);
100 if (roughnessHeight_ > 0.0)
103 const scalar c_1 = 1/(90 - 2.25) + roughnessConstant_;
104 static const scalar c_2 = 2.25/(90 - 2.25);
105 static const scalar c_3 = 2.0*
atan(1.0)/
log(90/2.25);
106 static const scalar c_4 = c_3*
log(2.25);
114 const scalar magUpara =
magUp[facei];
115 const scalar
Re = magUpara*
y[facei]/nuw[facei];
116 const scalar kappaRe = kappa_*
Re;
118 scalar yp = yPlusLam_;
119 const scalar ryPlusLam = 1.0/yp;
122 scalar yPlusLast = 0.0;
123 scalar dKsPlusdYPlus = roughnessHeight_/
y[facei];
126 dKsPlusdYPlus *= roughnessFactor_;
133 scalar KsPlus = yp*dKsPlusdYPlus;
138 scalar yPlusGPrime = 0.0;
142 const scalar t_1 = 1 + roughnessConstant_*KsPlus;
144 yPlusGPrime = roughnessConstant_*KsPlus/t_1;
146 else if (KsPlus > 2.25)
148 const scalar t_1 = c_1*KsPlus - c_2;
149 const scalar t_2 = c_3*
log(KsPlus) - c_4;
150 const scalar sint_2 =
sin(t_2);
151 const scalar logt_1 =
log(t_1);
154 (c_1*sint_2*KsPlus/t_1) + (c_3*logt_1*
cos(t_2));
157 scalar denom = 1.0 +
log(E_*yp) -
G - yPlusGPrime;
158 if (
mag(denom) > VSMALL)
160 yp = (kappaRe + yp*(1 - yPlusGPrime))/denom;
164 mag(ryPlusLam*(yp - yPlusLast)) > 0.0001
178 const scalar magUpara =
magUp[facei];
179 const scalar
Re = magUpara*
y[facei]/nuw[facei];
180 const scalar kappaRe = kappa_*
Re;
182 scalar yp = yPlusLam_;
183 const scalar ryPlusLam = 1.0/yp;
186 scalar yPlusLast = 0.0;
191 yp = (kappaRe + yp)/(1.0 +
log(E_*yp));
193 }
while (
mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 10);
303 writeEntry(
"value", os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
static word groupName(Name name, const word &group)
#define forAll(list, i)
Loop across all elements in list.
const dimensionedScalar G
Newtonian constant of gravitation.
A class for managing temporary objects.
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
dimensionedScalar sin(const dimensionedScalar &ds)
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions ...
rDeltaT dimensionedInternalField()
scalarField Re(const UList< complex > &cf)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
virtual tmp< scalarField > calcYPlus(const scalarField &magUp) const
Calculate yPLus.
static const word propertiesName
Default name of the turbulence properties dictionary.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
const char *const group
Group name for atomic constants.
const nearWallDist & y() const
Return the near wall distances.
dimensioned< scalar > mag(const dimensioned< Type > &)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
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 any number of values (e....
nutURoughWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
dimensionedScalar log(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
scalar roughnessFactor_
Scale factor.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
scalar roughnessConstant_
Constant.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Traits class for primitives.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
scalar roughnessHeight_
Height.
dimensionedScalar atan(const dimensionedScalar &ds)
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const volVectorField & U() const
Access function to velocity field.
dimensionedScalar cos(const dimensionedScalar &ds)