epsilonWallFunctionFvPatchScalarField.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::epsilonWallFunctionFvPatchScalarField
29 
30 Group
31  grpWallFunctions
32 
33 Description
34  This boundary condition provides a wall constraint on the turbulent kinetic
35  energy dissipation rate, i.e. \c epsilon, and the turbulent kinetic
36  energy production contribution, i.e. \c G, for low- and high-Reynolds
37  number 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  \endverbatim
54 
55 Usage
56  Example of the boundary condition specification:
57  \verbatim
58  <patchName>
59  {
60  // Mandatory entries (unmodifiable)
61  type epsilonWallFunction;
62 
63  // Optional entries (unmodifiable)
64  lowReCorrection false;
65  blending stepwise;
66  n 2.0;
67 
68  // Optional (inherited) entries
69  ...
70  }
71  \endverbatim
72 
73  where the entries mean:
74  \table
75  Property | Description | Type | Req'd | Dflt
76  type | Type name: epsilonWallFunction | word | yes | -
77  lowReCorrection | Flag: apply low-Re correction | bool | no | false
78  blending | Viscous/inertial sublayer blending method <!--
79  --> | word | no | stepwise
80  n | Binomial blending exponent | scalar | no | 2.0
81  \endtable
82 
83  The inherited entries are elaborated in:
84  - \link fixedValueFvPatchField.H \endlink
85  - \link nutWallFunctionFvPatchScalarField.H \endlink
86 
87  Options for the \c blending entry:
88  \verbatim
89  stepwise | Stepwise switch (discontinuous)
90  max | Maximum value switch (discontinuous)
91  binomial | Binomial blending (smooth)
92  exponential | Exponential blending (smooth)
93  \endverbatim
94 
95  wherein \c epsilon predictions for the viscous and inertial sublayers are
96  blended according to the following expressions:
97 
98  - \c stepwise (default):
99 
100  \f[
101  \epsilon = \epsilon_{vis} \qquad if \quad y^+ < y^+_{lam}
102  \f]
103  \f[
104  \epsilon = \epsilon_{log} \qquad if \quad y^+ >= y^+_{lam}
105  \f]
106 
107  where
108  \vartable
109  \epsilon | \f$\epsilon\f$ at \f$y^+\f$
110  \epsilon_{vis} | \f$\epsilon\f$ computed by viscous subl. assumptions
111  \epsilon_{log} | \f$\epsilon\f$ computed by inertial subl. assumptions
112  y^+ | estimated wall-normal height of the cell centre in wall units
113  y^+_{lam} | estimated intersection of the viscous and inertial sublayers
114  \endvartable
115 
116 
117  - \c max (PH:Eq. 27):
118 
119  \f[
120  \epsilon = max(\epsilon_{vis}, \epsilon_{log})
121  \f]
122 
123 
124  - \c binomial (ME:Eqs. 15-16):
125 
126  \f[
127  \epsilon = ((\epsilon_{vis})^n + (\epsilon_{log})^n)^{1/n}
128  \f]
129  where
130  \vartable
131  n | Binomial blending exponent
132  \endvartable
133 
134 
135  - \c exponential (PH:Eq. 32):
136 
137  \f[
138  \epsilon = \epsilon_{vis} \exp[-\Gamma] +\epsilon_{log} \exp[-1/\Gamma]
139  \f]
140  where (PH:p. 193)
141  \vartable
142  \Gamma_\epsilon | \f$\Gamma = 0.001 (y^+)^4 / (1.0 + y^+)\f$
143  \Gamma_G | \f$\Gamma = 0.01 (y^+)^4 / (1.0 + 5.0 y^+)\f$
144  \Gamma_\epsilon | Blending expression for \f$\epsilon\f$
145  \Gamma_G | Blending expression for \f$G\f$
146  \endvartable
147 
148 
149  \c G predictions for the viscous and inertial sublayers are blended
150  in a stepwise manner, and \c G below \f$y^+_{lam}\f$ (i.e. in the viscous
151  sublayer) is presumed to be zero.
152 
153 Note
154  - The coefficients \c Cmu, \c kappa, and \c E are obtained from
155  the specified \c nutWallFunction in order to ensure that each patch
156  possesses the same set of values for these coefficients.
157  - \c lowReCorrection operates with only \c stepwise blending treatment to
158  ensure the backward compatibility.
159  - If \c lowReCorrection is \c on, \c stepwise blending treatment is fully
160  active.
161  - If \c lowReCorrection is \c off, only the inertial sublayer prediction
162  is used in the wall function, hence high-Re mode operation.
163 
164 See also
165  - Foam::omegaWallFunctionFvPatchScalarField
166 
167 SourceFiles
168  epsilonWallFunctionFvPatchScalarField.C
169 
170 \*---------------------------------------------------------------------------*/
171 
172 #ifndef epsilonWallFunctionFvPatchScalarField_H
173 #define epsilonWallFunctionFvPatchScalarField_H
174 
175 #include "fixedValueFvPatchField.H"
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 namespace Foam
180 {
181 
182 class turbulenceModel;
183 
184 /*---------------------------------------------------------------------------*\
185  Class epsilonWallFunctionFvPatchScalarField Declaration
186 \*---------------------------------------------------------------------------*/
187 
188 class epsilonWallFunctionFvPatchScalarField
189 :
190  public fixedValueFvPatchField<scalar>
191 {
192  // Private Enumerations
193 
194  //- Options for the blending treatment of viscous and inertial sublayers
195  enum blendingType
196  {
197  STEPWISE,
198  MAX,
199  BINOMIAL,
200  EXPONENTIAL
201  };
202 
203  //- Names for blendingType
204  static const Enum<blendingType> blendingTypeNames;
205 
206 
207  // Private Data
208 
209  //- Blending treatment (default = blendingType::STEPWISE)
210  const enum blendingType blending_;
211 
212  //- Binomial blending exponent being used when
213  //- blendingType is blendingType::BINOMIAL (default = 2)
214  const scalar n_;
215 
216 
217 protected:
218 
219  // Protected Data
220 
221  //- Tolerance used in weighted calculations
222  static scalar tolerance_;
223 
224  //- Apply low-Re correction term (default = no)
225  const bool lowReCorrection_;
226 
227  //- Initialised flag
228  bool initialised_;
229 
230  //- Master patch ID
231  label master_;
232 
233  //- Local copy of turbulence G field
234  scalarField G_;
235 
236  //- Local copy of turbulence epsilon field
238 
239  //- List of averaging corner weights
240  List<List<scalar>> cornerWeights_;
241 
242 
243  // Protected Member Functions
244 
245  //- Set the master patch - master is responsible for updating all
246  //- wall function patches
247  virtual void setMaster();
248 
249  //- Create the averaging weights for cells which are bounded by
250  //- multiple wall function faces
251  virtual void createAveragingWeights();
252 
253  //- Helper function to return non-const access to an epsilon patch
255  (
256  const label patchi
257  );
258 
259  //- Main driver to calculate the turbulence fields
260  virtual void calculateTurbulenceFields
261  (
263  scalarField& G0,
265  );
266 
267  //- Calculate the epsilon and G
268  virtual void calculate
269  (
271  const List<scalar>& cornerWeights,
272  const fvPatch& patch,
273  scalarField& G,
275  );
276 
277  //- Return non-const access to the master patch ID
278  virtual label& master()
279  {
280  return master_;
281  }
282 
283 
284 public:
285 
286  //- Runtime type information
287  TypeName("epsilonWallFunction");
288 
289 
290  // Constructors
291 
292  //- Construct from patch and internal field
294  (
295  const fvPatch&,
297  );
298 
299  //- Construct from patch, internal field and dictionary
301  (
302  const fvPatch&,
304  const dictionary&
305  );
306 
307  //- Construct by mapping given
308  //- epsilonWallFunctionFvPatchScalarField
309  //- onto a new patch
311  (
313  const fvPatch&,
315  const fvPatchFieldMapper&
316  );
317 
318  //- Construct as copy
320  (
322  );
323 
324  //- Construct and return a clone
325  virtual tmp<fvPatchScalarField> clone() const
326  {
328  (
330  );
331  }
332 
333  //- Construct as copy setting internal field reference
335  (
338  );
339 
340  //- Construct and return a clone setting internal field reference
342  (
344  ) const
345  {
347  (
349  );
350  }
351 
352 
353  //- Destructor
354  virtual ~epsilonWallFunctionFvPatchScalarField() = default;
355 
356 
357  // Member Functions
358 
359  // Access
360 
361  //- Return non-const access to the master's G field
362  scalarField& G(bool init = false);
363 
364  //- Return non-const access to the master's epsilon field
365  scalarField& epsilon(bool init = false);
366 
367 
368  // Evaluation functions
369 
370  //- Update the coefficients associated with the patch field
371  virtual void updateCoeffs();
372 
373  //- Update the coefficients associated with the patch field
374  virtual void updateWeightedCoeffs(const scalarField& weights);
375 
376  //- Manipulate matrix
377  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
378 
379  //- Manipulate matrix with given weights
380  virtual void manipulateMatrix
381  (
382  fvMatrix<scalar>& matrix,
383  const scalarField& weights
384  );
385 
386 
387  // I-O
388 
389  //- Write
390  virtual void write(Ostream&) const;
391 };
392 
393 
394 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395 
396 } // End namespace Foam
397 
398 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
399 
400 #endif
401 
402 // ************************************************************************* //
Foam::epsilonWallFunctionFvPatchScalarField::initialised_
bool initialised_
Initialised flag.
Definition: epsilonWallFunctionFvPatchScalarField.H:299
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::Enum< blendingType >
Foam::epsilonWallFunctionFvPatchScalarField::master
virtual label & master()
Return non-const access to the master patch ID.
Definition: epsilonWallFunctionFvPatchScalarField.H:349
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::epsilonWallFunctionFvPatchScalarField::updateWeightedCoeffs
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
Definition: epsilonWallFunctionFvPatchScalarField.C:505
Foam::epsilonWallFunctionFvPatchScalarField::tolerance_
static scalar tolerance_
Tolerance used in weighted calculations.
Definition: epsilonWallFunctionFvPatchScalarField.H:293
Foam::epsilonWallFunctionFvPatchScalarField::createAveragingWeights
virtual void createAveragingWeights()
Definition: epsilonWallFunctionFvPatchScalarField.C:84
Foam::epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: epsilonWallFunctionFvPatchScalarField.C:305
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::constant::electromagnetic::epsilon0
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
Foam::epsilonWallFunctionFvPatchScalarField::epsilon
scalarField & epsilon(bool init=false)
Return non-const access to the master's epsilon field.
Definition: epsilonWallFunctionFvPatchScalarField.C:441
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::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::epsilonWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent kinetic energy dissipation rate,...
Definition: epsilonWallFunctionFvPatchScalarField.H:259
init
mesh init(true)
Foam::epsilonWallFunctionFvPatchScalarField::calculate
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
Definition: epsilonWallFunctionFvPatchScalarField.C:195
Foam::epsilonWallFunctionFvPatchScalarField::lowReCorrection_
const bool lowReCorrection_
Apply low-Re correction term (default = no)
Definition: epsilonWallFunctionFvPatchScalarField.H:296
Foam::epsilonWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: epsilonWallFunctionFvPatchScalarField.C:621
Foam::epsilonWallFunctionFvPatchScalarField::G
scalarField & G(bool init=false)
Return non-const access to the master's G field.
Definition: epsilonWallFunctionFvPatchScalarField.C:422
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::epsilonWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: epsilonWallFunctionFvPatchScalarField.C:459
Foam::epsilonWallFunctionFvPatchScalarField::G_
scalarField G_
Local copy of turbulence G field.
Definition: epsilonWallFunctionFvPatchScalarField.H:305
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::epsilonWallFunctionFvPatchScalarField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
Definition: epsilonWallFunctionFvPatchScalarField.C:562
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::epsilonWallFunctionFvPatchScalarField::calculateTurbulenceFields
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
Definition: epsilonWallFunctionFvPatchScalarField.C:162
Foam::epsilonWallFunctionFvPatchScalarField::epsilonPatch
virtual epsilonWallFunctionFvPatchScalarField & epsilonPatch(const label patchi)
Helper function to return non-const access to an epsilon patch.
Definition: epsilonWallFunctionFvPatchScalarField.C:145
Foam::epsilonWallFunctionFvPatchScalarField::epsilon_
scalarField epsilon_
Local copy of turbulence epsilon field.
Definition: epsilonWallFunctionFvPatchScalarField.H:308
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::epsilonWallFunctionFvPatchScalarField::TypeName
TypeName("epsilonWallFunction")
Runtime type information.
Foam::epsilonWallFunctionFvPatchScalarField::master_
label master_
Master patch ID.
Definition: epsilonWallFunctionFvPatchScalarField.H:302
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::epsilonWallFunctionFvPatchScalarField::setMaster
virtual void setMaster()
Definition: epsilonWallFunctionFvPatchScalarField.C:54
fixedValueFvPatchField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::epsilonWallFunctionFvPatchScalarField::clone
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
Definition: epsilonWallFunctionFvPatchScalarField.H:396
Foam::fvPatchField< scalar >::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
Foam::epsilonWallFunctionFvPatchScalarField::~epsilonWallFunctionFvPatchScalarField
virtual ~epsilonWallFunctionFvPatchScalarField()=default
Destructor.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::epsilonWallFunctionFvPatchScalarField::cornerWeights_
List< List< scalar > > cornerWeights_
List of averaging corner weights.
Definition: epsilonWallFunctionFvPatchScalarField.H:311
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54