nutURoughWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 
39 tmp<scalarField> nutURoughWallFunctionFvPatchScalarField::calcNut() const
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 fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
53  const tmp<scalarField> tnuw = turbModel.nu(patchi);
54  const scalarField& nuw = tnuw();
55 
56  // The flow velocity at the adjacent cell centre
57  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
58 
59  tmp<scalarField> tyPlus = calcYPlus(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  if (yPlus[facei] > yPlusLam_)
68  {
69  const scalar Re = magUp[facei]*y[facei]/nuw[facei] + ROOTVSMALL;
70  nutw[facei] = nuw[facei]*(sqr(yPlus[facei])/Re - 1);
71  }
72  }
73 
74  return tnutw;
75 }
76 
77 
78 tmp<scalarField> nutURoughWallFunctionFvPatchScalarField::calcYPlus
79 (
80  const scalarField& magUp
81 ) const
82 {
83  const label patchi = patch().index();
84 
85  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
86  (
88  (
90  internalField().group()
91  )
92  );
93  const scalarField& y = turbModel.y()[patchi];
94  const tmp<scalarField> tnuw = turbModel.nu(patchi);
95  const scalarField& nuw = tnuw();
96 
97  tmp<scalarField> tyPlus(new scalarField(patch().size(), 0.0));
98  scalarField& yPlus = tyPlus.ref();
99 
100  if (roughnessHeight_ > 0.0)
101  {
102  // Rough Walls
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);
107 
108  //if (KsPlusBasedOnYPlus_)
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  const scalar magUpara = magUp[facei];
115  const scalar Re = magUpara*y[facei]/nuw[facei];
116  const scalar kappaRe = kappa_*Re;
117 
118  scalar yp = yPlusLam_;
119  const scalar ryPlusLam = 1.0/yp;
120 
121  int iter = 0;
122  scalar yPlusLast = 0.0;
123  scalar dKsPlusdYPlus = roughnessHeight_/y[facei];
124 
125  // Additional tuning parameter - nominally = 1
126  dKsPlusdYPlus *= roughnessFactor_;
127 
128  do
129  {
130  yPlusLast = yp;
131 
132  // The non-dimensional roughness height
133  scalar KsPlus = yp*dKsPlusdYPlus;
134 
135  // The extra term in the law-of-the-wall
136  scalar G = 0.0;
137 
138  scalar yPlusGPrime = 0.0;
139 
140  if (KsPlus >= 90)
141  {
142  const scalar t_1 = 1 + roughnessConstant_*KsPlus;
143  G = log(t_1);
144  yPlusGPrime = roughnessConstant_*KsPlus/t_1;
145  }
146  else if (KsPlus > 2.25)
147  {
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);
152  G = logt_1*sint_2;
153  yPlusGPrime =
154  (c_1*sint_2*KsPlus/t_1) + (c_3*logt_1*cos(t_2));
155  }
156 
157  scalar denom = 1.0 + log(E_*yp) - G - yPlusGPrime;
158  if (mag(denom) > VSMALL)
159  {
160  yp = (kappaRe + yp*(1 - yPlusGPrime))/denom;
161  }
162  } while
163  (
164  mag(ryPlusLam*(yp - yPlusLast)) > 0.0001
165  && ++iter < 10
166  && yp > VSMALL
167  );
168 
169  yPlus[facei] = max(0.0, yp);
170  }
171  }
172  }
173  else
174  {
175  // Smooth Walls
176  forAll(yPlus, facei)
177  {
178  const scalar magUpara = magUp[facei];
179  const scalar Re = magUpara*y[facei]/nuw[facei];
180  const scalar kappaRe = kappa_*Re;
181 
182  scalar yp = yPlusLam_;
183  const scalar ryPlusLam = 1.0/yp;
184 
185  int iter = 0;
186  scalar yPlusLast = 0.0;
187 
188  do
189  {
190  yPlusLast = yp;
191  yp = (kappaRe + yp)/(1.0 + log(E_*yp));
192 
193  } while (mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 10);
194 
195  yPlus[facei] = max(0.0, yp);
196  }
197  }
198 
199  return tyPlus;
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
204 
206 (
207  const fvPatch& p,
209 )
210 :
212  roughnessHeight_(Zero),
213  roughnessConstant_(Zero),
214  roughnessFactor_(Zero)
215 {}
216 
217 
219 (
221  const fvPatch& p,
223  const fvPatchFieldMapper& mapper
224 )
225 :
226  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
227  roughnessHeight_(ptf.roughnessHeight_),
228  roughnessConstant_(ptf.roughnessConstant_),
229  roughnessFactor_(ptf.roughnessFactor_)
230 {}
231 
232 
234 (
235  const fvPatch& p,
237  const dictionary& dict
238 )
239 :
241  roughnessHeight_(readScalar(dict.lookup("roughnessHeight"))),
242  roughnessConstant_(readScalar(dict.lookup("roughnessConstant"))),
243  roughnessFactor_(readScalar(dict.lookup("roughnessFactor")))
244 {}
245 
246 
248 (
250 )
251 :
253  roughnessHeight_(rwfpsf.roughnessHeight_),
254  roughnessConstant_(rwfpsf.roughnessConstant_),
255  roughnessFactor_(rwfpsf.roughnessFactor_)
256 {}
257 
258 
260 (
263 )
264 :
266  roughnessHeight_(rwfpsf.roughnessHeight_),
267  roughnessConstant_(rwfpsf.roughnessConstant_),
268  roughnessFactor_(rwfpsf.roughnessFactor_)
269 {}
270 
271 
272 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
273 
275 {
276  const label patchi = patch().index();
277 
278  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
279  (
281  (
283  internalField().group()
284  )
285  );
286  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
287  tmp<scalarField> magUp = mag(Uw.patchInternalField() - Uw);
288 
289  return calcYPlus(magUp());
290 }
291 
292 
294 {
296  writeLocalEntries(os);
297  os.writeKeyword("roughnessHeight")
298  << roughnessHeight_ << token::END_STATEMENT << nl;
299  os.writeKeyword("roughnessConstant")
300  << roughnessConstant_ << token::END_STATEMENT << nl;
301  os.writeKeyword("roughnessFactor")
302  << roughnessFactor_ << token::END_STATEMENT << nl;
303  writeEntry("value", os);
304 }
305 
306 
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308 
310 (
313 );
314 
315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316 
317 } // End namespace Foam
318 
319 // ************************************************************************* //
const char *const group
Group name for atomic constants.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
fvPatchField< vector > fvPatchVectorField
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions ...
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const volVectorField & U() const
Access function to velocity field.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
scalar magUp
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
dimensionedScalar cos(const dimensionedScalar &ds)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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.
static const zero Zero
Definition: zero.H:91
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar sin(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
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 > &)
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition: PtrList.H:54
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451