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-2020 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 turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
44  (
46  (
48  internalField().group()
49  )
50  );
51  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
52  const scalarField magGradU(mag(Uw.snGrad()));
53  const tmp<scalarField> tnuw = turbModel.nu(patchi);
54  const scalarField& nuw = tnuw();
55 
56 
57  // Calculate new viscosity
58  tmp<scalarField> tnutw
59  (
60  max
61  (
62  scalar(0),
63  sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
64  )
65  );
66 
67  if (tolerance_ != 0.01)
68  {
69  // User-specified tolerance. Find out if current nut already satisfies
70  // eqns.
71 
72  // Run ut for one iteration
73  scalarField err;
74  tmp<scalarField> UTau(calcUTau(magGradU, 1, err));
75 
76  // Preserve nutw if the initial error (err) already lower than the
77  // tolerance.
78 
79  scalarField& nutw = tnutw.ref();
80  forAll(err, facei)
81  {
82  if (err[facei] < tolerance_)
83  {
84  nutw[facei] = this->operator[](facei);
85  }
86  }
87  }
88  return tnutw;
89 }
90 
91 
94 (
95  const scalarField& magGradU
96 ) const
97 {
98  scalarField err;
99  return calcUTau(magGradU, maxIter_, err);
100 }
101 
102 
105 (
106  const scalarField& magGradU,
107  const label maxIter,
108  scalarField& err
109 ) const
110 {
111  const label patchi = patch().index();
112 
113  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
114  (
116  (
118  internalField().group()
119  )
120  );
121  const scalarField& y = turbModel.y()[patchi];
122 
123  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
124  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
125 
126  const tmp<scalarField> tnuw = turbModel.nu(patchi);
127  const scalarField& nuw = tnuw();
128 
129  const scalarField& nutw = *this;
130 
131  tmp<scalarField> tuTau(new scalarField(patch().size(), Zero));
132  scalarField& uTau = tuTau.ref();
133 
134  err.setSize(uTau.size());
135  err = 0.0;
136 
137  forAll(uTau, facei)
138  {
139  scalar ut = sqrt((nutw[facei] + nuw[facei])*magGradU[facei]);
140  // Note: for exact restart seed with laminar viscosity only:
141  //scalar ut = sqrt(nuw[facei]*magGradU[facei]);
142 
143  if (ROOTVSMALL < ut)
144  {
145  int iter = 0;
146 
147  do
148  {
149  scalar kUu = min(kappa_*magUp[facei]/ut, 50);
150  scalar fkUu = exp(kUu) - 1 - kUu*(1 + 0.5*kUu);
151 
152  scalar f =
153  - ut*y[facei]/nuw[facei]
154  + magUp[facei]/ut
155  + 1/E_*(fkUu - 1.0/6.0*kUu*sqr(kUu));
156 
157  scalar df =
158  y[facei]/nuw[facei]
159  + magUp[facei]/sqr(ut)
160  + 1/E_*kUu*fkUu/ut;
161 
162  scalar uTauNew = ut + f/df;
163  err[facei] = mag((ut - uTauNew)/ut);
164  ut = uTauNew;
165 
166  //iterations_++;
167 
168  } while
169  (
170  ut > ROOTVSMALL
171  && err[facei] > tolerance_
172  && ++iter < maxIter
173  );
174 
175  uTau[facei] = max(0.0, ut);
176 
177  //invocations_++;
178  //if (iter > 1)
179  //{
180  // nontrivial_++;
181  //}
182  //if (iter >= maxIter_)
183  //{
184  // nonconvergence_++;
185  //}
186  }
187  }
188 
189  return tuTau;
190 }
191 
192 
194 (
195  Ostream& os
196 ) const
197 {
199 
200  os.writeEntryIfDifferent<label>("maxIter", 10, maxIter_);
201  os.writeEntryIfDifferent<scalar>("tolerance", 0.01, tolerance_);
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
206 
209 (
210  const fvPatch& p,
212 )
213 :
215  maxIter_(10),
216  tolerance_(0.01)
217  //invocations_(0),
218  //nontrivial_(0),
219  //nonconvergence_(0),
220  //iterations_(0)
221 {}
222 
223 
226 (
228  const fvPatch& p,
230  const fvPatchFieldMapper& mapper
231 )
232 :
233  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
234  maxIter_(ptf.maxIter_),
235  tolerance_(ptf.tolerance_)
236  //invocations_(0),
237  //nontrivial_(0),
238  //nonconvergence_(0),
239  //iterations_(0)
240 {}
241 
242 
245 (
246  const fvPatch& p,
248  const dictionary& dict
249 )
250 :
252  maxIter_(dict.getOrDefault<label>("maxIter", 10)),
253  tolerance_(dict.getOrDefault<scalar>("tolerance", 0.01))
254  //invocations_(0),
255  //nontrivial_(0),
256  //nonconvergence_(0),
257  //iterations_(0)
258 {}
259 
260 
263 (
265 )
266 :
268  maxIter_(wfpsf.maxIter_),
269  tolerance_(wfpsf.tolerance_)
270  //invocations_(wfpsf.invocations_),
271  //nontrivial_(wfpsf.nontrivial_),
272  //nonconvergence_(wfpsf.nonconvergence_),
273  //iterations_(wfpsf.iterations_)
274 {}
275 
276 
279 (
282 )
283 :
285  maxIter_(wfpsf.maxIter_),
286  tolerance_(wfpsf.tolerance_)
287  //invocations_(0),
288  //nontrivial_(0),
289  //nonconvergence_(0),
290  //iterations_(0)
291 {}
292 
293 
294 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
295 
298 {
299  //if (debug)
300  //{
301  // Info<< "nutUSpaldingWallFunctionFvPatchScalarField :"
302  // << " total invocations:"
303  // << returnReduce(invocations_, sumOp<label>())
304  // << " total iterations:"
305  // << returnReduce(iterations_, sumOp<label>())
306  // << " total non-convergence:"
307  // << returnReduce(nonconvergence_, sumOp<label>())
308  // << " total non-trivial:"
309  // << returnReduce(nontrivial_, sumOp<label>())
310  // << endl;
311  //}
312 }
313 
314 
315 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
316 
319 {
320  const label patchi = patch().index();
321 
322  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
323  (
325  (
327  internalField().group()
328  )
329  );
330  const scalarField& y = turbModel.y()[patchi];
331  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
332  const tmp<scalarField> tnuw = turbModel.nu(patchi);
333  const scalarField& nuw = tnuw();
334 
335  return y*calcUTau(mag(Uw.snGrad()))/nuw;
336 }
337 
338 
340 (
341  Ostream& os
342 ) const
343 {
345  writeLocalEntries(os);
346  writeEntry("value", os);
347 }
348 
349 
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351 
352 namespace Foam
353 {
355  (
358  );
359 }
360 
361 
362 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::nutUSpaldingWallFunctionFvPatchScalarField::tolerance_
const scalar tolerance_
Convergence tolerance.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.H:157
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:217
Foam::Ostream::writeEntryIfDifferent
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:244
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::nutUSpaldingWallFunctionFvPatchScalarField::~nutUSpaldingWallFunctionFvPatchScalarField
virtual ~nutUSpaldingWallFunctionFvPatchScalarField()
Destructor.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:297
Foam::constant::atomic::group
constexpr const char *const group
Group name for atomic constants.
Definition: atomicConstants.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::nutWallFunctionFvPatchScalarField::U
virtual const volVectorField & U(const turbulenceModel &turb) const
Definition: nutWallFunctionFvPatchScalarField.C:75
fvPatchFieldMapper.H
Foam::nutUSpaldingWallFunctionFvPatchScalarField::maxIter_
const label maxIter_
Max iterations in calcNut.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.H:154
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
magUp
scalar magUp
Definition: evaluateNearWall.H:10
Foam::nutUSpaldingWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:318
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::nutUSpaldingWallFunctionFvPatchScalarField::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:194
Foam::nutUSpaldingWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Uncomment in case of intrumentation.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:39
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::Field< scalar >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:225
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
uTau
scalar uTau
Definition: evaluateNearWall.H:14
Foam::nutUSpaldingWallFunctionFvPatchScalarField::write
virtual void write(Ostream &os) const
Write.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:340
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::nutUSpaldingWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutUSpaldingWallFunctionFvPatchScalarField.H:145
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::nutUSpaldingWallFunctionFvPatchScalarField::calcUTau
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Calculate the friction velocity.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:94
U
U
Definition: pEqn.H:72
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::nutUSpaldingWallFunctionFvPatchScalarField::nutUSpaldingWallFunctionFvPatchScalarField
nutUSpaldingWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: nutUSpaldingWallFunctionFvPatchScalarField.C:209
Foam::nutWallFunctionFvPatchScalarField::nutw
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
Definition: nutWallFunctionFvPatchScalarField.C:235
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
nutUSpaldingWallFunctionFvPatchScalarField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
turbulenceModel.H
Foam::nutWallFunctionFvPatchScalarField::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Definition: nutWallFunctionFvPatchScalarField.C:91
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
y
scalar y
Definition: LISASMDCalcMethod1.H:14