nutUSpaldingWallFunctionFvPatchScalarField.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-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"
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 magGradU(mag(Uw.snGrad()));
54
55 const tmp<scalarField> tnuw = turbModel.nu(patchi);
56 const scalarField& nuw = tnuw();
57
58
59 // Calculate new viscosity
61 (
62 max
63 (
64 scalar(0),
65 sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
66 )
67 );
68
69 if (tolerance_ != 0.01)
70 {
71 // User-specified tolerance. Find out if current nut already satisfies
72 // eqns.
73
74 // Run ut for one iteration
75 scalarField err;
76 tmp<scalarField> UTau(calcUTau(magGradU, 1, err));
77
78 // Preserve nutw if the initial error (err) already lower than the
79 // tolerance.
80
81 scalarField& nutw = tnutw.ref();
82 forAll(err, facei)
83 {
84 if (err[facei] < tolerance_)
85 {
86 nutw[facei] = this->operator[](facei);
87 }
88 }
89 }
90 return tnutw;
91}
92
93
96(
97 const scalarField& magGradU
98) const
99{
100 scalarField err;
101 return calcUTau(magGradU, maxIter_, err);
102}
103
104
107(
108 const scalarField& magGradU,
109 const label maxIter,
110 scalarField& err
111) const
112{
113 const label patchi = patch().index();
114
115 const auto& turbModel = db().lookupObject<turbulenceModel>
116 (
118 (
120 internalField().group()
121 )
122 );
123
124 const scalarField& y = turbModel.y()[patchi];
125
126 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
127 const scalarField magUp(mag(Uw.patchInternalField() - Uw));
128
129 const tmp<scalarField> tnuw = turbModel.nu(patchi);
130 const scalarField& nuw = tnuw();
131
132 const scalar kappa = wallCoeffs_.kappa();
133 const scalar E = wallCoeffs_.E();
134
135 const scalarField& nutw = *this;
136
137 auto tuTau = tmp<scalarField>::New(patch().size(), Zero);
138 auto& uTau = tuTau.ref();
139
140 err.setSize(uTau.size());
141 err = 0.0;
142
143 forAll(uTau, facei)
144 {
145 scalar ut = sqrt((nutw[facei] + nuw[facei])*magGradU[facei]);
146 // Note: for exact restart seed with laminar viscosity only:
147 //scalar ut = sqrt(nuw[facei]*magGradU[facei]);
148
149 if (ROOTVSMALL < ut)
150 {
151 int iter = 0;
152
153 do
154 {
155 const scalar kUu = min(kappa*magUp[facei]/ut, scalar(50));
156 const scalar fkUu = exp(kUu) - 1 - kUu*(1 + 0.5*kUu);
157
158 const scalar f =
159 - ut*y[facei]/nuw[facei]
160 + magUp[facei]/ut
161 + 1.0/E*(fkUu - 1.0/6.0*kUu*sqr(kUu));
162
163 const scalar df =
164 y[facei]/nuw[facei]
165 + magUp[facei]/sqr(ut)
166 + 1.0/E*kUu*fkUu/ut;
167
168 const scalar uTauNew = ut + f/df;
169 err[facei] = mag((ut - uTauNew)/ut);
170 ut = uTauNew;
171
172 //iterations_++;
173
174 } while
175 (
176 ut > ROOTVSMALL
177 && err[facei] > tolerance_
178 && ++iter < maxIter
179 );
180
181 uTau[facei] = max(scalar(0), ut);
182
183 //invocations_++;
184 //if (iter > 1)
185 //{
186 // nontrivial_++;
187 //}
188 //if (iter >= maxIter_)
189 //{
190 // nonconvergence_++;
191 //}
192 }
193 }
194
195 return tuTau;
196}
197
198
200(
201 Ostream& os
202) const
203{
204 os.writeEntryIfDifferent<label>("maxIter", 10, maxIter_);
205 os.writeEntryIfDifferent<scalar>("tolerance", 0.01, tolerance_);
206}
207
208
209// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
210
211Foam::nutUSpaldingWallFunctionFvPatchScalarField::
212nutUSpaldingWallFunctionFvPatchScalarField
213(
214 const fvPatch& p,
216)
217:
219 maxIter_(10),
220 tolerance_(0.01)
221 //invocations_(0),
222 //nontrivial_(0),
223 //nonconvergence_(0),
224 //iterations_(0)
225{}
226
227
228Foam::nutUSpaldingWallFunctionFvPatchScalarField::
229nutUSpaldingWallFunctionFvPatchScalarField
230(
232 const fvPatch& p,
234 const fvPatchFieldMapper& mapper
235)
236:
237 nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
238 maxIter_(ptf.maxIter_),
239 tolerance_(ptf.tolerance_)
240 //invocations_(0),
241 //nontrivial_(0),
242 //nonconvergence_(0),
243 //iterations_(0)
244{}
245
246
247Foam::nutUSpaldingWallFunctionFvPatchScalarField::
248nutUSpaldingWallFunctionFvPatchScalarField
249(
250 const fvPatch& p,
252 const dictionary& dict
253)
254:
256 maxIter_(dict.getOrDefault<label>("maxIter", 10)),
257 tolerance_(dict.getOrDefault<scalar>("tolerance", 0.01))
258 //invocations_(0),
259 //nontrivial_(0),
260 //nonconvergence_(0),
261 //iterations_(0)
262{}
263
264
265Foam::nutUSpaldingWallFunctionFvPatchScalarField::
266nutUSpaldingWallFunctionFvPatchScalarField
267(
269)
270:
272 maxIter_(wfpsf.maxIter_),
273 tolerance_(wfpsf.tolerance_)
274 //invocations_(wfpsf.invocations_),
275 //nontrivial_(wfpsf.nontrivial_),
276 //nonconvergence_(wfpsf.nonconvergence_),
277 //iterations_(wfpsf.iterations_)
278{}
279
280
281Foam::nutUSpaldingWallFunctionFvPatchScalarField::
282nutUSpaldingWallFunctionFvPatchScalarField
283(
286)
287:
289 maxIter_(wfpsf.maxIter_),
290 tolerance_(wfpsf.tolerance_)
291 //invocations_(0),
292 //nontrivial_(0),
293 //nonconvergence_(0),
294 //iterations_(0)
295{}
296
297
298// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
299
302{
303 //if (debug)
304 //{
305 // Info<< "nutUSpaldingWallFunctionFvPatchScalarField :"
306 // << " total invocations:"
307 // << returnReduce(invocations_, sumOp<label>())
308 // << " total iterations:"
309 // << returnReduce(iterations_, sumOp<label>())
310 // << " total non-convergence:"
311 // << returnReduce(nonconvergence_, sumOp<label>())
312 // << " total non-trivial:"
313 // << returnReduce(nontrivial_, sumOp<label>())
314 // << endl;
315 //}
316}
317
318
319// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
320
323{
324 const label patchi = patch().index();
325
326 const auto& turbModel = db().lookupObject<turbulenceModel>
327 (
329 (
331 internalField().group()
332 )
333 );
334
335 const scalarField& y = turbModel.y()[patchi];
336
337 const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
338
339 const tmp<scalarField> tnuw = turbModel.nu(patchi);
340 const scalarField& nuw = tnuw();
341
342 return y*calcUTau(mag(Uw.snGrad()))/nuw;
343}
344
345
347(
348 Ostream& os
349) const
350{
352 writeLocalEntries(os);
353 writeEntry("value", os);
354}
355
356
357// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
358
359namespace Foam
360{
362 (
365 );
366}
367
368
369// ************************************************************************* //
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.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:251
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 > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:237
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
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.
tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
virtual tmp< scalarField > calcNut() const
Uncomment in case of intrumentation.
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.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
U
Definition: pEqn.H:72
volScalarField & p
scalar magUp
scalar uTau
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)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333