nutkWallFunctionFvPatchScalarField.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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
33#include "wallFvPatch.H"
35
36// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37
39calcNut() const
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 scalarField& y = turbModel.y()[patchi];
53
54 const tmp<volScalarField> tk = turbModel.k();
55 const volScalarField& k = tk();
56
57 const tmp<scalarField> tnuw = turbModel.nu(patchi);
58 const scalarField& nuw = tnuw();
59
60 const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
61 const scalar kappa = wallCoeffs_.kappa();
62 const scalar E = wallCoeffs_.E();
63 const scalar yPlusLam = wallCoeffs_.yPlusLam();
64
65 auto tnutw = tmp<scalarField>::New(patch().size(), Zero);
66 auto& nutw = tnutw.ref();
67
68 forAll(nutw, facei)
69 {
70 const label celli = patch().faceCells()[facei];
71
72 const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
73
74 // Viscous sublayer contribution
75 const scalar nutVis = 0;
76
77 // Inertial sublayer contribution
78 const scalar nutLog =
79 nuw[facei]*(yPlus*kappa/log(max(E*yPlus, 1 + 1e-4)) - 1.0);
80
81 switch (blender_)
82 {
83 case blenderType::STEPWISE:
84 {
85 if (yPlus > yPlusLam)
86 {
87 nutw[facei] = nutLog;
88 }
89 else
90 {
91 nutw[facei] = nutVis;
92 }
93 break;
94 }
95
96 case blenderType::MAX:
97 {
98 // (PH:Eq. 27)
99 nutw[facei] = max(nutVis, nutLog);
100 break;
101 }
102
103 case blenderType::BINOMIAL:
104 {
105 // (ME:Eqs. 15-16)
106 nutw[facei] =
107 pow
108 (
109 pow(nutVis, n_) + pow(nutLog, n_),
110 scalar(1)/n_
111 );
112 break;
113 }
114
115 case blenderType::EXPONENTIAL:
116 {
117 // (PH:Eq. 31)
118 const scalar Gamma = 0.01*pow4(yPlus)/(1 + 5*yPlus);
119 const scalar invGamma = scalar(1)/(Gamma + ROOTVSMALL);
120
121 nutw[facei] = nutVis*exp(-Gamma) + nutLog*exp(-invGamma);
122 break;
123 }
124
125 case blenderType::TANH:
126 {
127 // (KAS:Eqs. 33-34)
128 const scalar phiTanh = tanh(pow4(0.1*yPlus));
129 const scalar b1 = nutVis + nutLog;
130 const scalar b2 =
131 pow(pow(nutVis, 1.2) + pow(nutLog, 1.2), 1.0/1.2);
132
133 nutw[facei] = phiTanh*b1 + (1 - phiTanh)*b2;
134 break;
135 }
136 }
137 }
138
139 return tnutw;
140}
141
142
144(
145 Ostream& os
146) const
147{
149}
150
151
152// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
153
155(
156 const fvPatch& p,
158)
159:
162{}
163
164
166(
168 const fvPatch& p,
170 const fvPatchFieldMapper& mapper
171)
172:
173 nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
175{}
176
177
179(
180 const fvPatch& p,
182 const dictionary& dict
183)
184:
186 wallFunctionBlenders(dict, blenderType::STEPWISE, scalar(4))
187{}
188
189
191(
193)
194:
197{}
198
199
201(
204)
205:
208{}
209
210
211// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212
214yPlus() const
215{
216 const label patchi = patch().index();
217
218 const auto& turbModel = db().lookupObject<turbulenceModel>
219 (
221 (
223 internalField().group()
224 )
225 );
226
227 const scalarField& y = turbModel.y()[patchi];
228
229 tmp<volScalarField> tk = turbModel.k();
230 const volScalarField& k = tk();
231 tmp<scalarField> tkwc = k.boundaryField()[patchi].patchInternalField();
232 const scalarField& kwc = tkwc();
233
234 tmp<scalarField> tnuw = turbModel.nu(patchi);
235 const scalarField& nuw = tnuw();
236
237 tmp<scalarField> tnuEff = turbModel.nuEff(patchi);
238 const scalarField& nuEff = tnuEff();
239
240 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
241 const scalarField magGradUw(mag(Uw.snGrad()));
242
243 const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
244 const scalar yPlusLam = wallCoeffs_.yPlusLam();
245
246 auto tyPlus = tmp<scalarField>::New(patch().size(), Zero);
247 auto& yPlus = tyPlus.ref();
248
249 forAll(yPlus, facei)
250 {
251 // inertial sublayer
252 yPlus[facei] = Cmu25*y[facei]*sqrt(kwc[facei])/nuw[facei];
253
254 if (yPlusLam > yPlus[facei])
255 {
256 // viscous sublayer
257 yPlus[facei] =
258 y[facei]*sqrt(nuEff[facei]*magGradUw[facei])/nuw[facei];
259 }
260 }
261
262 return tyPlus;
263}
264
265
267(
268 Ostream& os
269) const
270{
272 writeLocalEntries(os);
273 writeEntry("value", os);
274}
275
276
277// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278
279namespace Foam
280{
282 (
285 );
286}
287
288
289// ************************************************************************* //
scalar y
label k
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:229
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
This boundary condition provides a wall function for the turbulent viscosity (i.e....
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
A class for managing temporary objects.
Definition: tmp.H:65
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
The class wallFunctionBlenders is a base class that hosts common entries for various derived wall-fun...
void writeEntries(Ostream &) const
Write wall-function blending data as dictionary entries.
blenderType
Options for the blending treatment of viscous and inertial sublayers.
enum blenderType blender_
Blending treatment.
scalar kappa() const noexcept
Return the object: kappa.
scalar E() const noexcept
Return the object: E.
scalar yPlusLam() const noexcept
Return the object: yPlusLam.
scalar Cmu() const noexcept
Return the object: Cmu.
U
Definition: pEqn.H:72
volScalarField & p
scalar yPlus
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333