omegaWallFunctionFvPatchScalarField.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::omegaWallFunctionFvPatchScalarField
29 
30 Group
31  grpWallFunctions
32 
33 Description
34  This boundary condition provides a wall constraint on the specific
35  dissipation rate, i.e. \c omega, and the turbulent kinetic energy
36  production contribution, i.e. \c G, for low- and high-Reynolds number
37  turbulence models.
38 
39  Reference:
40  \verbatim
41  Binomial blending of the viscous and inertial sublayers (tag:ME):
42  Menter, F., & Esch, T. (2001).
43  Elements of industrial heat transfer prediction.
44  In Proceedings of the 16th Brazilian Congress of Mechanical
45  Engineering (COBEM), November 2001. vol. 20, p. 117-127.
46 
47  Exponential/Max blending of the viscous and inertial sublayers (tag:PH):
48  Popovac, M., & Hanjalić, K. (2007).
49  Compound wall treatment for RANS computation of complex
50  turbulent flows and heat transfer.
51  Flow, turbulence and combustion, 78(2), 177-202.
52  DOI:10.1007/s10494-006-9067-x
53 
54  Tanh blending of the viscous and inertial sublayers (tag:KAS):
55  Knopp, T., Alrutz, T., & Schwamborn, D. (2006).
56  A grid and flow adaptive wall-function method for RANS
57  turbulence modelling.
58  Journal of Computational Physics, 220(1), 19-40.
59  DOI:10.1016/j.jcp.2006.05.003
60  \endverbatim
61 
62 Usage
63  Example of the boundary condition specification:
64  \verbatim
65  <patchName>
66  {
67  // Mandatory entries (unmodifiable)
68  type omegaWallFunction;
69 
70  // Optional entries (unmodifiable)
71  beta1 0.075;
72  blending binomial2;
73  n 2.0;
74 
75  // Optional (inherited) entries
76  ...
77  }
78  \endverbatim
79 
80  \table
81  Property | Description | Type | Req'd | Dflt
82  type | Type name: omegaWallFunction | word | yes | -
83  beta1 | Model coefficient | scalar | no | 0.075
84  blending | Viscous/inertial sublayer blending method <!--
85  --> | word | no | binomial2
86  n | Binomial blending exponent | scalar | no | 2.0
87  \endtable
88 
89  The inherited entries are elaborated in:
90  - \link fixedValueFvPatchField.H \endlink
91  - \link nutWallFunctionFvPatchScalarField.H \endlink
92 
93  Options for the \c blending entry:
94  \verbatim
95  stepwise | Stepwise switch (discontinuous)
96  max | Maximum value switch (discontinuous)
97  binomial2 | Binomial blending (smooth) n = 2
98  binomial | Binomial blending (smooth)
99  exponential | Exponential blending (smooth)
100  tanh | Tanh blending (smooth)
101  \endverbatim
102 
103  wherein \c omega predictions for the viscous and inertial sublayers are
104  blended according to the following expressions:
105 
106  - \c stepwise:
107 
108  \f[
109  \omega = \omega_{log} \qquad if \quad y^+ > y^+_{lam}
110  \f]
111  \f[
112  \omega = \omega_{vis} \qquad if \quad y^+ <= y^+_{lam}
113  \f]
114 
115  where
116  \vartable
117  \omega | \f$\omega\f$ at \f$y^+\f$
118  \omega_{vis} | \f$\omega\f$ computed by using viscous sublayer assumptions
119  \omega_{log} |\f$\omega\f$ computed by using inertial sublayer assumptions
120  y^+ | estimated wall-normal height of the cell centre in wall units
121  y^+_{lam} | estimated intersection of the viscous and inertial sublayers
122  \endvartable
123 
124 
125  - \c max (PH:Eq. 27):
126 
127  \f[
128  \omega = max(\omega_{vis}, \omega_{log})
129  \f]
130 
131 
132  - \c binomial2 (ME:Eq. 15) (default):
133 
134  \f[
135  \omega = \sqrt{(\omega_{vis})^2 + (\omega_{log})^2}
136  \f]
137 
138 
139  - \c binomial:
140 
141  \f[
142  \omega = ((\omega_{vis})^n + (\omega_{log})^n)^{1/n}
143  \f]
144 
145  where
146  \vartable
147  n | Binomial blending exponent
148  \endvartable
149 
150 
151  - \c exponential (PH:Eq. 32):
152 
153  \f[
154  \omega = \omega_{vis} \exp[-\Gamma] + \omega_{log} \exp[-1/\Gamma]
155  \f]
156 
157  where (PH:Eq. 31)
158  \vartable
159  \Gamma | Blending expression
160  \Gamma | \f$0.01 (y^+)^4 / (1.0 + 5.0 y^+)\f$
161  \endvartable
162 
163 
164  - \c tanh (KAS:Eqs. 33-34):
165 
166  \f[
167  \omega = \phi \omega_{b1} + (1 - \phi)\omega_{b2}
168  \f]
169 
170  where
171  \vartable
172  \phi | \f$tanh((y^+/10)^4)\f$
173  \omega_{b1} | \f$\omega_{vis} + \omega_{log}\f$
174  \omega_{b2} | \f$(\omega_{vis}^{1.2} + \omega_{log}^1.2)^{1/1.2}\f$
175  \endvartable
176 
177 
178  \c G predictions for the viscous and inertial sublayers are blended
179  in a stepwise manner, and \c G below \f$y^+_{lam}\f$ (i.e. in the viscous
180  sublayer) is presumed to be zero.
181 
182 Note
183  - The coefficients \c Cmu, \c kappa, and \c E are obtained from
184  the specified \c nutWallFunction in order to ensure that each patch
185  possesses the same set of values for these coefficients.
186  - The reason why \c binomial2 and \c binomial blending methods exist at
187  the same time is to ensure the bitwise regression with the previous
188  versions since \c binomial2 and \c binomial with \c n=2 will yield
189  slightly different output due to the miniscule differences in the
190  implementation of the basic functions (i.e. \c pow, \c sqrt, \c sqr).
191 
192 See also
193  - Foam::epsilonWallFunctionFvPatchScalarField
194 
195 SourceFiles
196  omegaWallFunctionFvPatchScalarField.C
197 
198 \*---------------------------------------------------------------------------*/
199 
200 #ifndef omegaWallFunctionFvPatchScalarField_H
201 #define omegaWallFunctionFvPatchScalarField_H
202 
203 #include "fixedValueFvPatchField.H"
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 namespace Foam
208 {
209 
210 class turbulenceModel;
211 
212 /*---------------------------------------------------------------------------*\
213  Class omegaWallFunctionFvPatchScalarField Declaration
214 \*---------------------------------------------------------------------------*/
215 
216 class omegaWallFunctionFvPatchScalarField
217 :
218  public fixedValueFvPatchField<scalar>
219 {
220  // Private Enumerations
221 
222  //- Options for the blending treatment of viscous and inertial sublayers
223  enum blendingType
224  {
225  STEPWISE,
226  MAX,
227  BINOMIAL2,
228  BINOMIAL,
229  EXPONENTIAL,
230  TANH
231  };
232 
233  //- Names for blendingType
234  static const Enum<blendingType> blendingTypeNames;
235 
236 
237  // Private Data
238 
239  //- Blending treatment (default = blendingType::BINOMIAL2)
240  enum blendingType blending_;
241 
242  //- Blending exponent being used when
243  //- blendingType is blendingType::BINOMIAL (default = 2)
244  const scalar n_;
245 
246 
247 protected:
248 
249  // Protected Data
250 
251  //- Tolerance used in weighted calculations
252  static scalar tolerance_;
253 
254  //- Deprecated(2019-11) Blending switch
255  // \deprecated(2019-11) - use blending:: options
256  bool blended_;
257 
258  //- Initialised flag
259  bool initialised_;
260 
261  //- Master patch ID
262  label master_;
263 
264  //- beta1 coefficient
265  scalar beta1_;
266 
267  //- Local copy of turbulence G field
268  scalarField G_;
269 
270  //- Local copy of turbulence omega field
272 
273  //- List of averaging corner weights
274  List<List<scalar>> cornerWeights_;
275 
276 
277  // Protected Member Functions
278 
279  //- Set the master patch - master is responsible for updating all
280  //- wall function patches
281  virtual void setMaster();
282 
283  //- Create the averaging weights for cells which are bounded by
284  //- multiple wall function faces
285  virtual void createAveragingWeights();
286 
287  //- Helper function to return non-const access to an omega patch
289  (
290  const label patchi
291  );
292 
293  //- Main driver to calculate the turbulence fields
294  virtual void calculateTurbulenceFields
295  (
297  scalarField& G0,
298  scalarField& omega0
299  );
300 
301  //- Calculate the omega and G
302  virtual void calculate
303  (
305  const List<scalar>& cornerWeights,
306  const fvPatch& patch,
307  scalarField& G,
309  );
310 
311  //- Return non-const access to the master patch ID
312  virtual label& master()
313  {
314  return master_;
315  }
316 
317 
318 public:
319 
320  //- Runtime type information
321  TypeName("omegaWallFunction");
322 
323 
324  // Constructors
325 
326  //- Construct from patch and internal field
328  (
329  const fvPatch&,
331  );
332 
333  //- Construct from patch, internal field and dictionary
335  (
336  const fvPatch&,
338  const dictionary&
339  );
340 
341  //- Construct by mapping given
342  //- omegaWallFunctionFvPatchScalarField
343  //- onto a new patch
345  (
347  const fvPatch&,
350  );
351 
352  //- Construct as copy
354  (
356  );
357 
358  //- Construct and return a clone
359  virtual tmp<fvPatchScalarField> clone() const
360  {
362  (
364  );
365  }
366 
367  //- Construct as copy setting internal field reference
369  (
372  );
373 
374  //- Construct and return a clone setting internal field reference
376  (
378  ) const
379  {
381  (
383  );
384  }
385 
386  //- Destructor
388 
389 
390  // Member Functions
391 
392  // Access
393 
394  //- Return non-const access to the master's G field
395  scalarField& G(bool init = false);
396 
397  //- Return non-const access to the master's omega field
398  scalarField& omega(bool init = false);
399 
400 
401  // Evaluation functions
402 
403  //- Update the coefficients associated with the patch field
404  virtual void updateCoeffs();
405 
406  //- Update the coefficients associated with the patch field
407  virtual void updateWeightedCoeffs(const scalarField& weights);
408 
409  //- Manipulate matrix
410  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
411 
412  //- Manipulate matrix with given weights
413  virtual void manipulateMatrix
414  (
415  fvMatrix<scalar>& matrix,
416  const scalarField& weights
417  );
418 
419 
420  // I-O
421 
422  //- Write
423  virtual void write(Ostream&) const;
424 };
425 
426 
427 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
428 
429 } // End namespace Foam
430 
431 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
432 
433 #endif
434 
435 // ************************************************************************* //
Foam::omegaWallFunctionFvPatchScalarField::createAveragingWeights
virtual void createAveragingWeights()
Definition: omegaWallFunctionFvPatchScalarField.C:85
Foam::omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: omegaWallFunctionFvPatchScalarField.C:315
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::omegaWallFunctionFvPatchScalarField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
Definition: omegaWallFunctionFvPatchScalarField.C:591
Foam::omegaWallFunctionFvPatchScalarField::calculate
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.
Definition: omegaWallFunctionFvPatchScalarField.C:195
Foam::Enum< blendingType >
Foam::omegaWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the specific dissipation rate,...
Definition: omegaWallFunctionFvPatchScalarField.H:291
Foam::omegaWallFunctionFvPatchScalarField::updateWeightedCoeffs
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
Definition: omegaWallFunctionFvPatchScalarField.C:534
turbulence
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::omegaWallFunctionFvPatchScalarField::~omegaWallFunctionFvPatchScalarField
virtual ~omegaWallFunctionFvPatchScalarField()=default
Destructor.
Foam::omegaWallFunctionFvPatchScalarField::G_
scalarField G_
Local copy of turbulence G field.
Definition: omegaWallFunctionFvPatchScalarField.H:343
Foam::omegaWallFunctionFvPatchScalarField::beta1_
scalar beta1_
beta1 coefficient
Definition: omegaWallFunctionFvPatchScalarField.H:340
Foam::omegaWallFunctionFvPatchScalarField::cornerWeights_
List< List< scalar > > cornerWeights_
List of averaging corner weights.
Definition: omegaWallFunctionFvPatchScalarField.H:349
Foam::omegaWallFunctionFvPatchScalarField::master_
label master_
Master patch ID.
Definition: omegaWallFunctionFvPatchScalarField.H:337
Foam::fixedValueFvPatchField
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Definition: fixedValueFvPatchField.H:80
Foam::omegaWallFunctionFvPatchScalarField::omega_
scalarField omega_
Local copy of turbulence omega field.
Definition: omegaWallFunctionFvPatchScalarField.H:346
Foam::omegaWallFunctionFvPatchScalarField::omega
scalarField & omega(bool init=false)
Return non-const access to the master's omega field.
Definition: omegaWallFunctionFvPatchScalarField.C:470
Foam::omegaWallFunctionFvPatchScalarField::TypeName
TypeName("omegaWallFunction")
Runtime type information.
Foam::Field< scalar >
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:63
Foam::constant::electromagnetic::G0
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
Foam::omegaWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: omegaWallFunctionFvPatchScalarField.C:650
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
init
mesh init(true)
Foam::omegaWallFunctionFvPatchScalarField::initialised_
bool initialised_
Initialised flag.
Definition: omegaWallFunctionFvPatchScalarField.H:334
Foam::omegaWallFunctionFvPatchScalarField::omegaPatch
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
Definition: omegaWallFunctionFvPatchScalarField.C:145
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::omegaWallFunctionFvPatchScalarField::calculateTurbulenceFields
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence fields.
Definition: omegaWallFunctionFvPatchScalarField.C:162
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::omegaWallFunctionFvPatchScalarField::blended_
bool blended_
Deprecated(2019-11) Blending switch.
Definition: omegaWallFunctionFvPatchScalarField.H:331
Foam::omegaWallFunctionFvPatchScalarField::clone
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
Definition: omegaWallFunctionFvPatchScalarField.H:434
Foam::omegaWallFunctionFvPatchScalarField::setMaster
virtual void setMaster()
Definition: omegaWallFunctionFvPatchScalarField.C:55
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::omegaWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: omegaWallFunctionFvPatchScalarField.C:488
Foam::omegaWallFunctionFvPatchScalarField::master
virtual label & master()
Return non-const access to the master patch ID.
Definition: omegaWallFunctionFvPatchScalarField.H:387
fixedValueFvPatchField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::fvPatchField< scalar >::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::omegaWallFunctionFvPatchScalarField::tolerance_
static scalar tolerance_
Tolerance used in weighted calculations.
Definition: omegaWallFunctionFvPatchScalarField.H:327
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::omegaWallFunctionFvPatchScalarField::G
scalarField & G(bool init=false)
Return non-const access to the master's G field.
Definition: omegaWallFunctionFvPatchScalarField.C:451