nutWallFunctionFvPatchScalarField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016, 2019 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "wallFvPatch.H"
33 #include "turbulenceModel.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(nutWallFunctionFvPatchScalarField, 0);
41 }
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 const Foam::Enum
46 <
48 >
50 ({
51  { blendingType::STEPWISE , "stepwise" },
52  { blendingType::MAX , "max" },
53  { blendingType::BINOMIAL , "binomial" },
54  { blendingType::EXPONENTIAL, "exponential" }
55 });
56 
57 
58 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
59 
61 {
62  if (!isA<wallFvPatch>(patch()))
63  {
65  << "Invalid wall function specification" << nl
66  << " Patch type for patch " << patch().name()
67  << " must be wall" << nl
68  << " Current patch type is " << patch().type() << nl << endl
69  << abort(FatalError);
70  }
71 }
72 
73 
75 (
76  const turbulenceModel& turb
77 ) const
78 {
79  if (UName_ == word::null)
80  {
81  return turb.U();
82  }
83  else
84  {
85  return db().lookupObject<volVectorField>(UName_);
86  }
87 }
88 
89 
91 (
92  Ostream& os
93 ) const
94 {
95  os.writeEntry("blending", blendingTypeNames[blending_]);
96 
97  if (blending_ == blendingType::BINOMIAL)
98  {
99  os.writeEntry("n", n_);
100  }
101 
102  if (UName_ != word::null)
103  {
104  os.writeEntry("U", UName_);
105  }
106 
107  os.writeEntry("Cmu", Cmu_);
108  os.writeEntry("kappa", kappa_);
109  os.writeEntry("E", E_);
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
114 
116 (
117  const fvPatch& p,
119 )
120 :
121  fixedValueFvPatchScalarField(p, iF),
122  blending_(blendingType::STEPWISE),
123  n_(4.0),
124  UName_(word::null),
125  Cmu_(0.09),
126  kappa_(0.41),
127  E_(9.8),
128  yPlusLam_(yPlusLam(kappa_, E_))
129 {
130  checkType();
131 }
132 
133 
135 (
137  const fvPatch& p,
139  const fvPatchFieldMapper& mapper
140 )
141 :
142  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
143  blending_(ptf.blending_),
144  n_(ptf.n_),
145  UName_(ptf.UName_),
146  Cmu_(ptf.Cmu_),
147  kappa_(ptf.kappa_),
148  E_(ptf.E_),
149  yPlusLam_(ptf.yPlusLam_)
150 {
151  checkType();
152 }
153 
154 
156 (
157  const fvPatch& p,
159  const dictionary& dict
160 )
161 :
162  fixedValueFvPatchScalarField(p, iF, dict),
163  blending_
164  (
165  blendingTypeNames.getOrDefault
166  (
167  "blending",
168  dict,
169  blendingType::STEPWISE
170  )
171  ),
172  n_
173  (
174  dict.getCheckOrDefault<scalar>
175  (
176  "n",
177  4.0,
179  )
180  ),
181  UName_(dict.getOrDefault<word>("U", word::null)),
182  Cmu_(dict.getOrDefault<scalar>("Cmu", 0.09)),
183  kappa_
184  (
185  dict.getCheckOrDefault<scalar>("kappa", 0.41, scalarMinMax::ge(SMALL))
186  ),
187  E_(dict.getCheckOrDefault<scalar>("E", 9.8, scalarMinMax::ge(SMALL))),
188  yPlusLam_(yPlusLam(kappa_, E_))
189 {
190  checkType();
191 }
192 
193 
195 (
197 )
198 :
199  fixedValueFvPatchScalarField(wfpsf),
200  blending_(wfpsf.blending_),
201  n_(wfpsf.n_),
202  UName_(wfpsf.UName_),
203  Cmu_(wfpsf.Cmu_),
204  kappa_(wfpsf.kappa_),
205  E_(wfpsf.E_),
206  yPlusLam_(wfpsf.yPlusLam_)
207 {
208  checkType();
209 }
210 
211 
213 (
216 )
217 :
218  fixedValueFvPatchScalarField(wfpsf, iF),
219  blending_(wfpsf.blending_),
220  n_(wfpsf.n_),
221  UName_(wfpsf.UName_),
222  Cmu_(wfpsf.Cmu_),
223  kappa_(wfpsf.kappa_),
224  E_(wfpsf.E_),
225  yPlusLam_(wfpsf.yPlusLam_)
226 {
227  checkType();
228 }
229 
230 
231 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
232 
235 (
236  const turbulenceModel& turbModel,
237  const label patchi
238 )
239 {
240  return
241  refCast<const nutWallFunctionFvPatchScalarField>
242  (
243  turbModel.nut()().boundaryField()[patchi],
244  patchi
245  );
246 }
247 
248 
250 (
251  const scalar kappa,
252  const scalar E
253 )
254 {
255  scalar ypl = 11.0;
256 
257  for (label i = 0; i < 10; ++i)
258  {
259  ypl = log(max(E*ypl, 1.0))/kappa;
260  }
261 
262  return ypl;
263 }
264 
265 
267 (
268  const scalar nutVis,
269  const scalar nutLog,
270  const scalar yPlus
271 ) const
272 {
273  scalar nutw = 0.0;
274 
275  switch (blending_)
276  {
277  case blendingType::STEPWISE:
278  {
279  if (yPlus > yPlusLam_)
280  {
281  nutw = nutLog;
282  }
283  else
284  {
285  nutw = nutVis;
286  }
287  break;
288  }
289 
290  case blendingType::MAX:
291  {
292  // (PH:Eq. 27)
293  nutw = max(nutVis, nutLog);
294  break;
295  }
296 
297  case blendingType::BINOMIAL:
298  {
299  // (ME:Eqs. 15-16)
300  nutw =
301  pow
302  (
303  pow(nutVis, n_) + pow(nutLog, n_),
304  1.0/n_
305  );
306  break;
307  }
308 
309  case blendingType::EXPONENTIAL:
310  {
311  // (PH:Eq. 31)
312  const scalar Gamma = 0.01*pow4(yPlus)/(1.0 + 5.0*yPlus);
313  const scalar invGamma = 1.0/(Gamma + ROOTVSMALL);
314 
315  nutw = nutVis*exp(-Gamma) + nutLog*exp(-invGamma);
316  break;
317  }
318  }
319 
320  return nutw;
321 }
322 
323 
325 {
326  return yPlusLam_;
327 }
328 
329 
331 {
332  if (updated())
333  {
334  return;
335  }
336 
337  operator==(calcNut());
338 
339  fixedValueFvPatchScalarField::updateCoeffs();
340 }
341 
342 
344 (
345  Ostream& os
346 ) const
347 {
349  writeLocalEntries(os);
350  writeEntry("value", os);
351 }
352 
353 
354 // ************************************************************************* //
Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutWallFunctionFvPatchScalarField.C:116
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::nutWallFunctionFvPatchScalarField::blendingType
blendingType
Options for the blending treatment of viscous and inertial sublayers.
Definition: nutWallFunctionFvPatchScalarField.H:260
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::nutWallFunctionFvPatchScalarField::n_
const scalar n_
Definition: nutWallFunctionFvPatchScalarField.H:279
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::nutWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: nutWallFunctionFvPatchScalarField.C:330
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Foam::nutWallFunctionFvPatchScalarField::yPlusLam_
scalar yPlusLam_
Definition: nutWallFunctionFvPatchScalarField.H:297
nutWallFunctionFvPatchScalarField.H
wallFvPatch.H
Foam::nutWallFunctionFvPatchScalarField::U
virtual const volVectorField & U(const turbulenceModel &turb) const
Definition: nutWallFunctionFvPatchScalarField.C:75
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
fvPatchFieldMapper.H
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::nutWallFunctionFvPatchScalarField::kappa_
scalar kappa_
von Kármán constant
Definition: nutWallFunctionFvPatchScalarField.H:290
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::nutWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: nutWallFunctionFvPatchScalarField.C:344
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::nutWallFunctionFvPatchScalarField::blending_
enum blendingType blending_
Blending treatment (default = blendingType::STEPWISE)
Definition: nutWallFunctionFvPatchScalarField.H:275
Foam::dictionary::getCheckOrDefault
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:209
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::nutWallFunctionFvPatchScalarField::E_
scalar E_
Wall roughness parameter.
Definition: nutWallFunctionFvPatchScalarField.H:293
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::nutWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Empirical model coefficient.
Definition: nutWallFunctionFvPatchScalarField.H:287
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::nutWallFunctionFvPatchScalarField::UName_
word UName_
Name of velocity field.
Definition: nutWallFunctionFvPatchScalarField.H:284
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nutWallFunctionFvPatchScalarField::blend
scalar blend(const scalar nutVis, const scalar nutLog, const scalar yPlus) const
Return the blended nut according to the chosen blending treatment.
Definition: nutWallFunctionFvPatchScalarField.C:267
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::nutWallFunctionFvPatchScalarField::blendingTypeNames
static const Enum< blendingType > blendingTypeNames
Names for blendingType.
Definition: nutWallFunctionFvPatchScalarField.H:269
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::MinMax::ge
static MinMax< T > ge(const T &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:31
Foam::nutWallFunctionFvPatchScalarField::nutw
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
Definition: nutWallFunctionFvPatchScalarField.C:235
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
turbulenceModel.H
Foam::nutWallFunctionFvPatchScalarField::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Definition: nutWallFunctionFvPatchScalarField.C:91
Foam::nutWallFunctionFvPatchScalarField::yPlusLam
scalar yPlusLam() const
Return the estimated y+ at the two-sublayer intersection.
Definition: nutWallFunctionFvPatchScalarField.C:324
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::nutWallFunctionFvPatchScalarField::checkType
virtual void checkType()
Check the type of the patch.
Definition: nutWallFunctionFvPatchScalarField.C:60