nutUWallFunctionFvPatchScalarField.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 OpenFOAM Foundation
9  Copyright (C) 2019-2021 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 "turbulenceModel.H"
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
34 
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
40 {
41  const label patchi = patch().index();
42 
43  const auto& turbModel = db().lookupObject<turbulenceModel>
44  (
46  (
48  internalField().group()
49  )
50  );
51 
52  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
53  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
54 
55  const tmp<scalarField> tnuw = turbModel.nu(patchi);
56  const scalarField& nuw = tnuw();
57 
58  const scalar kappa = wallCoeffs_.kappa();
59  const scalar E = wallCoeffs_.E();
60  const scalar yPlusLam = wallCoeffs_.yPlusLam();
61 
63  const scalarField& yPlus = tyPlus();
64 
65  auto tnutw = tmp<scalarField>::New(patch().size(), Zero);
66  auto& nutw = tnutw.ref();
67 
68  forAll(yPlus, facei)
69  {
70  // Viscous sublayer contribution
71  const scalar nutVis = 0;
72 
73  // Inertial sublayer contribution
74  const scalar nutLog =
75  nuw[facei]
76  *(yPlus[facei]*kappa/log(max(E*yPlus[facei], 1 + 1e-4)) - 1.0);
77 
78  switch (blender_)
79  {
80  case blenderType::STEPWISE:
81  {
82  if (yPlus[facei] > yPlusLam)
83  {
84  nutw[facei] = nutLog;
85  }
86  else
87  {
88  nutw[facei] = nutVis;
89  }
90  break;
91  }
92 
93  case blenderType::MAX:
94  {
95  // (PH:Eq. 27)
96  nutw[facei] = max(nutVis, nutLog);
97  break;
98  }
99 
100  case blenderType::BINOMIAL:
101  {
102  // (ME:Eqs. 15-16)
103  nutw[facei] =
104  pow
105  (
106  pow(nutVis, n_) + pow(nutLog, n_),
107  scalar(1)/n_
108  );
109  break;
110  }
111 
112  case blenderType::EXPONENTIAL:
113  {
114  // (PH:Eq. 31)
115  const scalar Gamma =
116  0.01*pow4(yPlus[facei])/(1 + 5*yPlus[facei]);
117  const scalar invGamma = scalar(1)/(Gamma + ROOTVSMALL);
118 
119  nutw[facei] = nutVis*exp(-Gamma) + nutLog*exp(-invGamma);
120  break;
121  }
122 
123  case blenderType::TANH:
124  {
125  // (KAS:Eqs. 33-34)
126  const scalar phiTanh = tanh(pow4(0.1*yPlus[facei]));
127  const scalar b1 = nutVis + nutLog;
128  const scalar b2 =
129  pow(pow(nutVis, 1.2) + pow(nutLog, 1.2), 1.0/1.2);
130 
131  nutw[facei] = phiTanh*b1 + (1 - phiTanh)*b2;
132  break;
133  }
134  }
135  }
137  return tnutw;
138 }
139 
140 
143 (
144  const scalarField& magUp
145 ) const
146 {
147  const label patchi = patch().index();
148 
149  const auto& turbModel = db().lookupObject<turbulenceModel>
150  (
152  (
154  internalField().group()
155  )
156  );
157 
158  const scalarField& y = turbModel.y()[patchi];
159 
160  const tmp<scalarField> tnuw = turbModel.nu(patchi);
161  const scalarField& nuw = tnuw();
162 
163  const scalar kappa = wallCoeffs_.kappa();
164  const scalar E = wallCoeffs_.E();
165  const scalar yPlusLam = wallCoeffs_.yPlusLam();
166 
167  auto tyPlus = tmp<scalarField>::New(patch().size(), Zero);
168  auto& yPlus = tyPlus.ref();
169 
170  forAll(yPlus, facei)
171  {
172  const scalar kappaRe = kappa*magUp[facei]*y[facei]/nuw[facei];
173 
174  scalar yp = yPlusLam;
175  const scalar ryPlusLam = 1.0/yp;
176 
177  int iter = 0;
178  scalar yPlusLast = 0.0;
179 
180  do
181  {
182  yPlusLast = yp;
183  yp = (kappaRe + yp)/(1.0 + log(E*yp));
184 
185  } while (mag(ryPlusLam*(yp - yPlusLast)) > 0.01 && ++iter < 10 );
186 
187  yPlus[facei] = max(scalar(0), yp);
188  }
189 
190  return tyPlus;
191 }
192 
193 
195 (
196  Ostream& os
197 ) const
198 {
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
204 
206 (
207  const fvPatch& p,
209 )
210 :
213 {}
214 
215 
217 (
219  const fvPatch& p,
221  const fvPatchFieldMapper& mapper
222 )
223 :
224  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
226 {}
227 
228 
230 (
231  const fvPatch& p,
233  const dictionary& dict
234 )
235 :
237  wallFunctionBlenders(dict, blenderType::STEPWISE, scalar(4))
238 {}
239 
240 
242 (
244 )
245 :
247  wallFunctionBlenders(sawfpsf)
248 {}
249 
250 
252 (
253  const nutUWallFunctionFvPatchScalarField& sawfpsf,
255 )
256 :
259 {}
260 
261 
262 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
263 
266 {
267  const label patchi = patch().index();
268 
269  const auto& turbModel = db().lookupObject<turbulenceModel>
270  (
272  (
274  internalField().group()
275  )
276  );
277 
278  const scalarField& y = turbModel.y()[patchi];
279 
280  tmp<scalarField> tnuw = turbModel.nu(patchi);
281  const scalarField& nuw = tnuw();
282 
283  tmp<scalarField> tnuEff = turbModel.nuEff(patchi);
284  const scalarField& nuEff = tnuEff();
285 
286  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
287  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
288  const scalarField magGradUw(mag(Uw.snGrad()));
289 
290  const scalar yPlusLam = wallCoeffs_.yPlusLam();
291 
292  tmp<scalarField> tyPlus = calcYPlus(magUp);
293  scalarField& yPlus = tyPlus.ref();
294 
295  forAll(yPlus, facei)
296  {
297  if (yPlusLam > yPlus[facei])
298  {
299  // viscous sublayer
300  yPlus[facei] =
301  y[facei]*sqrt(nuEff[facei]*magGradUw[facei])/nuw[facei];
302  }
303  }
304 
305  return tyPlus;
306 }
307 
308 
310 (
311  Ostream& os
312 ) const
313 {
315  writeLocalEntries(os);
316  writeEntry("value", os);
317 }
318 
319 
320 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
321 
322 namespace Foam
323 {
325  (
327  nutUWallFunctionFvPatchScalarField
328  );
329 }
330 
331 
332 // ************************************************************************* //
scalar n_
Binomial blending exponent being used when blenderType is blenderType::BINOMIAL.
dimensionedScalar tanh(const dimensionedScalar &ds)
dictionary dict
scalar yPlusLam() const noexcept
Return the object: yPlusLam.
virtual const volVectorField & U(const turbulenceModel &turb) const
Helper to return the velocity field either from the turbulence model (default) or the mesh database...
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:174
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
This boundary condition provides a wall function for the turbulent viscosity (i.e. nut) based on velocity (i.e. U) for low- and high-Reynolds number applications.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:210
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Abstract base class for turbulence models (RAS, LES and laminar).
The class wallFunctionBlenders is a base class that hosts common entries for various derived wall-fun...
Macros for easy insertion into run-time selection tables.
scalar magUp
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
constexpr const char *const group
Group name for atomic constants.
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dimensionedScalar exp(const dimensionedScalar &ds)
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
tmp< scalarField > calcYPlus(const scalarField &magUp) const
Calculate yPlus.
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:203
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:182
nutUWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
OBJstream os(runTime.globalPath()/outputName)
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
scalar yPlus
void writeEntries(Ostream &) const
Write wall-function blending data as dictionary entries.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
dimensionedScalar pow4(const dimensionedScalar &ds)
const std::string patch
OpenFOAM patch number as a std::string.
scalar E() const noexcept
Return the object: E.
enum blenderType blender_
Blending treatment.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
scalar kappa() const noexcept
Return the object: kappa.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
blenderType
Options for the blending treatment of viscous and inertial sublayers.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157