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-2019 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 "turbulenceModel.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 turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
44  (
46  (
48  internalField().group()
49  )
50  );
51  const scalarField& y = turbModel.y()[patchi];
52  const tmp<scalarField> tnuw = turbModel.nu(patchi);
53  const scalarField& nuw = tnuw();
54 
55  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
56  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
57 
58  tmp<scalarField> tyPlus = yPlus(magUp);
59  scalarField& yPlus = tyPlus.ref();
60 
61  tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
62  scalarField& nutw = tnutw.ref();
63 
64  forAll(yPlus, facei)
65  {
66  const scalar Re = magUp[facei]*y[facei]/nuw[facei];
67 
68  if (sqr(yPlus[facei]) > Re)
69  {
70  nutw[facei] = nuw[facei]*(sqr(yPlus[facei])/Re - 1);
71  }
72  else
73  {
74  nutw[facei] = 0;
75  }
76  }
77 
78  return tnutw;
79 }
80 
81 
83 (
84  const scalarField& magUp
85 ) const
86 {
87  const label patchi = patch().index();
88 
89  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
90  (
92  (
94  internalField().group()
95  )
96  );
97  const scalarField& y = turbModel.y()[patchi];
98  const tmp<scalarField> tnuw = turbModel.nu(patchi);
99  const scalarField& nuw = tnuw();
100 
101  tmp<scalarField> tyPlus(new scalarField(patch().size(), 0.0));
102  scalarField& yPlus = tyPlus.ref();
103 
104  static const scalar c2 = 2.25/(90 - 2.25);
105  static const scalar c3 = 2.0*atan(1.0)/log(90/2.25);
106  static const scalar c4 = c3*log(2.25);
107 
108  // If KsPlus is based on YPlus the extra term added to the law
109  // of the wall will depend on yPlus
110  forAll(yPlus, facei)
111  {
112  if (Ks_[facei] > 0)
113  {
114  // Rough Walls
115 
116  const scalar c1 = 1/(90 - 2.25) + Cs_[facei];
117 
118  const scalar magUpara = magUp[facei];
119  const scalar Re = magUpara*y[facei]/nuw[facei];
120  const scalar kappaRe = kappa_*Re;
121 
122  scalar yp = yPlusLam_;
123  const scalar ryPlusLam = 1/yp;
124 
125  int iter = 0;
126  scalar yPlusLast = 0.0;
127 
128  const scalar dKsPlusdYPlus = Ks_[facei]/y[facei];
129 
130  do
131  {
132  yPlusLast = yp;
133 
134  // The non-dimensional roughness height
135  scalar KsPlus = yp*dKsPlusdYPlus;
136 
137  // The extra term in the law-of-the-wall
138  scalar yPlusGPrime = 0;
139  scalar E = E_;
140 
141  if (KsPlus >= 90)
142  {
143  const scalar t1 = 1 + Cs_[facei]*KsPlus;
144  E = E_/t1;
145  yPlusGPrime = Cs_[facei]*KsPlus/t1;
146  }
147  else if (KsPlus > 2.25)
148  {
149  const scalar t1 = c1*KsPlus - c2;
150  const scalar t2 = c3*log(KsPlus) - c4;
151  const scalar sint2 = sin(t2);
152  const scalar logt1 = log(t1);
153  E = E_/pow(t1, sint2);
154  yPlusGPrime =
155  (c1*sint2*KsPlus/t1) + (c3*logt1*cos(t2));
156  }
157 
158  const scalar yPlusMin = constant::mathematical::e/E;
159 
160  if (kappa_*yPlusMin > 1)
161  {
162  yp = max
163  (
164  (kappaRe + yp*(1 - yPlusGPrime))
165  /(1 + log(E*yp) - yPlusGPrime),
166  sqrt(Re)
167  );
168  }
169  else
170  {
171  if (log(E*yp) < kappa_*yp)
172  {
173  yp = max
174  (
175  (kappaRe + yp*(1 - yPlusGPrime))
176  /(1 + log(E*yp) - yPlusGPrime),
177  yPlusMin
178  );
179  }
180  else
181  {
182  yp = sqrt(Re);
183  }
184  }
185  } while(mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 20);
186 
187  yPlus[facei] = yp;
188  }
189  else
190  {
191  // Smooth Walls
192  const scalar Re = magUp[facei]*y[facei]/nuw[facei];
193  const scalar ryPlusLam = 1/yPlusLam_;
194 
195  int iter = 0;
196  scalar yp = yPlusLam_;
197  scalar yPlusLast = yp;
198 
199  do
200  {
201  yPlusLast = yp;
202  if (yp > yPlusLam_)
203  {
204  yp = (kappa_*Re + yp)/(1 + log(E_*yp));
205  }
206  else
207  {
208  yp = sqrt(Re);
209  }
210  } while(mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 20);
211 
212  yPlus[facei] = yp;
213  }
214  }
215 
216  return tyPlus;
217 }
218 
219 
220 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
221 
223 (
224  const fvPatch& p,
226 )
227 :
229  Ks_(p.size(), 0.0),
230  Cs_(p.size(), 0.0)
231 {}
232 
233 
235 (
236  const fvPatch& p,
238  const dictionary& dict
239 )
240 :
242  Ks_("Ks", dict, p.size()),
243  Cs_("Cs", dict, p.size())
244 {}
245 
246 
248 (
250  const fvPatch& p,
252  const fvPatchFieldMapper& mapper
253 )
254 :
255  nutUWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
256  Ks_(mapper(ptf.Ks_)),
257  Cs_(mapper(ptf.Cs_))
258 {}
259 
260 
262 (
264 )
265 :
267  Ks_(rwfpsf.Ks_),
268  Cs_(rwfpsf.Cs_)
269 {}
270 
271 
273 (
276 )
277 :
279  Ks_(rwfpsf.Ks_),
280  Cs_(rwfpsf.Cs_)
281 {}
282 
283 
284 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
285 
287 (
288  const fvPatchFieldMapper& m
289 )
290 {
291  nutUWallFunctionFvPatchScalarField::autoMap(m);
292  m(Ks_, Ks_);
293  m(Cs_, Cs_);
294 }
295 
296 
298 (
299  const fvPatchScalarField& ptf,
300  const labelList& addr
301 )
302 {
303  nutUWallFunctionFvPatchScalarField::rmap(ptf, addr);
304 
306  refCast<const nutURoughWallFunctionFvPatchScalarField>(ptf);
307 
308  Ks_.rmap(nrwfpsf.Ks_, addr);
309  Cs_.rmap(nrwfpsf.Cs_, addr);
310 }
311 
312 
314 {
316  writeLocalEntries(os);
317  writeEntry(os, "Cs", Cs_);
318  writeEntry(os, "Ks", Ks_);
319  writeEntry(os, "value", *this);
320 }
321 
322 
323 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
324 
326 (
329 );
330 
331 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
332 
333 } // End namespace Foam
334 
335 // ************************************************************************* //
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.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
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.
const volVectorField & U() const
Access function to velocity field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
const dimensionedScalar c2
Second radiation constant: default SI units: [m K].
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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:280
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
scalar magUp
scalar y
dimensionedScalar cos(const dimensionedScalar &ds)
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
const nearWallDist & y() const
Return the near wall distances.
virtual label size() const
Return size.
Definition: fvPatch.H:155
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar sin(const dimensionedScalar &ds)
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:194
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
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:381
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...
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
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.
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.