Loading...
Searching...
No Matches
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-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"
33#include "wallFvPatch.H"
35
36// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37
39calcNut() const
40{
41 const label patchi = patch().index();
42
43 const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
44 const scalar kappa = wallCoeffs_.kappa();
45 const scalar E = wallCoeffs_.E();
46 const scalar yPlusLam = wallCoeffs_.yPlusLam();
47
48 const auto& turbModel = db().lookupObject<turbulenceModel>
49 (
51 (
53 internalField().group()
54 )
55 );
56
57 const labelUList& faceCells = patch().faceCells();
58
59 const scalarField& y = turbModel.y()[patchi];
60
61 const tmp<volScalarField> tk = turbModel.k();
62 const volScalarField& k = tk();
63
64 // Viscous sublayer contribution
65 const tmp<scalarField> tnutVis = turbModel.nu(patchi);
66 const scalarField& nutVis = tnutVis();
67
68 // Calculate y-plus
69 const auto yPlus = [&](const label facei) -> scalar
70 {
71 return (Cmu25*y[facei]*sqrt(k[faceCells[facei]])/nutVis[facei]);
72 };
73
74 // Inertial sublayer contribution
75 const auto nutLog = [&](const label facei) -> scalar
76 {
77 const scalar yPlusFace = yPlus(facei);
78 return
79 (
80 nutVis[facei]*yPlusFace*kappa
81 / log(max(E*yPlusFace, 1 + 1e-4))
82 );
83 };
84
85 auto tnutw = tmp<scalarField>::New(patch().size(), Zero);
86 auto& nutw = tnutw.ref();
87
88 switch (blender_)
89 {
90 case blenderType::STEPWISE:
91 {
92 forAll(nutw, facei)
93 {
94 if (yPlus(facei) > yPlusLam)
95 {
96 nutw[facei] = nutLog(facei);
97 }
98 else
99 {
100 nutw[facei] = nutVis[facei];
101 }
102 }
103 break;
104 }
105
106 case blenderType::MAX:
107 {
108 forAll(nutw, facei)
109 {
110 // (PH:Eq. 27)
111 nutw[facei] = max(nutVis[facei], nutLog(facei));
112 }
113 break;
114 }
115
116 case blenderType::BINOMIAL:
117 {
118 forAll(nutw, facei)
119 {
120 // (ME:Eqs. 15-16)
121 nutw[facei] =
122 pow
123 (
124 pow(nutVis[facei], n_) + pow(nutLog(facei), n_),
125 scalar(1)/n_
126 );
127 }
128 break;
129 }
130
131 case blenderType::EXPONENTIAL:
132 {
133 forAll(nutw, facei)
134 {
135 // (PH:Eq. 31)
136 const scalar yPlusFace = yPlus(facei);
137 const scalar Gamma = 0.01*pow4(yPlusFace)/(1 + 5*yPlusFace);
138 const scalar invGamma = scalar(1)/(Gamma + ROOTVSMALL);
139
140 nutw[facei] =
141 nutVis[facei]*exp(-Gamma) + nutLog(facei)*exp(-invGamma);
142 }
143 break;
144 }
145
146 case blenderType::TANH:
147 {
148 forAll(nutw, facei)
149 {
150 // (KAS:Eqs. 33-34)
151 const scalar nutLogFace = nutLog(facei);
152 const scalar b1 = nutVis[facei] + nutLogFace;
153 const scalar b2 =
154 pow
155 (
156 pow(nutVis[facei], 1.2) + pow(nutLogFace, 1.2),
157 1.0/1.2
158 );
159 const scalar phiTanh = tanh(pow4(0.1*yPlus(facei)));
160
161 nutw[facei] = phiTanh*b1 + (1 - phiTanh)*b2;
162 }
163 break;
164 }
165 }
167 nutw -= nutVis;
168
169 return tnutw;
170}
171
172
174(
175 Ostream& os
176) const
179}
180
181
182// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
183
185(
186 const fvPatch& p,
200 const fvPatchFieldMapper& mapper
212 const dictionary& dict
221(
234)
235:
238{}
239
240
241// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
242
244yPlus() const
245{
246 const label patchi = patch().index();
247
248 const auto& turbModel = db().lookupObject<turbulenceModel>
249 (
251 (
253 internalField().group()
254 )
255 );
256
257 const scalarField& y = turbModel.y()[patchi];
258
259 tmp<volScalarField> tk = turbModel.k();
260 const volScalarField& k = tk();
261 tmp<scalarField> tkwc = k.boundaryField()[patchi].patchInternalField();
262 const scalarField& kwc = tkwc();
263
264 tmp<scalarField> tnutVis = turbModel.nu(patchi);
265 const scalarField& nutVis = tnutVis();
266
267 tmp<scalarField> tnuEff = turbModel.nuEff(patchi);
268 const scalarField& nuEff = tnuEff();
269
270 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
271 const scalarField magGradUw(mag(Uw.snGrad()));
272
273 const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
274 const scalar yPlusLam = wallCoeffs_.yPlusLam();
275
276 auto tyPlus = tmp<scalarField>::New(patch().size(), Zero);
277 auto& yPlus = tyPlus.ref();
278
279 forAll(yPlus, facei)
280 {
281 // inertial sublayer
282 yPlus[facei] = Cmu25*y[facei]*sqrt(kwc[facei])/nutVis[facei];
283
284 if (yPlusLam > yPlus[facei])
285 {
286 // viscous sublayer
287 yPlus[facei] =
288 y[facei]*sqrt(nuEff[facei]*magGradUw[facei])/nutVis[facei];
290 }
291
292 return tyPlus;
293}
294
295
297(
298 Ostream& os
299) const
300{
302 writeLocalEntries(os);
305
306
307// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308
309namespace Foam
310{
312 (
314 nutkWallFunctionFvPatchScalarField
315 );
316}
317
318
319// ************************************************************************* //
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:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
A FieldMapper for finite-volume patch fields.
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
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.
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.
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
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
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
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)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299