nutURoughWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "momentumTransportModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
40 {
41  const label patchi = patch().index();
42 
43  const momentumTransportModel& turbModel =
44  db().lookupObject<momentumTransportModel>
45  (
47  (
48  momentumTransportModel::typeName,
49  internalField().group()
50  )
51  );
52  const scalarField& y = turbModel.y()[patchi];
53  const tmp<scalarField> tnuw = turbModel.nu(patchi);
54  const scalarField& nuw = tnuw();
55 
56  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
57  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
58 
59  tmp<scalarField> tyPlus = yPlus(magUp);
60  scalarField& yPlus = tyPlus.ref();
61 
62  tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
63  scalarField& nutw = tnutw.ref();
64 
65  forAll(yPlus, facei)
66  {
67  const scalar Re = magUp[facei]*y[facei]/nuw[facei];
68 
69  if (sqr(yPlus[facei]) > Re)
70  {
71  nutw[facei] = nuw[facei]*(sqr(yPlus[facei])/Re - 1);
72  }
73  else
74  {
75  nutw[facei] = 0;
76  }
77  }
78 
79  return tnutw;
80 }
81 
82 
84 (
85  const scalarField& magUp
86 ) const
87 {
88  const label patchi = patch().index();
89 
90  const momentumTransportModel& turbModel =
91  db().lookupObject<momentumTransportModel>
92  (
94  (
95  momentumTransportModel::typeName,
96  internalField().group()
97  )
98  );
99  const scalarField& y = turbModel.y()[patchi];
100  const tmp<scalarField> tnuw = turbModel.nu(patchi);
101  const scalarField& nuw = tnuw();
102 
103  tmp<scalarField> tyPlus(new scalarField(patch().size(), 0.0));
104  scalarField& yPlus = tyPlus.ref();
105 
106  static const scalar c2 = 2.25/(90 - 2.25);
107  static const scalar c3 = 2.0*atan(1.0)/log(90/2.25);
108  static const scalar c4 = c3*log(2.25);
109 
110  // If KsPlus is based on YPlus the extra term added to the law
111  // of the wall will depend on yPlus
112  forAll(yPlus, facei)
113  {
114  if (Ks_[facei] > 0)
115  {
116  // Rough Walls
117 
118  const scalar c1 = 1/(90 - 2.25) + Cs_[facei];
119 
120  const scalar magUpara = magUp[facei];
121  const scalar Re = magUpara*y[facei]/nuw[facei];
122  const scalar kappaRe = kappa_*Re;
123 
124  scalar yp = yPlusLam_;
125  const scalar ryPlusLam = 1/yp;
126 
127  int iter = 0;
128  scalar yPlusLast = 0.0;
129 
130  const scalar dKsPlusdYPlus = Ks_[facei]/y[facei];
131 
132  do
133  {
134  yPlusLast = yp;
135 
136  // The non-dimensional roughness height
137  scalar KsPlus = yp*dKsPlusdYPlus;
138 
139  // The extra term in the law-of-the-wall
140  scalar yPlusGPrime = 0;
141  scalar E = E_;
142 
143  if (KsPlus >= 90)
144  {
145  const scalar t1 = 1 + Cs_[facei]*KsPlus;
146  E = E_/t1;
147  yPlusGPrime = Cs_[facei]*KsPlus/t1;
148  }
149  else if (KsPlus > 2.25)
150  {
151  const scalar t1 = c1*KsPlus - c2;
152  const scalar t2 = c3*log(KsPlus) - c4;
153  const scalar sint2 = sin(t2);
154  const scalar logt1 = log(t1);
155  E = E_/pow(t1, sint2);
156  yPlusGPrime =
157  (c1*sint2*KsPlus/t1) + (c3*logt1*cos(t2));
158  }
159 
160  const scalar yPlusMin = constant::mathematical::e/E;
161 
162  if (kappa_*yPlusMin > 1)
163  {
164  yp = max
165  (
166  (kappaRe + yp*(1 - yPlusGPrime))
167  /(1 + log(E*yp) - yPlusGPrime),
168  sqrt(Re)
169  );
170  }
171  else
172  {
173  if (log(E*yp) < kappa_*yp)
174  {
175  yp = max
176  (
177  (kappaRe + yp*(1 - yPlusGPrime))
178  /(1 + log(E*yp) - yPlusGPrime),
179  yPlusMin
180  );
181  }
182  else
183  {
184  yp = sqrt(Re);
185  }
186  }
187  } while(mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 20);
188 
189  yPlus[facei] = yp;
190  }
191  else
192  {
193  // Smooth Walls
194  const scalar Re = magUp[facei]*y[facei]/nuw[facei];
195  const scalar ryPlusLam = 1/yPlusLam_;
196 
197  int iter = 0;
198  scalar yp = yPlusLam_;
199  scalar yPlusLast = yp;
200 
201  do
202  {
203  yPlusLast = yp;
204  if (yp > yPlusLam_)
205  {
206  yp = (kappa_*Re + yp)/(1 + log(E_*yp));
207  }
208  else
209  {
210  yp = sqrt(Re);
211  }
212  } while(mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 20);
213 
214  yPlus[facei] = yp;
215  }
216  }
217 
218  return tyPlus;
219 }
220 
221 
222 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
223 
225 (
226  const fvPatch& p,
228 )
229 :
231  Ks_(p.size(), 0.0),
232  Cs_(p.size(), 0.0)
233 {}
234 
235 
237 (
238  const fvPatch& p,
240  const dictionary& dict
241 )
242 :
244  Ks_("Ks", dict, p.size()),
245  Cs_("Cs", dict, p.size())
246 {}
247 
248 
250 (
252  const fvPatch& p,
254  const fvPatchFieldMapper& mapper
255 )
256 :
257  nutUWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
258  Ks_(mapper(ptf.Ks_)),
259  Cs_(mapper(ptf.Cs_))
260 {}
261 
262 
264 (
267 )
268 :
270  Ks_(rwfpsf.Ks_),
271  Cs_(rwfpsf.Cs_)
272 {}
273 
274 
275 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
276 
278 (
279  const fvPatchFieldMapper& m
280 )
281 {
282  nutUWallFunctionFvPatchScalarField::autoMap(m);
283  m(Ks_, Ks_);
284  m(Cs_, Cs_);
285 }
286 
287 
289 (
290  const fvPatchScalarField& ptf,
291  const labelList& addr
292 )
293 {
294  nutUWallFunctionFvPatchScalarField::rmap(ptf, addr);
295 
297  refCast<const nutURoughWallFunctionFvPatchScalarField>(ptf);
298 
299  Ks_.rmap(nrwfpsf.Ks_, addr);
300  Cs_.rmap(nrwfpsf.Cs_, addr);
301 }
302 
303 
305 (
306  const fvPatchScalarField& ptf
307 )
308 {
309  nutUWallFunctionFvPatchScalarField::reset(ptf);
310 
312  refCast<const nutURoughWallFunctionFvPatchScalarField>(ptf);
313 
314  Ks_.reset(nrwfpsf.Ks_);
315  Cs_.reset(nrwfpsf.Cs_);
316 }
317 
318 
320 {
322  writeLocalEntries(os);
323  writeEntry(os, "Cs", Cs_);
324  writeEntry(os, "Ks", Ks_);
325  writeEntry(os, "value", *this);
326 }
327 
328 
329 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
330 
332 (
335 );
336 
337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338 
339 } // End namespace Foam
340 
341 // ************************************************************************* //
const char *const group
Group name for atomic constants.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions ...
dimensionedScalar log(const dimensionedScalar &ds)
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
const volScalarField::Boundary & y() const
Return the near wall distance.
Macros for easy insertion into run-time selection tables.
const dimensionedScalar c2
Second radiation constant: default SI units: [m K].
scalar magUp
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
scalar y
dimensionedScalar cos(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
dimensionedScalar sin(const dimensionedScalar &ds)
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:174
nutUWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Abstract base class for turbulence models (RAS, LES and laminar).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const volVectorField & U() const
Access function to velocity field.
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:362
nutURoughWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label patchi
dimensionedScalar atan(const dimensionedScalar &ds)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
void reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:432
Namespace for OpenFOAM.