nutWallFunctionFvPatchScalarField.H
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 Class
28  Foam::nutWallFunctionFvPatchScalarField
29 
30 Group
31  grpWallFunctions
32 
33 Description
34  The class \c nutWallFunction is a base class that parents the derived
35  boundary conditions which provide a wall constraint on various fields, such
36  as turbulent viscosity, i.e. \c nut, or turbulent kinetic energy dissipation
37  rate, i.e. \c epsilon, for low- and high-Reynolds number turbulence models.
38  The class is not an executable itself, yet a provider for common entries
39  to its derived boundary conditions.
40 
41  Reference:
42  \verbatim
43  Default model coefficients (tag:VM):
44  Versteeg, H. K., & Malalasekera, W. (2011).
45  An introduction to computational fluid dynamics: The finite
46  volume method. Harlow: Pearson Education.
47  Subsection "3.5.2 k-epsilon model".
48 
49  Binomial blending of the viscous and inertial sublayers (tag:ME):
50  Menter, F., & Esch, T. (2001).
51  Elements of industrial heat transfer prediction.
52  In Proceedings of the 16th Brazilian Congress of Mechanical
53  Engineering (COBEM), November 2001. vol. 20, p. 117-127.
54 
55  Exponential/Max blending of the viscous and inertial sublayers (tag:PH):
56  Popovac, M., & Hanjalić, K. (2007).
57  Compound wall treatment for RANS computation of complex
58  turbulent flows and heat transfer.
59  Flow, turbulence and combustion, 78(2), 177-202.
60  DOI:10.1007/s10494-006-9067-x
61  \endverbatim
62 
63 Usage
64  Example of the boundary condition specification:
65  \verbatim
66  <patchName>
67  {
68  // Mandatory and other optional entries
69  ...
70 
71  // Optional (inherited) entries
72  Cmu 0.09;
73  kappa 0.41;
74  E 9.8;
75  blending stepwise;
76  n 4.0;
77  U U;
78 
79  // Optional (inherited) entries
80  ...
81  }
82  \endverbatim
83 
84  where the entries mean:
85  \table
86  Property | Description | Type | Req'd | Dflt
87  Cmu | Empirical model coefficient | scalar | no | 0.09
88  kappa | von Kármán constant | scalar | no | 0.41
89  E | Wall roughness parameter | scalar | no | 9.8
90  blending | Viscous/inertial sublayer blending | word | no | stepwise
91  n | Binomial blending exponent | scalar | no | 2.0
92  U | Name of the velocity field | word | no | U
93  \endtable
94 
95  The inherited entries are elaborated in:
96  - \link fixedValueFvPatchFields.H \endlink
97 
98  Options for the \c blending entry:
99  \verbatim
100  stepwise | Stepwise switch (discontinuous)
101  max | Maximum value switch (discontinuous)
102  binomial | Binomial blending (smooth)
103  exponential | Exponential blending (smooth)
104  \endverbatim
105 
106  wherein \c nut predictions for the viscous and inertial sublayers are
107  blended according to the following expressions:
108 
109  - \c stepwise (default):
110 
111  \f[
112  \nu_t = {\nu_t}_{log} \qquad if \quad y^+ > y^+_{lam}
113  \f]
114 
115  \f[
116  \nu_t = {\nu_t}_{vis} \qquad if \quad y^+ <= y^+_{lam}
117  \f]
118 
119  where
120  \vartable
121  {\nu_t}_{vis} | \f$\nu_t\f$ prediction in the viscous sublayer
122  {\nu_t}_{log} | \f$\nu_t\f$ prediction in the inertial sublayer
123  y^+ | estimated wall-normal height of the cell centre in wall units
124  y^+_{lam} | estimated intersection of the viscous and inertial sublayers
125  \endvartable
126 
127 
128  - \c max (PH:Eq. 27):
129 
130  \f[
131  \nu_t = max({\nu_t}_{vis}, {\nu_t}_{log})
132  \f]
133 
134 
135  - \c binomial (ME:Eqs. 15-16):
136 
137  \f[
138  \nu_t = (({\nu_t}_{vis})^n + ({\nu_t}_{log})^n)^{1/n}
139  \f]
140  where
141  \vartable
142  n | Binomial blending exponent
143  \endvartable
144 
145 
146  - \c exponential (PH:Eq. 32):
147 
148  \f[
149  \nu_t = {\nu_t}_{vis} \exp[-\Gamma] + {\nu_t}_{log} \exp[-1/\Gamma]
150  \f]
151 
152  where (PH:Eq. 31)
153  \vartable
154  \Gamma | Blending expression
155  \Gamma | \f$0.01 (y^+)^4 / (1.0 + 5.0 y^+)\f$
156  \endvartable
157 
158 See also
159  - Foam::fixedValueFvPatchField
160 
161 SourceFiles
162  nutWallFunctionFvPatchScalarField.C
163 
164 \*---------------------------------------------------------------------------*/
165 
166 #ifndef nutWallFunctionFvPatchScalarField_H
167 #define nutWallFunctionFvPatchScalarField_H
168 
169 #include "fixedValueFvPatchFields.H"
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 namespace Foam
174 {
175 
176 class turbulenceModel;
177 
178 /*---------------------------------------------------------------------------*\
179  Class nutWallFunctionFvPatchScalarField Declaration
180 \*---------------------------------------------------------------------------*/
181 
182 class nutWallFunctionFvPatchScalarField
183 :
184  public fixedValueFvPatchScalarField
185 {
186 protected:
187 
188  // Protected Enumerations
189 
190  //- Options for the blending treatment of viscous and inertial sublayers
191  enum blendingType
192  {
193  STEPWISE,
194  MAX,
195  BINOMIAL,
196  EXPONENTIAL
197  };
198 
199  //- Names for blendingType
200  static const Enum<blendingType> blendingTypeNames;
201 
202 
203  // Protected Data
204 
205  //- Blending treatment (default = blendingType::STEPWISE)
206  const enum blendingType blending_;
207 
208  //- Binomial blending exponent being used when
209  //- blendingType is blendingType::BINOMIAL (default = 4)
210  const scalar n_;
211 
212  //- Name of velocity field
213  // Default is null (not specified) in which case the velocity is
214  // retrieved from the turbulence model
215  word UName_;
216 
217  //- Empirical model coefficient
218  scalar Cmu_;
219 
220  //- von Kármán constant
221  scalar kappa_;
222 
223  //- Wall roughness parameter
224  scalar E_;
225 
226  //- Estimated y+ value at the intersection
227  //- of the viscous and inertial sublayers
228  scalar yPlusLam_;
229 
230 
231  // Protected Member Functions
232 
233  //- Helper to return the velocity field either from the turbulence
234  //- model (default) or the mesh database
235  virtual const volVectorField& U(const turbulenceModel& turb) const;
236 
237  //- Check the type of the patch
238  virtual void checkType();
239 
240  //- Calculate the turbulent viscosity
241  virtual tmp<scalarField> calcNut() const = 0;
242 
243  //- Write local wall function variables
244  virtual void writeLocalEntries(Ostream&) const;
245 
246 
247 public:
248 
249  //- Runtime type information
250  TypeName("nutWallFunction");
251 
252 
253  // Constructors
254 
255  //- Construct from patch and internal field
257  (
258  const fvPatch&,
260  );
261 
262  //- Construct from patch, internal field and dictionary
264  (
265  const fvPatch&,
267  const dictionary&
268  );
269 
270  //- Construct by mapping given
271  //- nutWallFunctionFvPatchScalarField
272  //- onto a new patch
274  (
276  const fvPatch&,
278  const fvPatchFieldMapper&
279  );
280 
281  //- Construct as copy
283  (
285  );
286 
287  //- Construct as copy setting internal field reference
289  (
292  );
293 
294  // No clone methods - abstract class
295 
296 
297  // Member Functions
298 
299  //- Return Cmu
300  scalar Cmu() const
301  {
302  return Cmu_;
303  }
304 
305  //- Return kappa
306  scalar kappa() const
307  {
308  return kappa_;
309  }
310 
311  //- Return E
312  scalar E() const
313  {
314  return E_;
315  }
316 
317  //- Return the nut patchField for the given wall patch
319  (
320  const turbulenceModel& turbModel,
321  const label patchi
322  );
323 
324  //- Estimate the y+ at the intersection of the two sublayers
325  static scalar yPlusLam(const scalar kappa, const scalar E);
326 
327  //- Return the estimated y+ at the two-sublayer intersection
328  scalar yPlusLam() const;
329 
330  //- Return the blended nut according to the chosen blending treatment
331  scalar blend
332  (
333  const scalar nutVis,
334  const scalar nutLog,
335  const scalar yPlus
336  ) const;
337 
338  //- Calculate and return the yPlus at the boundary
339  // yPlus is the first-cell-centre height from boundary in wall units
340  virtual tmp<scalarField> yPlus() const = 0;
341 
342 
343  // Evaluation functions
344 
345  //- Update the coefficients associated with the patch field
346  virtual void updateCoeffs();
347 
348 
349  // I-O
350 
351  //- Write
352  virtual void write(Ostream&) const;
353 };
354 
355 
356 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
357 
358 } // End namespace Foam
359 
360 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361 
362 #endif
363 
364 // ************************************************************************* //
Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutWallFunctionFvPatchScalarField.C:116
Foam::nutWallFunctionFvPatchScalarField::blendingType
blendingType
Options for the blending treatment of viscous and inertial sublayers.
Definition: nutWallFunctionFvPatchScalarField.H:260
Foam::Enum< blendingType >
Foam::nutWallFunctionFvPatchScalarField::n_
const scalar n_
Definition: nutWallFunctionFvPatchScalarField.H:279
Foam::nutWallFunctionFvPatchScalarField::Cmu
scalar Cmu() const
Return Cmu.
Definition: nutWallFunctionFvPatchScalarField.H:369
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::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::nutWallFunctionFvPatchScalarField::yPlusLam_
scalar yPlusLam_
Definition: nutWallFunctionFvPatchScalarField.H:297
Foam::nutWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const =0
Calculate and return the yPlus at the boundary.
Foam::nutWallFunctionFvPatchScalarField::U
virtual const volVectorField & U(const turbulenceModel &turb) const
Definition: nutWallFunctionFvPatchScalarField.C:75
Foam::nutWallFunctionFvPatchScalarField::EXPONENTIAL
"Exponential blending (smooth)"
Definition: nutWallFunctionFvPatchScalarField.H:265
Foam::nutWallFunctionFvPatchScalarField::E
scalar E() const
Return E.
Definition: nutWallFunctionFvPatchScalarField.H:381
Foam::nutWallFunctionFvPatchScalarField::kappa_
scalar kappa_
von Kármán constant
Definition: nutWallFunctionFvPatchScalarField.H:290
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:63
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::calcNut
virtual tmp< scalarField > calcNut() const =0
Calculate the turbulent viscosity.
Foam::nutWallFunctionFvPatchScalarField::blending_
enum blendingType blending_
Blending treatment (default = blendingType::STEPWISE)
Definition: nutWallFunctionFvPatchScalarField.H:275
Foam::nutWallFunctionFvPatchScalarField::TypeName
TypeName("nutWallFunction")
Runtime type information.
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::nutWallFunctionFvPatchScalarField::BINOMIAL
"Binomial blending (smooth)"
Definition: nutWallFunctionFvPatchScalarField.H:264
Foam::nutWallFunctionFvPatchScalarField::E_
scalar E_
Wall roughness parameter.
Definition: nutWallFunctionFvPatchScalarField.H:293
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::nutWallFunctionFvPatchScalarField::MAX
"Maximum value switch (discontinuous)"
Definition: nutWallFunctionFvPatchScalarField.H:263
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
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::nutWallFunctionFvPatchScalarField::blendingTypeNames
static const Enum< blendingType > blendingTypeNames
Names for blendingType.
Definition: nutWallFunctionFvPatchScalarField.H:269
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
fixedValueFvPatchFields.H
Foam::nutWallFunctionFvPatchScalarField::kappa
scalar kappa() const
Return kappa.
Definition: nutWallFunctionFvPatchScalarField.H:375
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::nutWallFunctionFvPatchScalarField::STEPWISE
"Stepwise switch (discontinuous)"
Definition: nutWallFunctionFvPatchScalarField.H:262
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField< vector, fvPatchField, volMesh >
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
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