Loading...
Searching...
No Matches
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-2023 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"
32#include "volFields.H"
34
35// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36
39{
40 const label patchi = patch().index();
41
42 const scalar kappa = wallCoeffs_.kappa();
43 const scalar E = wallCoeffs_.E();
44 const scalar yPlusLam = wallCoeffs_.yPlusLam();
45
46 const auto& turbModel = db().lookupObject<turbulenceModel>
47 (
49 (
51 internalField().group()
52 )
53 );
54
55 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
56 const scalarField magUp(mag(Uw.patchInternalField() - Uw));
57
59 const scalarField& yPlus = tyPlus();
60
61 // Viscous sublayer contribution
62 const tmp<scalarField> tnutVis = turbModel.nu(patchi);
63 const scalarField& nutVis = tnutVis();
64
65 // Inertial sublayer contribution
66 const auto nutLog = [&](const label facei) -> scalar
67 {
68 const scalar yPlusFace = yPlus[facei];
69 return
70 (
71 nutVis[facei]*yPlusFace*kappa/log(max(E*yPlusFace, 1 + 1e-4))
72 );
73 };
74
75 auto tnutw = tmp<scalarField>::New(patch().size(), Zero);
76 auto& nutw = tnutw.ref();
77
78 switch (blender_)
79 {
80 case blenderType::STEPWISE:
81 {
82 forAll(nutw, facei)
83 {
84 if (yPlus[facei] > yPlusLam)
85 {
86 nutw[facei] = nutLog(facei);
87 }
88 else
89 {
90 nutw[facei] = nutVis[facei];
91 }
92 }
93 break;
94 }
95
96 case blenderType::MAX:
97 {
98 forAll(nutw, facei)
99 {
100 // (PH:Eq. 27)
101 nutw[facei] = max(nutVis[facei], nutLog(facei));
102 }
103 break;
104 }
105
106 case blenderType::BINOMIAL:
107 {
108 forAll(nutw, facei)
109 {
110 // (ME:Eqs. 15-16)
111 nutw[facei] =
112 pow
113 (
114 pow(nutVis[facei], n_) + pow(nutLog(facei), n_),
115 scalar(1)/n_
116 );
117 }
118 break;
119 }
120
121 case blenderType::EXPONENTIAL:
122 {
123 forAll(nutw, facei)
124 {
125 // (PH:Eq. 31)
126 const scalar Gamma =
127 0.01*pow4(yPlus[facei])/(1 + 5*yPlus[facei]);
128 const scalar invGamma = scalar(1)/(Gamma + ROOTVSMALL);
129
130 nutw[facei] =
131 nutVis[facei]*exp(-Gamma) + nutLog(facei)*exp(-invGamma);
132 }
133 break;
134 }
135
136 case blenderType::TANH:
137 {
138 forAll(nutw, facei)
139 {
140 // (KAS:Eqs. 33-34)
141 const scalar nutLogFace = nutLog(facei);
142 const scalar b1 = nutVis[facei] + nutLogFace;
143 const scalar b2 =
144 pow
145 (
146 pow(nutVis[facei], 1.2) + pow(nutLogFace, 1.2),
147 1.0/1.2
148 );
149 const scalar phiTanh = tanh(pow4(0.1*yPlus[facei]));
150
151 nutw[facei] = phiTanh*b1 + (1 - phiTanh)*b2;
152 }
153 break;
154 }
155 }
156
157 nutw -= nutVis;
158
159 return tnutw;
160}
161
162
165(
166 const scalarField& magUp
167) const
168{
169 const label patchi = patch().index();
170
171 const auto& turbModel = db().lookupObject<turbulenceModel>
172 (
174 (
176 internalField().group()
177 )
178 );
179
180 const scalarField& y = turbModel.y()[patchi];
181
182 const tmp<scalarField> tnuw = turbModel.nu(patchi);
183 const scalarField& nuw = tnuw();
184
185 const scalar kappa = wallCoeffs_.kappa();
186 const scalar E = wallCoeffs_.E();
187 const scalar yPlusLam = wallCoeffs_.yPlusLam();
188
189 auto tyPlus = tmp<scalarField>::New(patch().size(), Zero);
190 auto& yPlus = tyPlus.ref();
191
192 forAll(yPlus, facei)
193 {
194 const scalar kappaRe = kappa*magUp[facei]*y[facei]/nuw[facei];
195
196 scalar yp = yPlusLam;
197 const scalar ryPlusLam = 1.0/yp;
198
199 int iter = 0;
200 scalar yPlusLast = 0.0;
201
202 do
203 {
204 yPlusLast = yp;
205 yp = (kappaRe + yp)/(1.0 + log(E*yp));
206
207 } while (mag(ryPlusLam*(yp - yPlusLast)) > 0.01 && ++iter < 10 );
208
209 yPlus[facei] = max(scalar(0), yp);
210 }
211
212 return tyPlus;
213}
214
215
217(
218 Ostream& os
219) const
222}
223
224
225// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
226
228(
229 const fvPatch& p,
243 const fvPatchFieldMapper& mapper
255 const dictionary& dict
264(
277)
278:
281{}
282
283
284// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
285
288{
289 const label patchi = patch().index();
290
291 const auto& turbModel = db().lookupObject<turbulenceModel>
292 (
294 (
296 internalField().group()
297 )
298 );
299
300 const scalarField& y = turbModel.y()[patchi];
301
302 tmp<scalarField> tnuw = turbModel.nu(patchi);
303 const scalarField& nuw = tnuw();
304
305 tmp<scalarField> tnuEff = turbModel.nuEff(patchi);
306 const scalarField& nuEff = tnuEff();
307
308 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
309 const scalarField magUp(mag(Uw.patchInternalField() - Uw));
310 const scalarField magGradUw(mag(Uw.snGrad()));
311
312 const scalar yPlusLam = wallCoeffs_.yPlusLam();
313
314 tmp<scalarField> tyPlus = calcYPlus(magUp);
315 scalarField& yPlus = tyPlus.ref();
316
317 forAll(yPlus, facei)
318 {
319 if (yPlusLam > yPlus[facei])
320 {
321 // viscous sublayer
322 yPlus[facei] =
323 y[facei]*sqrt(nuEff[facei]*magGradUw[facei])/nuw[facei];
325 }
326
327 return tyPlus;
328}
329
330
332(
333 Ostream& os
334) const
335{
337 writeLocalEntries(os);
340
341
342// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
343
344namespace Foam
345{
347 (
349 nutUWallFunctionFvPatchScalarField
350 );
351}
352
353
354// ************************************************************************* //
scalar y
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:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
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.
nutUWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
tmp< scalarField > calcYPlus(const scalarField &magUp) const
Calculate yPlus.
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
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.
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
void writeEntries(Ostream &) const
Write wall-function blending data as dictionary entries.
wallFunctionBlenders()
Default construct with default coefficients.
blenderType
Options for the blending treatment of viscous and inertial sublayers.
@ STEPWISE
"Stepwise switch (discontinuous)"
enum blenderType blender_
Blending treatment.
scalar n_
Binomial blending exponent being used when blenderType is blenderType::BINOMIAL.
U
Definition pEqn.H:72
volScalarField & p
scalar yPlus
scalar magUp
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
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:40
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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:127
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299